Flow123d  jenkins-Flow123d-linux-release-multijob-198
field_fe.impl.hh
Go to the documentation of this file.
1 /*
2  * field_constant.impl.hh
3  *
4  * Created on: Dec 15, 2012
5  * Author: jb
6  */
7 
8 #ifndef FIELD_FE_IMPL_HH_
9 #define FIELD_FE_IMPL_HH_
10 
11 #include "fields/field_fe.hh"
12 #include "input/input_type.hh"
13 #include "quadrature/quadrature.hh"
14 #include "fem/fe_values.hh"
15 #include "fem/finite_element.hh"
16 
17 
18 
19 
20 /// Implementation.
21 
22 namespace it = Input::Type;
23 
24 
25 
26 
27 template <int spacedim, class Value>
28 FieldFE<spacedim, Value>::FieldFE( unsigned int n_comp)
29 : FieldAlgorithmBase<spacedim, Value>(n_comp),
30  dh_(nullptr),
31  data_(nullptr),
32  data_vec_(nullptr),
33  dof_indices(nullptr),
34  map1_(nullptr),
35  map2_(nullptr),
36  map3_(nullptr)
37 {}
38 
39 
40 template <int spacedim, class Value>
42  Mapping<1,3> *map1,
43  Mapping<2,3> *map2,
44  Mapping<3,3> *map3,
45  const Vec *data)
46 {
47  dh_ = dh;
48  map1_ = map1;
49  map2_ = map2;
50  map3_ = map3;
51  data_vec_ = data;
52  VecGetArray(*data_vec_, &data_);
53 
54  unsigned int ndofs = max(dh_->fe<1>()->n_dofs(), max(dh_->fe<2>()->n_dofs(), dh_->fe<3>()->n_dofs()));
55  dof_indices = new unsigned int[ndofs];
56 }
57 
58 
59 
60 /**
61  * Returns one value in one given point. ResultType can be used to avoid some costly calculation if the result is trivial.
62  */
63 template <int spacedim, class Value>
64 typename Value::return_type const & FieldFE<spacedim, Value>::value(const Point &p, const ElementAccessor<spacedim> &elm)
65 {
66  Point p_rel = p - elm.element()->node[0]->point();
67  DOFHandlerBase::CellIterator cell = dh_->mesh()->element(elm.idx());
68 
69  if (elm.dim() == 1) {
70  arma::mat::fixed<3,1> m1 = elm.element()->node[1]->point() - elm.element()->node[0]->point();
71  arma::mat::fixed<1,3> im1 = pinv(m1);
72 
73  Quadrature<1> q1(1);
74  q1.set_point(0, im1*p_rel);
75 
76  FEValues<1,3> fe_values1(*map1_, q1, *dh_->fe<1>(), update_values);
77  fe_values1.reinit(cell);
78 
79  dh_->get_dof_indices(cell, dof_indices);
80 
81  if (dh_->fe<1>()->is_scalar()) {
82  double value = 0;
83  for (int i=0; i<dh_->fe<1>()->n_dofs(); i++)
84  value += data_[dof_indices[i]]*fe_values1.shape_value(i, 0);
85  this->value_(0,0) = value;
86  }
87  else {
88  arma::vec3 value;
89  value.zeros();
90  for (int i=0; i<dh_->fe<1>()->n_dofs(); i++)
91  value += data_[dof_indices[i]]*fe_values1.shape_vector(i, 0);
92  for (int i=0; i<3; i++)
93  this->value_(i,0) = value(i);
94  }
95  }
96  else if (elm.dim() == 2) {
97  arma::mat::fixed<3,2> m2;
98  m2.col(0) = elm.element()->node[1]->point() - elm.element()->node[0]->point();
99  m2.col(1) = elm.element()->node[2]->point() - elm.element()->node[0]->point();
100  arma::mat::fixed<2,3> im2 = pinv(m2);
101 
102  Quadrature<2> q2(1);
103  q2.set_point(0, im2*p_rel);
104 
105  FEValues<2,3> fe_values2(*map2_, q2, *dh_->fe<2>(), update_values);
106  fe_values2.reinit(cell);
107 
108  dh_->get_dof_indices(cell, dof_indices);
109 
110  if (dh_->fe<2>()->is_scalar()) {
111  double value = 0;
112  for (int i=0; i<dh_->fe<2>()->n_dofs(); i++)
113  value += data_[dof_indices[i]]*fe_values2.shape_value(i, 0);
114  this->value_(0,0) = value;
115  }
116  else {
117  arma::vec3 value;
118  value.zeros();
119  for (int i=0; i<dh_->fe<2>()->n_dofs(); i++)
120  value += data_[dof_indices[i]]*fe_values2.shape_vector(i, 0);
121  for (int i=0; i<3; i++)
122  this->value_(i,0) = value(i);
123  } }
124  else {
125  arma::mat33 m3;
126  m3.col(0) = elm.element()->node[1]->point() - elm.element()->node[0]->point();
127  m3.col(1) = elm.element()->node[2]->point() - elm.element()->node[0]->point();
128  m3.col(2) = elm.element()->node[3]->point() - elm.element()->node[0]->point();
129  arma::mat33 im3 = inv(m3);
130 
131  Quadrature<3> q3(1);
132  q3.set_point(0, im3*p_rel);
133 
134  FEValues<3,3> fe_values3(*map3_, q3, *dh_->fe<3>(), update_values);
135  fe_values3.reinit(cell);
136 
137  dh_->get_dof_indices(cell, dof_indices);
138 
139  if (dh_->fe<3>()->is_scalar()) {
140  double value = 0;
141  for (int i=0; i<dh_->fe<3>()->n_dofs(); i++)
142  value += data_[dof_indices[i]]*fe_values3.shape_value(i, 0);
143  this->value_(0,0) = value;
144  }
145  else {
146  arma::vec3 value;
147  value.zeros();
148  for (int i=0; i<dh_->fe<3>()->n_dofs(); i++)
149  value += data_[dof_indices[i]]*fe_values3.shape_vector(i, 0);
150  for (int i=0; i<3; i++)
151  this->value_(i,0) = value(i);
152  }
153  }
154 
155  return this->r_value_;
156 }
157 
158 
159 
160 /**
161  * Returns std::vector of scalar values in several points at once.
162  */
163 template <int spacedim, class Value>
166 {
167  ASSERT_EQUAL( point_list.size(), value_list.size() );
168 
169  DOFHandlerBase::CellIterator cell = dh_->mesh()->element( elm.idx() );
170 
171  if (elm.dim() == 1) {
172  arma::mat::fixed<3,1> m1 = elm.element()->node[1]->point() - elm.element()->node[0]->point();
173  arma::mat::fixed<1,3> im1 = pinv(m1);
174 
175  dh_->get_dof_indices(cell, dof_indices);
176 
177  for (int k=0; k<point_list.size(); k++) {
178  Quadrature<1> q1(1);
179  Point p_rel = point_list[k] - elm.element()->node[0]->point();
180  q1.set_point(0, im1*p_rel);
181 
182  FEValues<1,3> fe_values1(*map1_, q1, *dh_->fe<1>(), update_values);
183  fe_values1.reinit(cell);
184 
185  Value envelope(value_list[k]);
186 
187  if (dh_->fe<1>()->is_scalar()) {
188  double value = 0;
189  for (int i=0; i<dh_->fe<1>()->n_dofs(); i++)
190  value += data_[dof_indices[i]]*fe_values1.shape_value(i, 0);
191  envelope(0,0) = value;
192  }
193  else {
194  arma::vec3 value;
195  value.zeros();
196  for (int i=0; i<dh_->fe<1>()->n_dofs(); i++)
197  value += data_[dof_indices[i]]*fe_values1.shape_vector(i, 0);
198  for (int i=0; i<3; i++)
199  envelope(i,0) = value(i);
200  }
201  }
202  }
203  else if (elm.dim() == 2) {
204  arma::mat::fixed<3,2> m2;
205  m2.col(0) = elm.element()->node[1]->point() - elm.element()->node[0]->point();
206  m2.col(1) = elm.element()->node[2]->point() - elm.element()->node[0]->point();
207  arma::mat::fixed<2,3> im2 = pinv(m2);
208 
209  dh_->get_dof_indices(cell, dof_indices);
210 
211  for (int k=0; k<point_list.size(); k++) {
212  Quadrature<2> q2(1);
213  Point p_rel = point_list[k] - elm.element()->node[0]->point();
214  q2.set_point(0, im2*p_rel);
215 
216  FEValues<2,3> fe_values2(*map2_, q2, *dh_->fe<2>(), update_values);
217  fe_values2.reinit(cell);
218 
219  Value envelope(value_list[k]);
220 
221  if (dh_->fe<2>()->is_scalar()) {
222  double value = 0;
223  for (int i=0; i<dh_->fe<2>()->n_dofs(); i++)
224  value += data_[dof_indices[i]]*fe_values2.shape_value(i, 0);
225  envelope(0,0) = value;
226  }
227  else {
228  arma::vec3 value;
229  value.zeros();
230  for (int i=0; i<dh_->fe<2>()->n_dofs(); i++)
231  value += data_[dof_indices[i]]*fe_values2.shape_vector(i, 0);
232  for (int i=0; i<3; i++)
233  envelope(i,0) = value(i);
234  }
235  }
236  }
237  else {
238  arma::mat33 m3;
239  m3.col(0) = elm.element()->node[1]->point() - elm.element()->node[0]->point();
240  m3.col(1) = elm.element()->node[2]->point() - elm.element()->node[0]->point();
241  m3.col(2) = elm.element()->node[3]->point() - elm.element()->node[0]->point();
242  arma::mat33 im3 = inv(m3);
243 
244  dh_->get_dof_indices(cell, dof_indices);
245 
246  for (int k=0; k<point_list.size(); k++) {
247  Quadrature<3> q3(1);
248  Point p_rel = point_list[k] - elm.element()->node[0]->point();
249  q3.set_point(0, im3*p_rel);
250 
251  FEValues<3,3> fe_values3(*map3_, q3, *dh_->fe<3>(), update_values);
252  fe_values3.reinit(cell);
253 
254  Value envelope(value_list[k]);
255 
256  if (dh_->fe<3>()->is_scalar()) {
257  double value = 0;
258  for (int i=0; i<dh_->fe<3>()->n_dofs(); i++)
259  value += data_[dof_indices[i]]*fe_values3.shape_value(i, 0);
260  envelope(0,0) = value;
261  }
262  else {
263  arma::vec3 value;
264  value.zeros();
265  for (int i=0; i<dh_->fe<3>()->n_dofs(); i++)
266  value += data_[dof_indices[i]]*fe_values3.shape_vector(i, 0);
267  for (int i=0; i<3; i++)
268  envelope(i,0) = value(i);
269  }
270  }
271  }
272 }
273 
274 
275 
276 template <int spacedim, class Value>
278 {
279  if (dof_indices != nullptr) delete[] dof_indices;
280 }
281 
282 
283 
284 #endif /* FIELD_FE_IMPL_HH_ */
Shape function values.
Definition: update_flags.hh:98
virtual ~FieldFE()
const Element * element() const
Definition: accessors.hh:76
virtual void value_list(const std::vector< Point > &point_list, const ElementAccessor< spacedim > &elm, std::vector< typename Value::return_type > &value_list)
Class FEValues calculates finite element data on the actual cells such as shape function values...
Node ** node
Definition: elements.h:90
void set_fe_data(const DOFHandlerMultiDim *dh, Mapping< 1, 3 > *map1, Mapping< 2, 3 > *map2, Mapping< 3, 3 > *map3, const Vec *data)
#define ASSERT_EQUAL(a, b)
Definition: global_defs.h:136
Provides the numbering of the finite element degrees of freedom on the computational mesh...
Definition: dofhandler.hh:251
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
Basic definitions of numerical quadrature rules.
Space< spacedim >::Point Point
const arma::vec::fixed< spacedim > shape_vector(const unsigned int function_no, const unsigned int point_no)
Return the value of the function_no-th shape function at the point_no-th quadrature point...
Definition: fe_values.cc:160
void set_point(const unsigned int i, const arma::vec::fixed< dim > &p)
Sets individual quadrature point coordinates.
Definition: quadrature.hh:152
FieldFE(unsigned int n_comp=0)
Implementation.
const double shape_value(const unsigned int function_no, const unsigned int point_no)
Return the value of the function_no-th shape function at the point_no-th quadrature point...
Definition: fe_values.cc:148
void reinit(ElementFullIter &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:244
Calculates finite element data on the actual cell.
Definition: fe_values.hh:370
Abstract class for description of finite elements.
unsigned int dim() const
Definition: accessors.hh:73
arma::vec3 & point()
Definition: nodes.hh:80
unsigned int idx() const
Definition: accessors.hh:103