Flow123d  release_3.0.0-859-g84487d0
intersectionquadrature.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file intersection.cc
15  * @brief
16  */
17 
18 #include "mesh/intersection.hh"
19 #include "mesh/mesh.h"
20 #include "mesh/accessors.hh"
21 
22 #include <boost/tokenizer.hpp>
23 #include "boost/lexical_cast.hpp"
24 #include <armadillo>
25 
26 // inicializovat objekt, cist zbytek tokenu z tok a naplnit map a shift pro master a slave
27 // viz dokumentace k Armadillu
29  const IntersectionLocal *isec)
30 : dim(isec->n_points() - 1),
31  master(ele_master), slave(ele_slave),
32  master_map(master->dim(), dim), slave_map(slave->dim(), dim),
33  master_shift(master->dim()), slave_shift(slave->dim())
34 {
35  ///otestuje se jestli dimenze masteru je mensi nez dimenze slave - chybova hlaska (vyjimka - throw)
36  ///pocet pointu=dim+1
37  if (master->dim() > slave->dim()) {
38  cout << "Exception: master->dim() > slave->dim()" << endl;
39  //throw((char*) "master->dim > slave->dim");
40  }
41 
43 
44  arma::vec master_tmp(master_shift), slave_tmp(slave_shift);
45  // cyklus pres body pruniku
46  for (unsigned int i = 1; i < (dim + 1); ++i) {
47  intersection_point_to_vectors(isec->get_point(i),master_tmp, slave_tmp);
48  master_tmp -= master_shift;
49  slave_tmp -= slave_shift;
50 
51  master_map.col(i-1) = master_tmp;
52  slave_map.col(i-1) = slave_tmp;
53  }
54 }
55 
56 
57 
59  {return master->dim();}
60 
61 
62 
64  {return slave->dim();}
65 
66 
67 
68 void Intersection::intersection_point_to_vectors(const IntersectionPoint *point, arma::vec &vec1, arma::vec &vec2)
69 {
70  const vector<double> &coord_el1 = point->el1_coord();
71  OLD_ASSERT_EQUAL(coord_el1.size() , vec1.n_elem);
72  vec1=arma::vec(coord_el1);
73 
74  const vector<double> &coord_el2 = point->el2_coord();
75  OLD_ASSERT_EQUAL(coord_el2.size() , vec2.n_elem);
76  vec2=arma::vec(coord_el2);
77 }
78 
79 
80 arma::vec Intersection::map_to_master(const arma::vec &point) const
81 {
82  //dim = dimenze intersec elementu
83  OLD_ASSERT(( point.n_elem == dim ),"Map to slave: point.n_elem(%d) != dim(%d) \n", point.n_elem, dim);
84  int result_dim = master->dim();
85  arma::vec result(result_dim+1);
86  result(0)=1.0;
87  result.subvec(1, result_dim) = (master_map * point + master_shift);
88  return result;
89 }
90 
91 arma::vec Intersection::map_to_slave(const arma::vec &point) const
92 {
93  OLD_ASSERT(( point.n_elem == dim ),"Map to slave: point.n_elem(%d) != dim(%d) \n", point.n_elem, dim);
94  int result_dim = slave->dim();
95  arma::vec result(result_dim+1);
96  result(0)=1.0;
97  result.subvec(1, result_dim) = (slave_map * point + slave_shift);
98  return result;
99 }
100 
102 
103  static const double factorial[4] = {1.0, 1.0, 2.0, 6.0};
104  return (master.measure() * det(master_map) / factorial[dim]);
105 }
106 
107 
arma::vec map_to_master(const arma::vec &point) const
ElementAccessor< 3 > slave
arma::vec map_to_slave(const arma::vec &point) const
Class represents intersection of two elements.
Class represents an intersection point of simplex<N> and simplex<M>. It contains barycentric coordina...
ElementAccessor< 3 > master
unsigned int dim
dimenze pruniku
unsigned int dim() const
Definition: elements.h:124
#define OLD_ASSERT(...)
Definition: global_defs.h:131
unsigned int master_dim()
dimension of the master element
arma::Mat< double > master_map
matrix part of linear transform from reference element of intersection to reference element of master...
unsigned int slave_dim()
dimension of the slave element
arma::vec master_shift
shift vector of the linear transform
void intersection_point_to_vectors(const IntersectionPoint *point, arma::vec &vec1, arma::vec &vec2)
double measure() const
Computes the measure of the element.
Definition: accessors.hh:254
Intersection(ElementAccessor< 3 > ele_master, ElementAccessor< 3 > ele_slave, const IntersectionLocal *isec)
arma::Mat< double > slave_map
double intersection_true_size() const
#define OLD_ASSERT_EQUAL(a, b)
Definition: global_defs.h:133