Flow123d
intersection.cc
Go to the documentation of this file.
1 /*
2  * intersection.cc
3  *
4  * Created on: May 24, 2011
5  * Author: jb
6  */
7 
8 #include "mesh/intersection.hh"
9 #include "mesh/mesh.h"
10 
11 #include <boost/tokenizer.hpp>
12 #include "boost/lexical_cast.hpp"
13 #include <armadillo>
14 
15 // inicializovat objekt, cist zbytek tokenu z tok a naplnit map a shift pro master a slave
16 // viz dokumentace k Armadillu
18  const IntersectionLocal *isec)
19 : dim(isec->n_points() - 1),
20  master(ele_master), slave(ele_slave),
21  master_map(master->dim(), dim), slave_map(slave->dim(), dim),
22  master_shift(master->dim()), slave_shift(slave->dim())
23 {
24  ///otestuje se jestli dimenze masteru je mensi nez dimenze slave - chybova hlaska (vyjimka - throw)
25  ///pocet pointu=dim+1
26  if (master->dim() > slave->dim()) {
27  cout << "Exception: master->dim() > slave->dim()" << endl;
28  //throw((char*) "master->dim > slave->dim");
29  }
30 
32 
33  arma::vec master_tmp(master_shift), slave_tmp(slave_shift);
34  // cyklus pres body pruniku
35  for (unsigned int i = 1; i < (dim + 1); ++i) {
36  intersection_point_to_vectors(isec->get_point(i),master_tmp, slave_tmp);
37  master_tmp -= master_shift;
38  slave_tmp -= slave_shift;
39 
40  master_map.col(i-1) = master_tmp;
41  slave_map.col(i-1) = slave_tmp;
42  }
43 }
44 
45 
46 
48  {return master->dim();}
49 
50 
51 
53  {return slave->dim();}
54 
55 
56 
57 void Intersection::intersection_point_to_vectors(const IntersectionPoint *point, arma::vec &vec1, arma::vec &vec2)
58 {
59  const vector<double> &coord_el1 = point->el1_coord();
60  ASSERT_EQUAL(coord_el1.size() , vec1.n_elem);
61  vec1=arma::vec(coord_el1);
62 
63  const vector<double> &coord_el2 = point->el2_coord();
64  ASSERT_EQUAL(coord_el2.size() , vec2.n_elem);
65  vec2=arma::vec(coord_el2);
66 }
67 
68 
69 arma::vec Intersection::map_to_master(const arma::vec &point) const
70 {
71  //dim = dimenze intersec elementu
72  ASSERT(( point.n_elem == dim ),"Map to slave: point.n_elem(%d) != dim(%d) \n", point.n_elem, dim);
73  int result_dim = master->dim();
74  arma::vec result(result_dim+1);
75  result(0)=1.0;
76  result.subvec(1, result_dim) = (master_map * point + master_shift);
77  return result;
78 }
79 
80 arma::vec Intersection::map_to_slave(const arma::vec &point) const
81 {
82  ASSERT(( point.n_elem == dim ),"Map to slave: point.n_elem(%d) != dim(%d) \n", point.n_elem, dim);
83  int result_dim = slave->dim();
84  arma::vec result(result_dim+1);
85  result(0)=1.0;
86  //DBGMSG("s dim: %d dim:%d\n", master->dim, )
87  result.subvec(1, result_dim) = (slave_map * point + slave_shift);
88  return result;
89 }
90 
92 
93  static const double factorial[4] = {1.0, 1.0, 2.0, 6.0};
94  return (master->measure() * det(master_map) / factorial[dim]);
95 }
96 
97