Flow123d  DF_patch_fe_data_tables-398e8e7
patch_point_values.hh
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 patch_point_values.hh
15  * @brief Store finite element data on the actual patch
16  * such as shape function values, gradients, Jacobian
17  * of the mapping from the reference cell etc.
18  * @author David Flanderka
19  */
20 
21 #ifndef PATCH_POINT_VALUES_HH_
22 #define PATCH_POINT_VALUES_HH_
23 
24 #include <Eigen/Dense>
25 
26 #include "fem/eigen_tools.hh"
27 #include "fem/dh_cell_accessor.hh"
28 #include "fem/element_values.hh"
30 
31 
32 template<unsigned int spacedim> class PatchFEValues;
33 template<unsigned int spacedim> class ElOp;
34 template<unsigned int dim> class BulkValues;
35 template<unsigned int dim> class SideValues;
36 using Scalar = double;
39 
40 
41 
42 
43 /// Type for conciseness
44 using ReinitFunction = std::function<void(std::vector<ElOp<3>> &, TableDbl &, TableInt &)>;
45 
46 
47 namespace FeBulk {
48  /// enum of bulk operation
49  enum BulkOps
50  {
51  /// operations evaluated on elements
52  opElCoords, ///< coordinations of all nodes of element
53  opJac, ///< Jacobian of element
54  opInvJac, ///< inverse Jacobian
55  opJacDet, ///< determinant of Jacobian
56  /// operation executed expansion to quadrature point (value of element > values on quadrature points)
57  opExpansionCoords, ///< expands coordinates
58  opExpansionJac, ///< expands Jacobian
59  opExpansionInvJac, ///< expands inverse Jacobian
60  opExpansionJacDet, ///< expands Jacobian determinant
61  /// operations evaluated on quadrature points
62  opCoords, ///< coordinations of quadrature point
63  opWeights, ///< weight of quadrature point
64  opJxW ///< JxW value of quadrature point
65  };
66 }
67 
68 
69 namespace FeSide {
70  /// enum of side operation
71  enum SideOps
72  {
73  /// operations evaluated on elements
74  opElCoords, ///< coordinations of all nodes of element
75  opElJac, ///< Jacobian of element
76  opElInvJac, ///< inverse Jacobian of element
77  /// operations evaluated on sides
78  opSideCoords, ///< coordinations of all nodes of side
79  opSideJac, ///< Jacobian of element
80  opSideJacDet, ///< determinant of Jacobian of side
81  /// operation executed expansion to quadrature point (value of element / side > values on quadrature points)
82  opExpansionElCoords, ///< expands coordinates on element
83  opExpansionElJac, ///< expands Jacobian on element
84  opExpansionElInvJac, ///< expands inverse Jacobian on element
85  opExpansionSideCoords, ///< expands coordinates on side
86  opExpansionSideJac, ///< expands Jacobian on side
87  opExpansionSideJacDet, ///< expands Jacobian determinant on side
88  /// operations evaluated on quadrature points
89  opCoords, ///< coordinations of quadrature point
90  opWeights, ///< weight of quadrature point
91  opJxW, ///< JxW value of quadrature point
92  opNormalVec ///< normal vector of quadrature point
93  };
94 }
95 
96 
97 /**
98  * @brief Class for storing FE data of quadrature points.
99  *
100  * Store data of bulk or side quadrature points of one dimension.
101  */
102 template<unsigned int spacedim = 3>
104 {
105 public:
106  /**
107  * Constructor
108  *
109  * Set dimension
110  */
112  : dim_(dim), n_rows_(0), elements_map_(300, 0), points_map_(300, 0) {}
113 
114  /**
115  * Initialize object, set number of columns (quantities) in tables.
116  *
117  * Number of columns of int_vals_ table is passed by argument, number of columns
118  * of other tables is given by n_rows_ value.
119  */
120  void initialize(uint int_cols) {
121  this->reset();
122 
123  point_vals_.resize(n_rows_);
124  int_vals_.resize(int_cols);
125  }
126 
127  /// Reset number of rows (points and elements)
128  inline void reset() {
129  n_points_ = 0;
130  n_elems_ = 0;
131  }
132 
133  /// Getter of dim_
134  inline uint dim() const {
135  return dim_;
136  }
137 
138  /// Getter of n_rows_
139  inline uint n_rows() const {
140  return n_rows_;
141  }
142 
143  /// Getter of n_elems_
144  inline uint n_elems() const {
145  return n_elems_;
146  }
147 
148  /// Getter of n_points_
149  inline uint n_points() const {
150  return n_points_;
151  }
152 
153  /// Return quadrature
155  return quad_;
156  }
157 
158  /// Adds the number of rows equal to n_added, returns index of first of them
159  /// Temporary method allow acces to old structure PatchFEValues::DimPatchFEValues
160  inline uint add_rows(uint n_added) {
161  uint old_size = n_rows_;
162  n_rows_ += n_added;
163  return old_size;
164  }
165 
166  /// Resize data tables
170  }
171 
172  /// Register element, add to point_vals_ table
173  uint register_element(arma::mat coords, uint element_patch_idx) {
174  uint res_column = operations_[FeBulk::BulkOps::opElCoords].result_row();
175  for (uint i_col=0; i_col<coords.n_cols; ++i_col)
176  for (uint i_row=0; i_row<coords.n_rows; ++i_row) {
177  point_vals_(res_column)(n_elems_) = coords(i_row, i_col);
178  ++res_column;
179  }
180 
181  elements_map_[element_patch_idx] = n_elems_;
182  return n_elems_++;
183  }
184 
185  /// Register side, add to point_vals_ table
186  uint register_side(arma::mat elm_coords, arma::mat side_coords) {
187  uint res_column = operations_[FeSide::SideOps::opElCoords].result_row();
188  for (uint i_col=0; i_col<elm_coords.n_cols; ++i_col)
189  for (uint i_row=0; i_row<elm_coords.n_rows; ++i_row) {
190  point_vals_(res_column)(n_elems_) = elm_coords(i_row, i_col);
191  ++res_column;
192  }
193 
194  res_column = operations_[FeSide::SideOps::opSideCoords].result_row();
195  for (uint i_col=0; i_col<side_coords.n_cols; ++i_col)
196  for (uint i_row=0; i_row<side_coords.n_rows; ++i_row) {
197  point_vals_(res_column)(n_elems_) = side_coords(i_row, i_col);
198  ++res_column;
199  }
200 
201  return n_elems_++;
202  }
203 
204  /// Register bulk point, add to int_vals_ table
205  uint register_bulk_point(uint elem_table_row, uint value_patch_idx, uint elem_idx) {
206  int_vals_(0)(n_points_) = value_patch_idx;
207  int_vals_(1)(n_points_) = elem_table_row;
208  int_vals_(2)(n_points_) = elem_idx;
209 
210  points_map_[value_patch_idx] = n_points_;
211  return n_points_++;
212  }
213 
214  /// Register side point, add to int_vals_ table
215  uint register_side_point(uint elem_table_row, uint value_patch_idx, uint elem_idx, uint side_idx) {
216  int_vals_(0)(n_points_) = value_patch_idx;
217  int_vals_(1)(n_points_) = elem_table_row;
218  int_vals_(2)(n_points_) = elem_idx;
219  int_vals_(3)(n_points_) = side_idx;
220 
221  points_map_[value_patch_idx] = n_points_;
222  return n_points_++;
223  }
224 
225  /// Add accessor to operations_ vector
226  ElOp<spacedim> &make_new_op(std::initializer_list<uint> shape, ReinitFunction reinit_f, std::vector<uint> input_ops_vec) {
227  ElOp<spacedim> op_accessor(this->dim_, shape, this->n_rows_, reinit_f, input_ops_vec);
228  this->n_rows_ += op_accessor.n_comp();
229  operations_.push_back(op_accessor);
230  return operations_[operations_.size()-1];
231  }
232 
233  /// Add accessor to operations_ vector
234  ElOp<spacedim> &make_expansion(ElOp<spacedim> &el_op, std::initializer_list<uint> shape, ReinitFunction reinit_f) {
235  ElOp<spacedim> op_accessor(this->dim_, shape, el_op.result_row(), reinit_f);
236  // shape passed from el_op throws:
237  // C++ exception with description "std::bad_alloc" thrown in the test body.
238  operations_.push_back(op_accessor);
239  return operations_[operations_.size()-1];
240  }
241 
242  /// Add accessor to operations_ vector
243  ElOp<spacedim> &make_fe_op(std::initializer_list<uint> shape, ReinitFunction reinit_f, std::vector<uint> input_ops_vec, uint n_dofs) {
244  ElOp<spacedim> op_accessor(this->dim_, shape, this->n_rows_, reinit_f, input_ops_vec, n_dofs);
245  this->n_rows_ += op_accessor.n_comp() * n_dofs;
246  operations_.push_back(op_accessor);
247  return operations_[operations_.size()-1];
248  }
249 
250 
251  /// Reinit data.
252  void reinit_patch() {
253  if (n_elems_ == 0) return; // skip if tables are empty
254  // precompute data on point_vals_ table
255  for (uint i=0; i<operations_.size(); ++i)
256  operations_[i].reinit_function(operations_, point_vals_, int_vals_);
257 
258 // // copy data from el_vals_ to point_vals_
259 // std::vector<uint> copied; // list of columns that will be copied
260 // for (uint i_op=0; i_op<operations_.size(); ++i_op)
261 // if (operations_[i_op].copy_vals()) {
262 // uint res_col = operations_[i_op].result_row();
263 // for (uint i_col=res_col; i_col<res_col+operations_[i_op].n_comp(); ++i_col)
264 // copied.push_back(i_col);
265 // }
266 // for (uint i_pt=0; i_pt<n_points_; ++i_pt) {
267 // uint el_table_idx = int_vals_(1)(i_pt);
268 // for (uint i_q=0; i_q<copied.size(); ++i_q)
269 // point_vals_(copied[i_q])(i_pt) = el_vals_(copied[i_q])(el_table_idx);
270 // }
271 //
272 // // precompute data on point_vals_ table
273 // for (uint i=0; i<operations_.size(); ++i)
274 // operations_[i].reinit_points(operations_, point_vals_);
275  }
276 
277  inline Scalar scalar_val(uint result_row, uint point_idx) const {
278  return point_vals_(result_row)(elements_map_[point_idx]);
279  }
280 
281  inline Vector vector_val(uint result_row, uint point_idx) const {
282  Vector val;
283  for (uint i=0; i<3; ++i)
284  val(i) = point_vals_(result_row+i)(elements_map_[point_idx]);
285  return val;
286  }
287 
288  inline Tensor tensor_val(uint result_row, uint point_idx) const {
289  Tensor val;
290  for (uint i=0; i<3; ++i)
291  for (uint j=0; j<3; ++j)
292  val(i,j) = point_vals_(result_row+3*i+j)(elements_map_[point_idx]);
293  return val;
294  }
295 
296  /// Temporary development method
297  void print(bool points, bool ints) const {
298  std::cout << "** Dimension: " << dim_ << std::endl;
299  if (points) {
300  std::cout << "Point vals: " << point_vals_.rows() << " - " << point_vals_.cols() << std::endl;
301  for (uint i_row=0; i_row<n_points_; ++i_row) {
302  for (uint i_col=0; i_col<n_rows_; ++i_col)
303  std::cout << point_vals_(i_col)(i_row) << " ";
304  std::cout << std::endl;
305  }
306  std::cout << std::endl;
307  }
308  if (ints) {
309  std::cout << "Int vals: " << int_vals_.rows() << " - " << int_vals_.cols() << std::endl;
310  for (uint i_row=0; i_row<n_points_; ++i_row) {
311  for (uint i_col=0; i_col<3; ++i_col)
312  std::cout << int_vals_(i_col)(i_row) << " ";
313  std::cout << std::endl;
314  }
315  std::cout << std::endl;
316  }
317  std::cout << "*****************" << std::endl;
318  }
319 
320 protected:
321  /**
322  * Store data of bulk or side quadrature points of one dimension
323  *
324  * Number of columns is given by n_rows_, number of used rows by n_points_.
325  */
327  /**
328  * Hold integer values of quadrature points of previous table.
329  *
330  * Table contains following columns:
331  * 0: Index of quadrature point on patch
332  * 1: Row of element/side in point_vals_ table in registration step (before expansion)
333  * 2: Element idx in Mesh
334  * 3: Index of side in element (column is allocated only for side point table)
335  * Number of used rows is given by n_points_.
336  */
338 
339  /// Vector of all defined operations
341 
342  uint dim_; ///< Dimension
343  uint n_rows_; ///< Number of columns of \p point_vals table
344  uint n_points_; ///< Number of points in patch
345  uint n_elems_; ///< Number of elements in patch
346  Quadrature *quad_; ///< Quadrature of given dimension and order passed in constructor.
347 
348  std::vector<uint> elements_map_; ///< Map of element patch indices to el_vals_ table
349  std::vector<uint> points_map_; ///< Map of point patch indices to point_vals_ and int_vals_ tables
350 
351  friend class PatchFEValues<spacedim>;
352  friend class ElOp<spacedim>;
353  template<unsigned int dim>
354  friend class BulkValues;
355  template<unsigned int dim>
356  friend class SideValues;
357  template<unsigned int dim>
358  friend class JoinValues;
359 };
360 
361 /**
362  * @brief Class represents FE operations.
363  */
364 template<unsigned int spacedim = 3>
365 class ElOp {
366 public:
367  /// Constructor
368  ElOp(uint dim, std::initializer_list<uint> shape, uint result_row, ReinitFunction reinit_f, std::vector<uint> input_ops = {}, uint n_dofs = 1)
370  {}
371 
372  /// Number of components computed from shape_ vector
373  inline uint n_comp() const {
374  if (shape_.size() == 1) return shape_[0];
375  else return shape_[0] * shape_[1];
376  }
377 
378  /// Getter of dimension
379  inline uint dim() const {
380  return dim_;
381  }
382 
383  /// Getter of result_row_
384  inline uint result_row() const {
385  return result_row_;
386  }
387 
388  /// Getter of input_ops_
389  inline const std::vector<uint> &input_ops() const {
390  return input_ops_;
391  }
392 
393  /// Getter of shape_
394  inline const std::vector<uint> &shape() const {
395  return shape_;
396  }
397 
398  /// Call reinit function on element table if function is defined
399  inline void reinit_function(std::vector<ElOp<spacedim>> &operations, TableDbl &data_table, TableInt &int_table) {
400  reinit_func(operations, data_table, int_table);
401  }
402 
403  /// Return map referenced Eigen::Matrix of given dimension
404  template<unsigned int dim1, unsigned int dim2>
405  Eigen::Map<Eigen::Matrix<ArrayDbl, dim1, dim2>> value(TableDbl &op_results, uint i_dof = 0) const {
406  return Eigen::Map<Eigen::Matrix<ArrayDbl, dim1, dim2>>(op_results.data() + result_row_ + i_dof * n_comp(), dim1, dim2);
407  }
408 
409 
410 protected:
411  uint dim_; ///< Dimension
412  std::vector<uint> shape_; ///< Shape of stored data (size of vector or number of rows and cols of matrix)
413  uint result_row_; ///< First row to scalar, vector or matrix result
414  std::vector<uint> input_ops_; ///< Indices of operations in PatchPointValues::operations_ vector on which ElOp is depended
415  uint n_dofs_; ///< Number of DOFs of FE operations (or 1 in case of element operations)
416 
417  ReinitFunction reinit_func; ///< Pointer to patch reinit function of element data table specialized by operation
418 };
419 
420 
421 /// Defines common functionality of reinit operations.
423  // empty base operation
424  static inline void op_base(FMT_UNUSED std::vector<ElOp<3>> &operations, FMT_UNUSED TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
425  // empty
426  }
427 
428  // expansion operation
429  static inline void expand_data(ElOp<3> &op, TableDbl &op_results, TableInt &el_table) {
430  uint row_begin = op.result_row();
431  uint row_end = row_begin + op.n_comp();
432  uint size = op_results(row_begin).rows();
433  for (int i_pt=size-1; i_pt>=0; --i_pt) {
434  uint el_table_idx = el_table(1)(i_pt);
435  for (uint i_q=row_begin; i_q<row_end; ++i_q)
436  op_results(i_q)(i_pt) = op_results(i_q)(el_table_idx);
437  }
438  }
439 };
440 
441 /// Defines reinit operations on bulk points.
442 struct bulk_reinit {
443  // element operations
444  template<unsigned int dim>
445  static inline void elop_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
446  // result matrix(spacedim, dim), input matrix(spacedim, dim+1)
447  auto &op = operations[FeBulk::BulkOps::opJac];
448  auto jac_value = op.value<3, dim>(op_results);
449  auto coords_value = operations[ op.input_ops()[0] ].value<3, dim+1>(op_results);
450  jac_value = eigen_tools::jacobian<3,dim>(coords_value);
451  }
452  template<unsigned int dim>
453  static inline void elop_inv_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
454  // result matrix(spacedim, dim), input matrix(spacedim, dim+1)
455  auto &op = operations[FeBulk::BulkOps::opInvJac];
456  auto inv_jac_value = op.value<dim, 3>(op_results);
457  auto jac_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
458  inv_jac_value = eigen_tools::inverse<3, dim>(jac_value);
459  }
460  template<unsigned int dim>
461  static inline void elop_jac_det(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
462  // result double, input matrix(spacedim, dim)
463  auto &op = operations[FeBulk::BulkOps::opJacDet];
464  auto &jac_det_value = op_results(op.result_row());
465  auto jac_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
466  jac_det_value = eigen_tools::determinant<Eigen::Matrix<ArrayDbl, 3, dim>>(jac_value);
467  }
468 
469  // expansion operations
470  static inline void expd_coords(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
471  auto &op = operations[FeBulk::BulkOps::opExpansionCoords];
472  common_reinit::expand_data(op, op_results, el_table);
473  }
474  static inline void expd_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
475  auto &op = operations[FeBulk::BulkOps::opExpansionJac];
476  common_reinit::expand_data(op, op_results, el_table);
477  }
478  static inline void expd_inv_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
479  auto &op = operations[FeBulk::BulkOps::opExpansionInvJac];
480  common_reinit::expand_data(op, op_results, el_table);
481  }
482  static inline void expd_jac_det(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
483  auto &op = operations[FeBulk::BulkOps::opExpansionJacDet];
484  common_reinit::expand_data(op, op_results, el_table);
485  }
486 
487  // point operations
488  static inline void ptop_coords(FMT_UNUSED std::vector<ElOp<3>> &operations, FMT_UNUSED TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
489  // Implement
490  }
491  static inline void ptop_weights(std::vector<ElOp<3>> &operations, TableDbl &op_results, std::vector<double> point_weights) {
492  auto &op = operations[FeBulk::BulkOps::opWeights];
493  ArrayDbl &result_row = op_results( op.result_row() );
494  auto n_points = point_weights.size();
495  for (uint i=0; i<result_row.rows(); ++i)
496  result_row(i) = point_weights[i%n_points];
497  }
498  static inline void ptop_JxW(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
499  auto &op = operations[FeBulk::BulkOps::opJxW];
500  ArrayDbl &weights_row = op_results( operations[op.input_ops()[0]].result_row() );
501  ArrayDbl &jac_det_row = op_results( operations[op.input_ops()[1]].result_row() );
502  ArrayDbl &result_row = op_results( op.result_row() );
503  result_row = jac_det_row * weights_row;
504  }
505  static inline void ptop_scalar_shape(std::vector<ElOp<3>> &operations, TableDbl &op_results,
506  std::vector< std::vector<double> > shape_values, uint scalar_shape_op_idx) {
507  auto &op = operations[scalar_shape_op_idx];
508  uint n_points = shape_values.size();
509 
510  for (uint i_row=0; i_row<shape_values[0].size(); ++i_row) {
511  ArrayDbl &result_row = op_results( op.result_row()+i_row );
512  for (uint i_pt=0; i_pt<result_row.rows(); ++i_pt)
513  result_row(i_pt) = shape_values[i_pt % n_points][i_row];
514  }
515  }
516  template<unsigned int dim>
517  static inline void ptop_scalar_shape_grads(std::vector<ElOp<3>> &operations, TableDbl &op_results,
518  std::vector< std::vector<arma::mat> > ref_shape_grads, uint scalar_shape_grads_op_idx) {
519  auto &op = operations[scalar_shape_grads_op_idx];
520  uint n_points = ref_shape_grads.size();
521  uint n_dofs = ref_shape_grads[0].size();
522 
523  Eigen::Vector<ArrayDbl, Eigen::Dynamic> ref_shape_grads_expd;
524  ref_shape_grads_expd.resize(dim*n_dofs);
525  for (uint i=0; i<ref_shape_grads_expd.rows(); ++i)
526  ref_shape_grads_expd(i).resize(op_results(0).rows());
527 
528  for (uint i_dof=0; i_dof<n_dofs; ++i_dof)
529  for (uint i_c=0; i_c<dim; ++i_c) {
530  ArrayDbl &shape_grad_row = ref_shape_grads_expd(i_dof*dim+i_c);
531  for (uint i_pt=0; i_pt<shape_grad_row.rows(); ++i_pt)
532  shape_grad_row(i_pt) = ref_shape_grads[i_pt % n_points][i_dof][i_c];
533  }
534 
535  auto inv_jac_value = operations[ op.input_ops()[0] ].value<dim, 3>(op_results);
536  for (uint i_dof=0; i_dof<n_dofs; ++i_dof) {
537  auto shape_grad_value = op.value<3, 1>(op_results, i_dof);
538  Eigen::Map<Eigen::Matrix<ArrayDbl, dim, 1>> ref_shape_grads_dof_value(ref_shape_grads_expd.data() + dim*i_dof, dim, 1);
539  shape_grad_value = inv_jac_value.transpose() * ref_shape_grads_dof_value;
540  }
541  }
542 };
543 
544 
545 
546 /// Defines reinit operations on side points.
547 struct side_reinit {
548  // element operations
549  template<unsigned int dim>
550  static inline void elop_el_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
551  // result matrix(spacedim, dim), input matrix(spacedim, dim+1)
552  auto &op = operations[FeSide::SideOps::opElJac];
553  auto jac_value = op.value<3, dim>(op_results);
554  auto coords_value = operations[ op.input_ops()[0] ].value<3, dim+1>(op_results);
555  jac_value = eigen_tools::jacobian<3,dim>(coords_value);
556  }
557  template<unsigned int dim>
558  static inline void elop_el_inv_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
559  auto &op = operations[FeSide::SideOps::opElInvJac];
560  auto inv_jac_value = op.value<dim, 3>(op_results);
561  auto jac_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
562  inv_jac_value = eigen_tools::inverse<3,dim>(jac_value);
563  }
564  template<unsigned int dim>
565  static inline void elop_sd_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
566  // result matrix(spacedim, dim), input matrix(spacedim, dim+1)
567  auto &op = operations[FeSide::SideOps::opSideJac];
568  auto jac_value = op.value<3, dim-1>(op_results);
569  auto coords_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
570  jac_value = eigen_tools::jacobian<3, dim-1>(coords_value);
571  }
572  template<unsigned int dim>
573  static inline void elop_sd_jac_det(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
574  // result double, input matrix(spacedim, dim)
575  auto &op = operations[FeSide::SideOps::opSideJacDet];
576  ArrayDbl &det_value = op_results( op.result_row() );
577  auto jac_value = operations[ op.input_ops()[0] ].value<3, dim-1>(op_results);
578  det_value = eigen_tools::determinant<Eigen::Matrix<ArrayDbl, 3, dim-1>>(jac_value);
579  }
580 
581  // expansion operations
582  static inline void expd_el_coords(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
583  auto &op = operations[FeSide::SideOps::opExpansionElCoords];
584  common_reinit::expand_data(op, op_results, el_table);
585  }
586  static inline void expd_el_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
587  auto &op = operations[FeSide::SideOps::opExpansionElJac];
588  common_reinit::expand_data(op, op_results, el_table);
589  }
590  static inline void expd_el_inv_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
591  auto &op = operations[FeSide::SideOps::opExpansionElInvJac];
592  common_reinit::expand_data(op, op_results, el_table);
593  }
594 
595  static inline void expd_sd_coords(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
596  auto &op = operations[FeSide::SideOps::opExpansionSideCoords];
597  common_reinit::expand_data(op, op_results, el_table);
598  }
599  static inline void expd_sd_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
600  auto &op = operations[FeSide::SideOps::opExpansionSideJac];
601  common_reinit::expand_data(op, op_results, el_table);
602  }
603  static inline void expd_sd_jac_det(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
604  auto &op = operations[FeSide::SideOps::opExpansionSideJacDet];
605  common_reinit::expand_data(op, op_results, el_table);
606  }
607 
608  // Point operations
609  static inline void ptop_coords(FMT_UNUSED std::vector<ElOp<3>> &operations, FMT_UNUSED TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
610  // Implement
611  }
612  static inline void ptop_weights(std::vector<ElOp<3>> &operations, TableDbl &op_results, std::vector<double> point_weights) {
613  auto &op = operations[FeSide::SideOps::opWeights];
614  ArrayDbl &result_row = op_results( op.result_row() );
615  auto n_points = point_weights.size();
616  for (uint i=0; i<result_row.rows(); ++i)
617  result_row(i) = point_weights[i%n_points];
618  }
619  static inline void ptop_JxW(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
620  auto &op = operations[FeSide::SideOps::opJxW];
621  ArrayDbl &weights_row = op_results( operations[op.input_ops()[0]].result_row() );
622  ArrayDbl &jac_det_row = op_results( operations[op.input_ops()[1]].result_row() );
623  ArrayDbl &result_row = op_results( op.result_row() );
624  result_row = jac_det_row * weights_row;
625  }
626  template<unsigned int dim>
627  static inline void ptop_normal_vec(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
628  auto &op = operations[FeSide::SideOps::opNormalVec];
629  auto normal_value = op.value<3, 1>(op_results);
630  auto inv_jac_mat_value = operations[ op.input_ops()[0] ].value<dim, 3>(op_results);
631  normal_value = inv_jac_mat_value.transpose() * RefElement<dim>::normal_vector_array( el_table(3) );
632  }
633  static inline void ptop_scalar_shape(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table,
634  std::vector< std::vector< std::vector<double> > > shape_values, uint scalar_shape_op_idx) {
635  auto &op = operations[scalar_shape_op_idx];
636  uint n_points = shape_values[0].size();
637 
638  for (uint i_row=0; i_row<shape_values[0][0].size(); ++i_row) {
639  ArrayDbl &result_row = op_results( op.result_row()+i_row );
640  for (uint i_pt=0; i_pt<result_row.rows(); ++i_pt)
641  result_row(i_pt) = shape_values[el_table(3)(i_pt)][i_pt % n_points][i_row];
642  }
643  }
644  template<unsigned int dim>
645  static inline void ptop_scalar_shape_grads(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table,
646  std::vector< std::vector< std::vector<arma::mat> > > ref_shape_grads, uint scalar_shape_grads_op_idx) {
647  auto &op = operations[scalar_shape_grads_op_idx];
648  uint n_points = ref_shape_grads[0].size();
649  uint n_dofs = ref_shape_grads[0][0].size();
650 
651  Eigen::Vector<ArrayDbl, Eigen::Dynamic> ref_shape_grads_expd;
652  ref_shape_grads_expd.resize(dim*n_dofs);
653  for (uint i=0; i<ref_shape_grads_expd.rows(); ++i)
654  ref_shape_grads_expd(i).resize(op_results(0).rows());
655 
656  for (uint i_dof=0; i_dof<n_dofs; ++i_dof)
657  for (uint i_c=0; i_c<dim; ++i_c) {
658  ArrayDbl &shape_grad_row = ref_shape_grads_expd(i_dof*dim+i_c);
659  for (uint i_pt=0; i_pt<shape_grad_row.rows(); ++i_pt)
660  shape_grad_row(i_pt) = ref_shape_grads[el_table(3)(i_pt)][i_pt % n_points][i_dof][i_c];
661  }
662 
663  auto inv_jac_value = operations[ op.input_ops()[0] ].value<dim, 3>(op_results);
664  for (uint i_dof=0; i_dof<n_dofs; ++i_dof) {
665  auto shape_grad_value = op.value<3, 1>(op_results, i_dof);
666  Eigen::Map<Eigen::Matrix<ArrayDbl, dim, 1>> ref_shape_grads_dof_value(ref_shape_grads_expd.data() + dim*i_dof, dim, 1);
667  shape_grad_value = inv_jac_value.transpose() * ref_shape_grads_dof_value;
668  }
669  }
670 };
671 
672 
673 
674 namespace FeBulk {
675 
676  /// Bulk data specialization, order of item in operations_ vector corresponds to the BulkOps enum
677  template<unsigned int spacedim = 3>
678  class PatchPointValues : public ::PatchPointValues<spacedim> {
679  public:
680  /// Constructor
682  : ::PatchPointValues<spacedim>(dim) {
683  this->quad_ = new QGauss(dim, 2*quad_order);
684  switch (dim) {
685  case 1:
686  init<1>();
687  break;
688  case 2:
689  init<2>();
690  break;
691  case 3:
692  init<3>();
693  break;
694  }
695  }
696 
697  private:
698  /// Initialize operations vector
699  template<unsigned int dim>
700  void init() {
701  // First step: adds element values operations
702  auto &el_coords = this->make_new_op( {spacedim, this->dim_+1}, &common_reinit::op_base, {} );
703 
704  auto &el_jac = this->make_new_op( {spacedim, this->dim_}, &bulk_reinit::elop_jac<dim>, {BulkOps::opElCoords} );
705 
706  auto &el_inv_jac = this->make_new_op( {this->dim_, spacedim}, &bulk_reinit::elop_inv_jac<dim>, {BulkOps::opJac} );
707 
708  auto &el_jac_det = this->make_new_op( {1}, &bulk_reinit::elop_jac_det<dim>, {BulkOps::opJac} );
709 
710  // Second step: adds expand operations (element values to point values)
711  this->make_expansion( el_coords, {spacedim, this->dim_+1}, &bulk_reinit::expd_coords );
712 
713  this->make_expansion( el_jac, {spacedim, this->dim_}, &bulk_reinit::expd_jac );
714 
715  this->make_expansion( el_inv_jac, {this->dim_, spacedim}, &bulk_reinit::expd_inv_jac );
716 
717  this->make_expansion( el_jac_det, {1}, &bulk_reinit::expd_jac_det );
718 
719  // Third step: adds point values operations
720  /*auto &pt_coords =*/ this->make_new_op( {spacedim}, &bulk_reinit::ptop_coords, {} );
721 
722  // use lambda reinit function
723  std::vector<double> point_weights = this->quad_->get_weights();
724  auto lambda_weights = [point_weights](std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
725  bulk_reinit::ptop_weights(operations, op_results, point_weights);
726  };
727  /*auto &weights =*/ this->make_new_op( {1}, lambda_weights, {} );
728 
729  /*auto &JxW =*/ this->make_new_op( {1}, &bulk_reinit::ptop_JxW, {BulkOps::opWeights, BulkOps::opJacDet} );
730  }
731  };
732 
733 } // closing namespace FeBulk
734 
735 
736 
737 namespace FeSide {
738 
739 /// Bulk Side specialization, order of item in operations_ vector corresponds to the SideOps enum
740  template<unsigned int spacedim = 3>
741  class PatchPointValues : public ::PatchPointValues<spacedim> {
742  public:
743  /// Constructor
745  : ::PatchPointValues<spacedim>(dim) {
746  this->quad_ = new QGauss(dim-1, 2*quad_order);
747  switch (dim) {
748  case 1:
749  init<1>();
750  break;
751  case 2:
752  init<2>();
753  break;
754  case 3:
755  init<3>();
756  break;
757  }
758  }
759 
760  private:
761  /// Initialize operations vector
762  template<unsigned int dim>
763  void init() {
764  // First step: adds element values operations
765  auto &el_coords = this->make_new_op( {spacedim, this->dim_+1}, &common_reinit::op_base, {} );
766 
767  auto &el_jac = this->make_new_op( {spacedim, this->dim_}, &side_reinit::elop_el_jac<dim>, {SideOps::opElCoords} );
768 
769  auto &el_inv_jac = this->make_new_op( {this->dim_, spacedim}, &side_reinit::elop_el_inv_jac<dim>, {SideOps::opElJac} );
770 
771  auto &sd_coords = this->make_new_op( {spacedim, this->dim_}, &common_reinit::op_base, {} );
772 
773  auto &sd_jac = this->make_new_op( {spacedim, this->dim_-1}, &side_reinit::elop_sd_jac<dim>, {SideOps::opSideCoords} );
774 
775  auto &sd_jac_det = this->make_new_op( {1}, &side_reinit::elop_sd_jac_det<dim>, {SideOps::opSideJac} );
776 
777  // Second step: adds expand operations (element values to point values)
778  this->make_expansion( el_coords, {spacedim, this->dim_+1}, &side_reinit::expd_el_coords );
779 
780  this->make_expansion( el_jac, {spacedim, this->dim_}, &side_reinit::expd_el_jac );
781 
782  this->make_expansion( el_inv_jac, {this->dim_, spacedim}, &side_reinit::expd_el_inv_jac );
783 
784  this->make_expansion( sd_coords, {spacedim, this->dim_}, &side_reinit::expd_sd_coords );
785 
786  this->make_expansion( sd_jac, {spacedim, this->dim_-1}, &side_reinit::expd_sd_jac );
787 
788  this->make_expansion( sd_jac_det, {1}, &side_reinit::expd_sd_jac_det );
789 
790  // Third step: adds point values operations
791  /*auto &coords =*/ this->make_new_op( {spacedim}, &side_reinit::ptop_coords, {} );
792 
793  // use lambda reinit function
794  std::vector<double> point_weights = this->quad_->get_weights();
795  auto lambda_weights = [point_weights](std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
796  side_reinit::ptop_weights(operations, op_results, point_weights);
797  };
798  /*auto &weights =*/ this->make_new_op( {1}, lambda_weights, {} );
799 
801 
802  /*auto &normal_vec =*/ this->make_new_op( {spacedim}, &side_reinit::ptop_normal_vec<dim>, {SideOps::opElInvJac} );
803  }
804  };
805 
806 } // closing namespace FeSide
807 
808 
809 
810 #endif /* PATCH_POINT_VALUES_HH_ */
Class represents FE operations.
std::vector< uint > input_ops_
Indices of operations in PatchPointValues::operations_ vector on which ElOp is depended.
uint dim_
Dimension.
std::vector< uint > shape_
Shape of stored data (size of vector or number of rows and cols of matrix)
uint result_row() const
Getter of result_row_.
const std::vector< uint > & shape() const
Getter of shape_.
const std::vector< uint > & input_ops() const
Getter of input_ops_.
uint result_row_
First row to scalar, vector or matrix result.
ReinitFunction reinit_func
Pointer to patch reinit function of element data table specialized by operation.
ElOp(uint dim, std::initializer_list< uint > shape, uint result_row, ReinitFunction reinit_f, std::vector< uint > input_ops={}, uint n_dofs=1)
Constructor.
uint n_dofs_
Number of DOFs of FE operations (or 1 in case of element operations)
uint dim() const
Getter of dimension.
Eigen::Map< Eigen::Matrix< ArrayDbl, dim1, dim2 > > value(TableDbl &op_results, uint i_dof=0) const
Return map referenced Eigen::Matrix of given dimension.
uint n_comp() const
Number of components computed from shape_ vector.
void reinit_function(std::vector< ElOp< spacedim >> &operations, TableDbl &data_table, TableInt &int_table)
Call reinit function on element table if function is defined.
Bulk data specialization, order of item in operations_ vector corresponds to the BulkOps enum.
PatchPointValues(uint dim, uint quad_order)
Constructor.
void init()
Initialize operations vector.
Bulk Side specialization, order of item in operations_ vector corresponds to the SideOps enum.
PatchPointValues(uint dim, uint quad_order)
Constructor.
void init()
Initialize operations vector.
Class for storing FE data of quadrature points.
ElOp< spacedim > & make_expansion(ElOp< spacedim > &el_op, std::initializer_list< uint > shape, ReinitFunction reinit_f)
Add accessor to operations_ vector.
uint add_rows(uint n_added)
Temporary method allow acces to old structure PatchFEValues::DimPatchFEValues.
void initialize(uint int_cols)
void resize_tables(uint n_points)
Resize data tables.
uint n_rows() const
Getter of n_rows_.
std::vector< uint > points_map_
Map of point patch indices to point_vals_ and int_vals_ tables.
uint register_side_point(uint elem_table_row, uint value_patch_idx, uint elem_idx, uint side_idx)
Register side point, add to int_vals_ table.
uint n_rows_
Number of columns of point_vals table.
uint register_element(arma::mat coords, uint element_patch_idx)
Register element, add to point_vals_ table.
Tensor tensor_val(uint result_row, uint point_idx) const
uint register_side(arma::mat elm_coords, arma::mat side_coords)
Register side, add to point_vals_ table.
std::vector< ElOp< spacedim > > operations_
Vector of all defined operations.
Quadrature * get_quadrature() const
Return quadrature.
uint n_elems() const
Getter of n_elems_.
void reset()
Reset number of rows (points and elements)
uint n_points() const
Getter of n_points_.
std::vector< uint > elements_map_
Map of element patch indices to el_vals_ table.
uint dim() const
Getter of dim_.
void print(bool points, bool ints) const
Temporary development method.
ElOp< spacedim > & make_fe_op(std::initializer_list< uint > shape, ReinitFunction reinit_f, std::vector< uint > input_ops_vec, uint n_dofs)
Add accessor to operations_ vector.
Vector vector_val(uint result_row, uint point_idx) const
void reinit_patch()
Reinit data.
uint n_points_
Number of points in patch.
Quadrature * quad_
Quadrature of given dimension and order passed in constructor.
Scalar scalar_val(uint result_row, uint point_idx) const
uint register_bulk_point(uint elem_table_row, uint value_patch_idx, uint elem_idx)
Register bulk point, add to int_vals_ table.
uint n_elems_
Number of elements in patch.
ElOp< spacedim > & make_new_op(std::initializer_list< uint > shape, ReinitFunction reinit_f, std::vector< uint > input_ops_vec)
Add accessor to operations_ vector.
Symmetric Gauss-Legendre quadrature formulae on simplices.
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
const std::vector< double > & get_weights() const
Return a reference to the whole array of weights.
Definition: quadrature.hh:111
static Eigen::Vector< ArrayDbl, dim > normal_vector_array(Eigen::Array< uint, Eigen::Dynamic, 1 > loc_side_idx_array)
Definition: ref_element.cc:279
Store finite element data on the actual patch such as shape function values, gradients,...
Eigen::Vector< ArrayDbl, Eigen::Dynamic > TableDbl
Definition: eigen_tools.hh:33
Eigen::Array< double, Eigen::Dynamic, 1 > ArrayDbl
Definitions of Eigen structures.
Definition: eigen_tools.hh:31
Eigen::Vector< ArrayInt, Eigen::Dynamic > TableInt
Definition: eigen_tools.hh:34
Class ElementValues calculates data related to transformation of reference cell to actual cell (Jacob...
unsigned int uint
ArmaMat< double, N, M > mat
Definition: armor.hh:936
BulkOps
enum of bulk operation
@ opExpansionJac
expands Jacobian
@ opCoords
operations evaluated on quadrature points
@ opInvJac
inverse Jacobian
@ opWeights
weight of quadrature point
@ opJac
Jacobian of element.
@ opExpansionInvJac
expands inverse Jacobian
@ opExpansionCoords
operation executed expansion to quadrature point (value of element > values on quadrature points)
@ opElCoords
operations evaluated on elements
@ opJxW
JxW value of quadrature point.
SideOps
enum of side operation
@ opExpansionElCoords
operation executed expansion to quadrature point (value of element / side > values on quadrature poin...
@ opExpansionElJac
expands Jacobian on element
@ opNormalVec
normal vector of quadrature point
@ opSideCoords
operations evaluated on sides
@ opExpansionSideCoords
expands coordinates on side
@ opJxW
JxW value of quadrature point.
@ opCoords
operations evaluated on quadrature points
@ opExpansionElInvJac
expands inverse Jacobian on element
@ opExpansionSideJac
expands Jacobian on side
@ opElCoords
operations evaluated on elements
@ opElJac
Jacobian of element.
@ opWeights
weight of quadrature point
@ opSideJac
Jacobian of element.
Eigen::Matrix< ArrayDbl, spacedim, dim > jacobian(const Eigen::Matrix< ArrayDbl, spacedim, dim+1 > &coords)
Definition: eigen_tools.hh:228
void resize_table(typename Eigen::Vector< ET, Eigen::Dynamic > &table, uint size)
Resize vector of Eigen::Array to given size.
Definition: eigen_tools.hh:41
ArrayDbl determinant(const T &M)
Calculates determinant of a rectangular matrix.
std::function< void(std::vector< ElOp< 3 > > &, TableDbl &, TableInt &)> ReinitFunction
Type for conciseness.
double Scalar
#define FMT_UNUSED
Definition: posix.h:75
Definitions of particular quadrature rules on simplices.
Defines reinit operations on bulk points.
static void expd_jac_det(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
static void ptop_scalar_shape_grads(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, std::vector< std::vector< arma::mat > > ref_shape_grads, uint scalar_shape_grads_op_idx)
static void ptop_coords(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void ptop_scalar_shape(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, std::vector< std::vector< double > > shape_values, uint scalar_shape_op_idx)
static void elop_jac_det(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void elop_inv_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void elop_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void ptop_weights(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, std::vector< double > point_weights)
static void expd_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
static void ptop_JxW(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void expd_coords(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
static void expd_inv_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
Defines common functionality of reinit operations.
static void expand_data(ElOp< 3 > &op, TableDbl &op_results, TableInt &el_table)
static void op_base(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED TableDbl &op_results, FMT_UNUSED TableInt &el_table)
Defines reinit operations on side points.
static void ptop_normal_vec(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
static void elop_el_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void ptop_weights(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, std::vector< double > point_weights)
static void expd_sd_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
static void elop_sd_jac_det(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void expd_sd_coords(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
static void ptop_scalar_shape(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table, std::vector< std::vector< std::vector< double > > > shape_values, uint scalar_shape_op_idx)
static void expd_el_inv_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
static void expd_sd_jac_det(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
static void ptop_JxW(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void expd_el_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
static void expd_el_coords(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
static void elop_el_inv_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void ptop_coords(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void elop_sd_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void ptop_scalar_shape_grads(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table, std::vector< std::vector< std::vector< arma::mat > > > ref_shape_grads, uint scalar_shape_grads_op_idx)