Flow123d  build_with_4.0.3-95067a1
mapping_p1.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 mapping_p1.cc
15  * @brief Class MappingP1 implements the affine transformation of
16  * the unit cell onto the actual cell.
17  * @author Jan Stebel
18  */
19 
20 #include "fem/mapping_p1.hh"
21 #include "quadrature/quadrature.hh"
22 
23 
24 
25 using namespace std;
26 
27 
28 
29 template<unsigned int dim, unsigned int spacedim>
31 {
32  UpdateFlags f = flags;
33 
34  if (flags & update_normal_vectors)
36 
37  if ((flags & update_volume_elements) |
38  (flags & update_JxW_values) |
39  (flags & update_inverse_jacobians))
40  f |= update_jacobians;
41 
42  return f;
43 }
44 
45 
46 
47 template<unsigned int dim, unsigned int spacedim>
49 {
50  ElementMap coords;
51  for (unsigned int i=0; i<dim+1; i++)
52  coords.col(i) = *elm.node(i);
53  return coords;
54 }
55 
56 
57 template<unsigned int dim, unsigned int spacedim>
58 auto MappingP1<dim,spacedim>::jacobian(const ElementMap &coords) -> arma::mat::fixed<spacedim,dim>
59 {
60  arma::mat::fixed<spacedim,dim> jac;
61  for (unsigned int i=0; i<spacedim; i++)
62  for (unsigned int j=0; j<dim; j++)
63  jac(i,j) = coords(i,j+1) - coords(i,0);
64  return jac;
65 }
66 
67 
68 template<unsigned int dim, unsigned int spacedim>
70 {
71  arma::mat::fixed<3, dim> A = map.cols(1,dim);
72  for(unsigned int i=0; i < dim; i++ ) {
73  A.col(i) -= map.col(0);
74  }
75 
76  arma::mat::fixed<dim, dim> AtA = A.t()*A;
77  arma::vec::fixed<dim> Atb = A.t()*(point - map.col(0));
78  arma::vec::fixed<dim+1> bary_coord;
79  bary_coord.rows(1, dim) = arma::solve(AtA, Atb);
80  bary_coord( 0 ) = 1.0 - arma::sum( bary_coord.rows(1,dim) );
81  return bary_coord;
82 }
83 
84 template<unsigned int dim, unsigned int spacedim>
86 {
87  return map * point;
88 }
89 
90 template<unsigned int dim, unsigned int spacedim>
92  return RefElement<dim>::clip(barycentric);
93 }
94 
95 template <unsigned int dim, unsigned int spacedim>
97 {
98  arma::vec projection = project_real_to_unit(point, element_map(elm));
99  return (projection.min() >= -BoundingBox::epsilon);
100 }
101 
102 
103 
104 template class MappingP1<0,3>; // Only for compilation of DG transport assemble methods, do not use this instance!
105 template class MappingP1<1,3>;
106 template class MappingP1<2,3>;
107 template class MappingP1<3,3>;
108 
109 
110 
111 
static const double epsilon
stabilization parameter
Definition: bounding_box.hh:65
Affine mapping between reference and actual cell.
Definition: mapping_p1.hh:59
arma::vec::fixed< dim+1 > BaryPoint
Definition: mapping_p1.hh:62
arma::vec::fixed< spacedim > RealPoint
Definition: mapping_p1.hh:63
static BaryPoint clip_to_element(BaryPoint &barycentric)
Definition: mapping_p1.cc:91
static arma::mat::fixed< spacedim, dim > jacobian(const ElementMap &coords)
Definition: mapping_p1.cc:58
static bool contains_point(arma::vec point, ElementAccessor< 3 > elm)
Test if element contains given point.
Definition: mapping_p1.cc:96
static RealPoint project_unit_to_real(const BaryPoint &point, const ElementMap &map)
Definition: mapping_p1.cc:85
arma::mat::fixed< spacedim, dim+1 > ElementMap
Definition: mapping_p1.hh:64
static BaryPoint project_real_to_unit(const RealPoint &point, const ElementMap &map)
Definition: mapping_p1.cc:69
static ElementMap element_map(ElementAccessor< 3 > elm)
Definition: mapping_p1.cc:48
static UpdateFlags update_each(UpdateFlags flags)
Determines which additional quantities have to be computed.
Definition: mapping_p1.cc:30
static BaryPoint clip(const BaryPoint &barycentric)
Definition: ref_element.cc:332
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell.
ArmaVec< double, N > vec
Definition: armor.hh:885
Basic definitions of numerical quadrature rules.
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
Definition: update_flags.hh:68
@ update_volume_elements
Determinant of the Jacobian.
@ update_normal_vectors
Normal vectors.
@ update_JxW_values
Transformed quadrature weights.
@ update_jacobians
Volume element.
@ update_inverse_jacobians
Volume element.