Flow123d  DF_patch_fe_data_tables-65e7202
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  /**
49  * Enumeration of element bulk operations
50  *
51  * Operations are stored in fix order. Order in enum is equal to order
52  * in PatchPointVale::operations_ vector. FE operations are added dynamically
53  * by request of user.
54  */
55  enum BulkOps
56  {
57  /// operations evaluated on elements
58  opElCoords, ///< coordinations of all nodes of element
59  opJac, ///< Jacobian of element
60  opInvJac, ///< inverse Jacobian
61  opJacDet, ///< determinant of Jacobian
62  /// operation executed expansion to quadrature point (value of element > values on quadrature points)
63  opExpansionCoords, ///< expands coordinates
64  opExpansionJac, ///< expands Jacobian
65  opExpansionInvJac, ///< expands inverse Jacobian
66  opExpansionJacDet, ///< expands Jacobian determinant
67  /// operations evaluated on quadrature points
68  opCoords, ///< coordinations of quadrature point
69  opWeights, ///< weight of quadrature point
70  opJxW ///< JxW value of quadrature point
71  };
72 }
73 
74 
75 namespace FeSide {
76  /**
77  * Enumeration of element side operations
78  *
79  * Operations are stored in fix order. Order in enum is equal to order
80  * in PatchPointVale::operations_ vector. FE operations are added dynamically
81  * by request of user.
82  */
83  enum SideOps
84  {
85  /// operations evaluated on elements
86  opElCoords, ///< coordinations of all nodes of element
87  opElJac, ///< Jacobian of element
88  opElInvJac, ///< inverse Jacobian of element
89  /// operations evaluated on sides
90  opSideCoords, ///< coordinations of all nodes of side
91  opSideJac, ///< Jacobian of element
92  opSideJacDet, ///< determinant of Jacobian of side
93  /// operation executed expansion to quadrature point (value of element / side > values on quadrature points)
94  opExpansionElCoords, ///< expands coordinates on element
95  opExpansionElJac, ///< expands Jacobian on element
96  opExpansionElInvJac, ///< expands inverse Jacobian on element
97  opExpansionSideCoords, ///< expands coordinates on side
98  opExpansionSideJac, ///< expands Jacobian on side
99  opExpansionSideJacDet, ///< expands Jacobian determinant on side
100  /// operations evaluated on quadrature points
101  opCoords, ///< coordinations of quadrature point
102  opWeights, ///< weight of quadrature point
103  opJxW, ///< JxW value of quadrature point
104  opNormalVec ///< normal vector of quadrature point
105  };
106 }
107 
108 
109 /**
110  * v Class for storing FE data of quadrature points on one patch.
111  *
112  * Store data of bulk or side quadrature points of one dimension.
113  */
114 template<unsigned int spacedim = 3>
116 {
117 public:
118  /**
119  * Constructor
120  *
121  * @param dim Set dimension
122  */
124  : dim_(dim), n_rows_(0), elements_map_(300, 0), points_map_(300, 0) {}
125 
126  /**
127  * Initialize object, set number of columns (quantities) in tables.
128  *
129  * Number of columns of int_vals_ table is passed by argument \p int_cols, number of columns
130  * of other tables is given by n_rows_ value.
131  */
132  void initialize(uint int_cols) {
133  this->reset();
134 
135  point_vals_.resize(n_rows_);
136  int_vals_.resize(int_cols);
137  }
138 
139  /// Reset number of columns (points and elements)
140  inline void reset() {
141  n_points_ = 0;
142  n_elems_ = 0;
143  }
144 
145  /// Getter for dim_
146  inline uint dim() const {
147  return dim_;
148  }
149 
150  /// Getter for n_rows_
151  inline uint n_rows() const {
152  return n_rows_;
153  }
154 
155  /// Getter for n_elems_
156  inline uint n_elems() const {
157  return n_elems_;
158  }
159 
160  /// Getter for n_points_
161  inline uint n_points() const {
162  return n_points_;
163  }
164 
165  /// Getter for quadrature
167  return quad_;
168  }
169 
170  /// Resize data tables. Method is called before reinit of patch.
174  }
175 
176  /**
177  * Register element, add to point_vals_ table
178  *
179  * @param coords Coordinates of element nodes.
180  * @param element_patch_idx Index of element on patch.
181  */
182  uint register_element(arma::mat coords, uint element_patch_idx) {
183  uint res_column = operations_[FeBulk::BulkOps::opElCoords].result_row();
184  for (uint i_col=0; i_col<coords.n_cols; ++i_col)
185  for (uint i_row=0; i_row<coords.n_rows; ++i_row) {
186  point_vals_(res_column)(n_elems_) = coords(i_row, i_col);
187  ++res_column;
188  }
189 
190  elements_map_[element_patch_idx] = n_elems_;
191  return n_elems_++;
192  }
193 
194  /**
195  * Register side, add to point_vals_ table
196  *
197  * @param coords Coordinates of element nodes.
198  * @param side_coords Coordinates of side nodes.
199  */
200  uint register_side(arma::mat elm_coords, arma::mat side_coords) {
201  uint res_column = operations_[FeSide::SideOps::opElCoords].result_row();
202  for (uint i_col=0; i_col<elm_coords.n_cols; ++i_col)
203  for (uint i_row=0; i_row<elm_coords.n_rows; ++i_row) {
204  point_vals_(res_column)(n_elems_) = elm_coords(i_row, i_col);
205  ++res_column;
206  }
207 
208  res_column = operations_[FeSide::SideOps::opSideCoords].result_row();
209  for (uint i_col=0; i_col<side_coords.n_cols; ++i_col)
210  for (uint i_row=0; i_row<side_coords.n_rows; ++i_row) {
211  point_vals_(res_column)(n_elems_) = side_coords(i_row, i_col);
212  ++res_column;
213  }
214 
215  return n_elems_++;
216  }
217 
218  /**
219  * Register bulk point, add to int_vals_ table
220  *
221  * @param elem_table_row Index of element in temporary element table.
222  * @param value_patch_idx Index of point in ElementCacheMap.
223  * @param elem_idx Index of element in Mesh.
224  */
225  uint register_bulk_point(uint elem_table_row, uint value_patch_idx, uint elem_idx) {
226  int_vals_(0)(n_points_) = value_patch_idx;
227  int_vals_(1)(n_points_) = elem_table_row;
228  int_vals_(2)(n_points_) = elem_idx;
229 
230  points_map_[value_patch_idx] = n_points_;
231  return n_points_++;
232  }
233 
234  /**
235  * Register side point, add to int_vals_ table
236  *
237  * @param elem_table_row Index of side in temporary element table.
238  * @param value_patch_idx Index of point in ElementCacheMap.
239  * @param elem_idx Index of element in Mesh.
240  * @param side_idx Index of side on element.
241  */
242  uint register_side_point(uint elem_table_row, uint value_patch_idx, uint elem_idx, uint side_idx) {
243  int_vals_(0)(n_points_) = value_patch_idx;
244  int_vals_(1)(n_points_) = elem_table_row;
245  int_vals_(2)(n_points_) = elem_idx;
246  int_vals_(3)(n_points_) = side_idx;
247 
248  points_map_[value_patch_idx] = n_points_;
249  return n_points_++;
250  }
251 
252  /**
253  * Adds accessor to of new operation operations_ vector
254  *
255  * @param shape Shape of function output
256  * @param reinit_f Reinitialize function
257  * @param input_ops_vec Indices of input operations in operations_ vector.
258  */
259  ElOp<spacedim> &make_new_op(std::initializer_list<uint> shape, ReinitFunction reinit_f, std::vector<uint> input_ops_vec) {
260  ElOp<spacedim> op_accessor(this->dim_, shape, this->n_rows_, reinit_f, input_ops_vec);
261  this->n_rows_ += op_accessor.n_comp();
262  operations_.push_back(op_accessor);
263  return operations_[operations_.size()-1];
264  }
265 
266  /**
267  * Adds accessor of expansion operation to operations_ vector
268  *
269  * @param el_op Source operation of expansion.
270  * @param shape Shape of function output
271  * @param reinit_f Reinitialize function
272  */
273  ElOp<spacedim> &make_expansion(ElOp<spacedim> &el_op, std::initializer_list<uint> shape, ReinitFunction reinit_f) {
274  ElOp<spacedim> op_accessor(this->dim_, shape, el_op.result_row(), reinit_f);
275  // shape passed from el_op throws:
276  // C++ exception with description "std::bad_alloc" thrown in the test body.
277  operations_.push_back(op_accessor);
278  return operations_[operations_.size()-1];
279  }
280 
281  /**
282  * Adds accessor of FE operation and adds operation dynamically to operations_ vector
283  *
284  * @param shape Shape of function output
285  * @param reinit_f Reinitialize function
286  * @param input_ops_vec Indices of input operations in operations_ vector.
287  * @param n_dofs Number of DOFs
288  */
289  ElOp<spacedim> &make_fe_op(std::initializer_list<uint> shape, ReinitFunction reinit_f, std::vector<uint> input_ops_vec, uint n_dofs) {
290  ElOp<spacedim> op_accessor(this->dim_, shape, this->n_rows_, reinit_f, input_ops_vec, n_dofs);
291  this->n_rows_ += op_accessor.n_comp() * n_dofs;
292  operations_.push_back(op_accessor);
293  return operations_[operations_.size()-1];
294  }
295 
296 
297  /**
298  * Reinitializes patch data.
299  *
300  * Calls reinit functions defined on each operations.
301  */
302  void reinit_patch() {
303  if (n_elems_ == 0) return; // skip if tables are empty
304  for (uint i=0; i<operations_.size(); ++i)
305  operations_[i].reinit_function(operations_, point_vals_, int_vals_);
306  }
307 
308  /**
309  * Returns scalar output value given by index of first row and index of quadrature point.
310  *
311  * @param result_row Row of operation in point_vals_ data table
312  * @param point_idx Index of quadrature point in ElementCacheMap
313  */
314  inline Scalar scalar_val(uint result_row, uint point_idx) const {
315  return point_vals_(result_row)(points_map_[point_idx]);
316  }
317 
318  /**
319  * Returns vector output value given by index of first row and index of quadrature point.
320  *
321  * @param result_row First row of operation in point_vals_ data table
322  * @param point_idx Index of quadrature point in ElementCacheMap
323  */
324  inline Vector vector_val(uint result_row, uint point_idx) const {
325  Vector val;
326  for (uint i=0; i<3; ++i)
327  val(i) = point_vals_(result_row+i)(points_map_[point_idx]);
328  return val;
329  }
330 
331  /**
332  * Returns tensor output value given by index of first row and index of quadrature point.
333  *
334  * @param result_row First row of operation in point_vals_ data table
335  * @param point_idx Index of quadrature point in ElementCacheMap
336  */
337  inline Tensor tensor_val(uint result_row, uint point_idx) const {
338  Tensor val;
339  for (uint i=0; i<3; ++i)
340  for (uint j=0; j<3; ++j)
341  val(i,j) = point_vals_(result_row+3*i+j)(points_map_[point_idx]);
342  return val;
343  }
344 
345  /**
346  * Performs output of data tables to stream.
347  *
348  * Development method.
349  * @param points Allows switched off output of point table,
350  * @param ints Allows switched off output of int (connectivity to elements) table,
351  */
352  void print_data_tables(ostream& stream, bool points, bool ints) const {
353  if (points) {
354  stream << "Point vals: " << point_vals_.rows() << " - " << point_vals_.cols() << std::endl;
355  for (uint i_row=0; i_row<n_points_; ++i_row) {
356  for (uint i_col=0; i_col<n_rows_; ++i_col)
357  stream << point_vals_(i_col)(i_row) << " ";
358  stream << std::endl;
359  }
360  stream << std::endl;
361  }
362  if (ints) {
363  stream << "Int vals: " << int_vals_.rows() << " - " << int_vals_.cols() << std::endl;
364  for (uint i_row=0; i_row<n_points_; ++i_row) {
365  for (uint i_col=0; i_col<3; ++i_col)
366  stream << int_vals_(i_col)(i_row) << " ";
367  stream << std::endl;
368  }
369  stream << std::endl;
370  }
371  }
372 
373  /**
374  * Performs table of fixed operations to stream.
375  *
376  * Development method.
377  * @param bulk_side Needs set 0 (bulk) or 1 (side) for correct output of operation names.
378  */
379  void print_operations(ostream& stream, uint bulk_side) const {
381  {
382  { "el_coords", "jacobian", "inv_jac", "jac_det", "exp_coords", "exp_jacobian", "exp_in_jac", "exp_jac_det", "pt_coords", "weights",
383  "JxW", "", "", "", "", "" },
384  { "el_coords", "el_jac", "el_inv_jac", "side_coords", "side_jac", "side_jac_det", "exp_el_coords", "exp_el_jac", "exp_el_inv_jac",
385  "exp_side_coords", "exp_side_jac", "exp_side_jac_det", "pt_coords", "weights", "JxW", "normal_vec", "", "", "", "", "" }
386  };
387  stream << std::setfill(' ') << " Operation" << setw(12) << "" << "Shape" << setw(2) << ""
388  << "Result row" << setw(2) << "" << "n DOFs" << setw(2) << "" << "Input operations" << endl;
389  for (uint i=0; i<operations_.size(); ++i) {
390  stream << " " << std::left << setw(20) << op_names[bulk_side][i] << "" << " " << setw(6) << operations_[i].format_shape() << "" << " "
391  << setw(11) << operations_[i].result_row() << "" << " " << setw(7) << operations_[i].n_dofs() << "" << " ";
392  auto &input_ops = operations_[i].input_ops();
393  for (auto i_o : input_ops) stream << op_names[bulk_side][i_o] << " ";
394  stream << std::endl;
395  }
396  }
397 
398 protected:
399  /**
400  * Store data of bulk or side quadrature points of one dimension
401  *
402  * Number of columns is given by n_rows_, number of used rows by n_points_.
403  */
405  /**
406  * Hold integer values of quadrature points of previous table.
407  *
408  * Table contains following columns:
409  * 0: Index of quadrature point on patch
410  * 1: Row of element/side in point_vals_ table in registration step (before expansion)
411  * 2: Element idx in Mesh
412  * 3: Index of side in element (column is allocated only for side point table)
413  * Number of used rows is given by n_points_.
414  */
416 
417  /// Vector of all defined operations
419 
420  uint dim_; ///< Dimension
421  uint n_rows_; ///< Number of columns of \p point_vals table
422  uint n_points_; ///< Number of points in patch
423  uint n_elems_; ///< Number of elements in patch
424  Quadrature *quad_; ///< Quadrature of given dimension and order passed in constructor.
425 
426  std::vector<uint> elements_map_; ///< Map of element patch indices to el_vals_ table
427  std::vector<uint> points_map_; ///< Map of point patch indices to point_vals_ and int_vals_ tables
428 
429  friend class PatchFEValues<spacedim>;
430  friend class ElOp<spacedim>;
431  template<unsigned int dim>
432  friend class BulkValues;
433  template<unsigned int dim>
434  friend class SideValues;
435  template<unsigned int dim>
436  friend class JoinValues;
437 };
438 
439 /**
440  * @brief Class represents element or FE operations.
441  */
442 template<unsigned int spacedim = 3>
443 class ElOp {
444 public:
445  /**
446  * Constructor
447  *
448  * Set all data members.
449  */
450  ElOp(uint dim, std::initializer_list<uint> shape, uint result_row, ReinitFunction reinit_f, std::vector<uint> input_ops = {}, uint n_dofs = 1)
452  {}
453 
454  /**
455  * Return number of operation components
456  *
457  * Value is computed from shape_ vector
458  */
459  inline uint n_comp() const {
460  if (shape_.size() == 1) return shape_[0];
461  else return shape_[0] * shape_[1];
462  }
463 
464  /// Getter for dimension
465  inline uint dim() const {
466  return dim_;
467  }
468 
469  /// Getter for result_row_
470  inline uint result_row() const {
471  return result_row_;
472  }
473 
474  /// Getter for n_dofs_
475  inline uint n_dofs() const {
476  return n_dofs_;
477  }
478 
479  /// Getter for input_ops_
480  inline const std::vector<uint> &input_ops() const {
481  return input_ops_;
482  }
483 
484  /// Getter for shape_
485  inline const std::vector<uint> &shape() const {
486  return shape_;
487  }
488 
489  /**
490  * Format shape to string
491  *
492  * Method is used in output development method.
493  */
494  inline std::string format_shape() const {
495  stringstream ss;
496  ss << shape_[0];
497  if (shape_.size() > 1) ss << "x" << shape_[1];
498  return ss.str();
499  }
500 
501  /// Call reinit function on element table if function is defined
502  inline void reinit_function(std::vector<ElOp<spacedim>> &operations, TableDbl &data_table, TableInt &int_table) {
503  reinit_func(operations, data_table, int_table);
504  }
505 
506  /// Return map referenced Eigen::Matrix of given dimension
507  template<unsigned int dim1, unsigned int dim2>
508  Eigen::Map<Eigen::Matrix<ArrayDbl, dim1, dim2>> value(TableDbl &op_results, uint i_dof = 0) const {
509  return Eigen::Map<Eigen::Matrix<ArrayDbl, dim1, dim2>>(op_results.data() + result_row_ + i_dof * n_comp(), dim1, dim2);
510  }
511 
512 
513 protected:
514  uint dim_; ///< Dimension
515  std::vector<uint> shape_; ///< Shape of stored data (size of vector or number of rows and cols of matrix)
516  uint result_row_; ///< First row to scalar, vector or matrix result
517  std::vector<uint> input_ops_; ///< Indices of operations in PatchPointValues::operations_ vector on which ElOp is depended
518  uint n_dofs_; ///< Number of DOFs of FE operations (or 1 in case of element operations)
519 
520  ReinitFunction reinit_func; ///< Pointer to patch reinit function of element data table specialized by operation
521 };
522 
523 
524 /// Defines common functionality of reinit operations.
526  // empty base operation
527  static inline void op_base(FMT_UNUSED std::vector<ElOp<3>> &operations, FMT_UNUSED TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
528  // empty
529  }
530 
531  // expansion operation
532  static inline void expand_data(ElOp<3> &op, TableDbl &op_results, TableInt &el_table) {
533  uint row_begin = op.result_row();
534  uint row_end = row_begin + op.n_comp();
535  uint size = op_results(row_begin).rows();
536  for (int i_pt=size-1; i_pt>=0; --i_pt) {
537  uint el_table_idx = el_table(1)(i_pt);
538  for (uint i_q=row_begin; i_q<row_end; ++i_q)
539  op_results(i_q)(i_pt) = op_results(i_q)(el_table_idx);
540  }
541  }
542 };
543 
544 /// Defines reinit operations on bulk points.
545 struct bulk_reinit {
546  // element operations
547  template<unsigned int dim>
548  static inline void elop_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
549  // result matrix(spacedim, dim), input matrix(spacedim, dim+1)
550  auto &op = operations[FeBulk::BulkOps::opJac];
551  auto jac_value = op.value<3, dim>(op_results);
552  auto coords_value = operations[ op.input_ops()[0] ].value<3, dim+1>(op_results);
553  jac_value = eigen_tools::jacobian<3,dim>(coords_value);
554  }
555  template<unsigned int dim>
556  static inline void elop_inv_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
557  // result matrix(spacedim, dim), input matrix(spacedim, dim+1)
558  auto &op = operations[FeBulk::BulkOps::opInvJac];
559  auto inv_jac_value = op.value<dim, 3>(op_results);
560  auto jac_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
561  inv_jac_value = eigen_tools::inverse<3, dim>(jac_value);
562  }
563  template<unsigned int dim>
564  static inline void elop_jac_det(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
565  // result double, input matrix(spacedim, dim)
566  auto &op = operations[FeBulk::BulkOps::opJacDet];
567  auto &jac_det_value = op_results(op.result_row());
568  auto jac_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
569  jac_det_value = eigen_tools::determinant<Eigen::Matrix<ArrayDbl, 3, dim>>(jac_value).array().abs();
570  }
571 
572  // expansion operations
573  static inline void expd_coords(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
574  auto &op = operations[FeBulk::BulkOps::opExpansionCoords];
575  common_reinit::expand_data(op, op_results, el_table);
576  }
577  static inline void expd_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
578  auto &op = operations[FeBulk::BulkOps::opExpansionJac];
579  common_reinit::expand_data(op, op_results, el_table);
580  }
581  static inline void expd_inv_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
582  auto &op = operations[FeBulk::BulkOps::opExpansionInvJac];
583  common_reinit::expand_data(op, op_results, el_table);
584  }
585  static inline void expd_jac_det(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
586  auto &op = operations[FeBulk::BulkOps::opExpansionJacDet];
587  common_reinit::expand_data(op, op_results, el_table);
588  }
589 
590  // point operations
591  static inline void ptop_coords(FMT_UNUSED std::vector<ElOp<3>> &operations, FMT_UNUSED TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
592  // Implement
593  }
594  static inline void ptop_weights(std::vector<ElOp<3>> &operations, TableDbl &op_results, std::vector<double> point_weights) {
595  auto &op = operations[FeBulk::BulkOps::opWeights];
596  ArrayDbl &result_row = op_results( op.result_row() );
597  auto n_points = point_weights.size();
598  for (uint i=0; i<result_row.rows(); ++i)
599  result_row(i) = point_weights[i%n_points];
600  }
601  static inline void ptop_JxW(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
602  auto &op = operations[FeBulk::BulkOps::opJxW];
603  ArrayDbl &weights_row = op_results( operations[op.input_ops()[0]].result_row() );
604  ArrayDbl &jac_det_row = op_results( operations[op.input_ops()[1]].result_row() );
605  ArrayDbl &result_row = op_results( op.result_row() );
606  result_row = jac_det_row * weights_row;
607  }
608  static inline void ptop_scalar_shape(std::vector<ElOp<3>> &operations, TableDbl &op_results,
609  std::vector< std::vector<double> > shape_values, uint scalar_shape_op_idx) {
610  auto &op = operations[scalar_shape_op_idx];
611  uint n_points = shape_values.size();
612 
613  for (uint i_row=0; i_row<shape_values[0].size(); ++i_row) {
614  ArrayDbl &result_row = op_results( op.result_row()+i_row );
615  for (uint i_pt=0; i_pt<result_row.rows(); ++i_pt)
616  result_row(i_pt) = shape_values[i_pt % n_points][i_row];
617  }
618  }
619  template<unsigned int dim>
620  static inline void ptop_scalar_shape_grads(std::vector<ElOp<3>> &operations, TableDbl &op_results,
621  std::vector< std::vector<arma::mat> > ref_shape_grads, uint scalar_shape_grads_op_idx) {
622  auto &op = operations[scalar_shape_grads_op_idx];
623  uint n_points = ref_shape_grads.size();
624  uint n_dofs = ref_shape_grads[0].size();
625 
626  Eigen::Vector<ArrayDbl, Eigen::Dynamic> ref_shape_grads_expd;
627  ref_shape_grads_expd.resize(dim*n_dofs);
628  for (uint i=0; i<ref_shape_grads_expd.rows(); ++i)
629  ref_shape_grads_expd(i).resize(op_results(0).rows());
630 
631  for (uint i_dof=0; i_dof<n_dofs; ++i_dof)
632  for (uint i_c=0; i_c<dim; ++i_c) {
633  ArrayDbl &shape_grad_row = ref_shape_grads_expd(i_dof*dim+i_c);
634  for (uint i_pt=0; i_pt<shape_grad_row.rows(); ++i_pt)
635  shape_grad_row(i_pt) = ref_shape_grads[i_pt % n_points][i_dof][i_c];
636  }
637 
638  auto inv_jac_value = operations[ op.input_ops()[0] ].value<dim, 3>(op_results);
639  for (uint i_dof=0; i_dof<n_dofs; ++i_dof) {
640  auto shape_grad_value = op.value<3, 1>(op_results, i_dof);
641  Eigen::Map<Eigen::Matrix<ArrayDbl, dim, 1>> ref_shape_grads_dof_value(ref_shape_grads_expd.data() + dim*i_dof, dim, 1);
642  shape_grad_value = inv_jac_value.transpose() * ref_shape_grads_dof_value;
643  }
644  }
645 };
646 
647 
648 
649 /// Defines reinit operations on side points.
650 struct side_reinit {
651  // element operations
652  template<unsigned int dim>
653  static inline void elop_el_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
654  // result matrix(spacedim, dim), input matrix(spacedim, dim+1)
655  auto &op = operations[FeSide::SideOps::opElJac];
656  auto jac_value = op.value<3, dim>(op_results);
657  auto coords_value = operations[ op.input_ops()[0] ].value<3, dim+1>(op_results);
658  jac_value = eigen_tools::jacobian<3,dim>(coords_value);
659  }
660  template<unsigned int dim>
661  static inline void elop_el_inv_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
662  auto &op = operations[FeSide::SideOps::opElInvJac];
663  auto inv_jac_value = op.value<dim, 3>(op_results);
664  auto jac_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
665  inv_jac_value = eigen_tools::inverse<3,dim>(jac_value);
666  }
667  template<unsigned int dim>
668  static inline void elop_sd_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
669  // result matrix(spacedim, dim), input matrix(spacedim, dim+1)
670  auto &op = operations[FeSide::SideOps::opSideJac];
671  auto jac_value = op.value<3, dim-1>(op_results);
672  auto coords_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
673  jac_value = eigen_tools::jacobian<3, dim-1>(coords_value);
674  }
675  template<unsigned int dim>
676  static inline void elop_sd_jac_det(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
677  // result double, input matrix(spacedim, dim)
678  auto &op = operations[FeSide::SideOps::opSideJacDet];
679  ArrayDbl &det_value = op_results( op.result_row() );
680  auto jac_value = operations[ op.input_ops()[0] ].value<3, dim-1>(op_results);
681  det_value = eigen_tools::determinant<Eigen::Matrix<ArrayDbl, 3, dim-1>>(jac_value).array().abs();
682  }
683 
684  // expansion operations
685  static inline void expd_el_coords(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
686  auto &op = operations[FeSide::SideOps::opExpansionElCoords];
687  common_reinit::expand_data(op, op_results, el_table);
688  }
689  static inline void expd_el_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
690  auto &op = operations[FeSide::SideOps::opExpansionElJac];
691  common_reinit::expand_data(op, op_results, el_table);
692  }
693  static inline void expd_el_inv_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
694  auto &op = operations[FeSide::SideOps::opExpansionElInvJac];
695  common_reinit::expand_data(op, op_results, el_table);
696  }
697 
698  static inline void expd_sd_coords(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
699  auto &op = operations[FeSide::SideOps::opExpansionSideCoords];
700  common_reinit::expand_data(op, op_results, el_table);
701  }
702  template<unsigned int dim>
703  static inline void expd_sd_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
704  auto &op = operations[FeSide::SideOps::opExpansionSideJac];
705  common_reinit::expand_data(op, op_results, el_table);
706  }
707  static inline void expd_sd_jac_det(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
708  auto &op = operations[FeSide::SideOps::opExpansionSideJacDet];
709  common_reinit::expand_data(op, op_results, el_table);
710  }
711 
712  // Point operations
713  static inline void ptop_coords(FMT_UNUSED std::vector<ElOp<3>> &operations, FMT_UNUSED TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
714  // Implement
715  }
716  static inline void ptop_weights(std::vector<ElOp<3>> &operations, TableDbl &op_results, std::vector<double> point_weights) {
717  auto &op = operations[FeSide::SideOps::opWeights];
718  ArrayDbl &result_row = op_results( op.result_row() );
719  auto n_points = point_weights.size();
720  for (uint i=0; i<result_row.rows(); ++i)
721  result_row(i) = point_weights[i%n_points];
722  }
723  static inline void ptop_JxW(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
724  auto &op = operations[FeSide::SideOps::opJxW];
725  ArrayDbl &weights_row = op_results( operations[op.input_ops()[0]].result_row() );
726  ArrayDbl &jac_det_row = op_results( operations[op.input_ops()[1]].result_row() );
727  ArrayDbl &result_row = op_results( op.result_row() );
728  result_row = jac_det_row * weights_row;
729  }
730  template<unsigned int dim>
731  static inline void ptop_normal_vec(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
732  auto &op = operations[FeSide::SideOps::opNormalVec];
733  auto normal_value = op.value<3, 1>(op_results);
734  auto inv_jac_mat_value = operations[ op.input_ops()[0] ].value<dim, 3>(op_results);
735  normal_value = inv_jac_mat_value.transpose() * RefElement<dim>::normal_vector_array( el_table(3) );
736 
737  ArrayDbl norm_vec;
738  norm_vec.resize(normal_value(0).rows());
739  Eigen::VectorXd A(3);
740 
741  for (uint i=0; i<normal_value(0).rows(); ++i) {
742  A(0) = normal_value(0)(i);
743  A(1) = normal_value(1)(i);
744  A(2) = normal_value(2)(i);
745  norm_vec(i) = A.norm();
746  }
747  for (uint i=0; i<3; ++i) {
748  normal_value(i) /= norm_vec;
749  }
750  }
751  static inline void ptop_scalar_shape(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table,
752  std::vector< std::vector< std::vector<double> > > shape_values, uint scalar_shape_op_idx) {
753  auto &op = operations[scalar_shape_op_idx];
754  uint n_points = shape_values[0].size();
755 
756  for (uint i_row=0; i_row<shape_values[0][0].size(); ++i_row) {
757  ArrayDbl &result_row = op_results( op.result_row()+i_row );
758  for (uint i_pt=0; i_pt<result_row.rows(); ++i_pt)
759  result_row(i_pt) = shape_values[el_table(3)(i_pt)][i_pt % n_points][i_row];
760  }
761  }
762  template<unsigned int dim>
763  static inline void ptop_scalar_shape_grads(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table,
764  std::vector< std::vector< std::vector<arma::mat> > > ref_shape_grads, uint scalar_shape_grads_op_idx) {
765  auto &op = operations[scalar_shape_grads_op_idx];
766  uint n_points = ref_shape_grads[0].size();
767  uint n_dofs = ref_shape_grads[0][0].size();
768 
769  Eigen::Vector<ArrayDbl, Eigen::Dynamic> ref_shape_grads_expd;
770  ref_shape_grads_expd.resize(dim*n_dofs);
771  for (uint i=0; i<ref_shape_grads_expd.rows(); ++i)
772  ref_shape_grads_expd(i).resize(op_results(0).rows());
773 
774  for (uint i_dof=0; i_dof<n_dofs; ++i_dof)
775  for (uint i_c=0; i_c<dim; ++i_c) {
776  ArrayDbl &shape_grad_row = ref_shape_grads_expd(i_dof*dim+i_c);
777  for (uint i_pt=0; i_pt<shape_grad_row.rows(); ++i_pt)
778  shape_grad_row(i_pt) = ref_shape_grads[el_table(3)(i_pt)][i_pt % n_points][i_dof][i_c];
779  }
780 
781  auto inv_jac_value = operations[ op.input_ops()[0] ].value<dim, 3>(op_results);
782  for (uint i_dof=0; i_dof<n_dofs; ++i_dof) {
783  auto shape_grad_value = op.value<3, 1>(op_results, i_dof);
784  Eigen::Map<Eigen::Matrix<ArrayDbl, dim, 1>> ref_shape_grads_dof_value(ref_shape_grads_expd.data() + dim*i_dof, dim, 1);
785  shape_grad_value = inv_jac_value.transpose() * ref_shape_grads_dof_value;
786  }
787  }
788 };
789 
790 
791 // template specialization
792 template<>
793 inline void side_reinit::elop_sd_jac<1>(FMT_UNUSED std::vector<ElOp<3>> &operations, FMT_UNUSED TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
794 }
795 
796 template<>
797 inline void side_reinit::elop_sd_jac_det<1>(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
798  auto &op = operations[FeSide::SideOps::opSideJacDet];
799  ArrayDbl &result_vec = op_results( op.result_row() );
800  for (uint i=0;i<result_vec.size(); ++i) {
801  result_vec(i) = 1.0;
802  }
803 }
804 
805 template<>
806 inline void side_reinit::expd_sd_jac<1>(FMT_UNUSED std::vector<ElOp<3>> &operations, FMT_UNUSED TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
807 }
808 
809 
810 
811 namespace FeBulk {
812 
813  /// Bulk data specialization, order of item in operations_ vector corresponds to the BulkOps enum
814  template<unsigned int spacedim = 3>
815  class PatchPointValues : public ::PatchPointValues<spacedim> {
816  public:
817  /// Constructor
819  : ::PatchPointValues<spacedim>(dim) {
820  this->quad_ = new QGauss(dim, 2*quad_order);
821  switch (dim) {
822  case 1:
823  init<1>();
824  break;
825  case 2:
826  init<2>();
827  break;
828  case 3:
829  init<3>();
830  break;
831  }
832  }
833 
834  private:
835  /// Initialize operations vector
836  template<unsigned int dim>
837  void init() {
838  // First step: adds element values operations
839  auto &el_coords = this->make_new_op( {spacedim, this->dim_+1}, &common_reinit::op_base, {} );
840 
841  auto &el_jac = this->make_new_op( {spacedim, this->dim_}, &bulk_reinit::elop_jac<dim>, {BulkOps::opElCoords} );
842 
843  auto &el_inv_jac = this->make_new_op( {this->dim_, spacedim}, &bulk_reinit::elop_inv_jac<dim>, {BulkOps::opJac} );
844 
845  auto &el_jac_det = this->make_new_op( {1}, &bulk_reinit::elop_jac_det<dim>, {BulkOps::opJac} );
846 
847  // Second step: adds expand operations (element values to point values)
848  this->make_expansion( el_coords, {spacedim, this->dim_+1}, &bulk_reinit::expd_coords );
849 
850  this->make_expansion( el_jac, {spacedim, this->dim_}, &bulk_reinit::expd_jac );
851 
852  this->make_expansion( el_inv_jac, {this->dim_, spacedim}, &bulk_reinit::expd_inv_jac );
853 
854  this->make_expansion( el_jac_det, {1}, &bulk_reinit::expd_jac_det );
855 
856  // Third step: adds point values operations
857  /*auto &pt_coords =*/ this->make_new_op( {spacedim}, &bulk_reinit::ptop_coords, {} );
858 
859  // use lambda reinit function
860  std::vector<double> point_weights = this->quad_->get_weights();
861  auto lambda_weights = [point_weights](std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
862  bulk_reinit::ptop_weights(operations, op_results, point_weights);
863  };
864  /*auto &weights =*/ this->make_new_op( {1}, lambda_weights, {} );
865 
866  /*auto &JxW =*/ this->make_new_op( {1}, &bulk_reinit::ptop_JxW, {BulkOps::opWeights, BulkOps::opJacDet} );
867  }
868  };
869 
870 } // closing namespace FeBulk
871 
872 
873 
874 namespace FeSide {
875 
876 /// Bulk Side specialization, order of item in operations_ vector corresponds to the SideOps enum
877  template<unsigned int spacedim = 3>
878  class PatchPointValues : public ::PatchPointValues<spacedim> {
879  public:
880  /// Constructor
882  : ::PatchPointValues<spacedim>(dim) {
883  this->quad_ = new QGauss(dim-1, 2*quad_order);
884  switch (dim) {
885  case 1:
886  init<1>();
887  break;
888  case 2:
889  init<2>();
890  break;
891  case 3:
892  init<3>();
893  break;
894  }
895  }
896 
897  private:
898  /// Initialize operations vector
899  template<unsigned int dim>
900  void init() {
901  // First step: adds element values operations
902  auto &el_coords = this->make_new_op( {spacedim, this->dim_+1}, &common_reinit::op_base, {} );
903 
904  auto &el_jac = this->make_new_op( {spacedim, this->dim_}, &side_reinit::elop_el_jac<dim>, {SideOps::opElCoords} );
905 
906  auto &el_inv_jac = this->make_new_op( {this->dim_, spacedim}, &side_reinit::elop_el_inv_jac<dim>, {SideOps::opElJac} );
907 
908  auto &sd_coords = this->make_new_op( {spacedim, this->dim_}, &common_reinit::op_base, {} );
909 
910  auto &sd_jac = this->make_new_op( {spacedim, this->dim_-1}, &side_reinit::elop_sd_jac<dim>, {SideOps::opSideCoords} );
911 
912  auto &sd_jac_det = this->make_new_op( {1}, &side_reinit::elop_sd_jac_det<dim>, {SideOps::opSideJac} );
913 
914  // Second step: adds expand operations (element values to point values)
915  this->make_expansion( el_coords, {spacedim, this->dim_+1}, &side_reinit::expd_el_coords );
916 
917  this->make_expansion( el_jac, {spacedim, this->dim_}, &side_reinit::expd_el_jac );
918 
919  this->make_expansion( el_inv_jac, {this->dim_, spacedim}, &side_reinit::expd_el_inv_jac );
920 
921  this->make_expansion( sd_coords, {spacedim, this->dim_}, &side_reinit::expd_sd_coords );
922 
923  this->make_expansion( sd_jac, {spacedim, this->dim_-1}, &side_reinit::expd_sd_jac<dim> );
924 
925  this->make_expansion( sd_jac_det, {1}, &side_reinit::expd_sd_jac_det );
926 
927  // Third step: adds point values operations
928  /*auto &coords =*/ this->make_new_op( {spacedim}, &side_reinit::ptop_coords, {} );
929 
930  // use lambda reinit function
931  std::vector<double> point_weights = this->quad_->get_weights();
932  auto lambda_weights = [point_weights](std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
933  side_reinit::ptop_weights(operations, op_results, point_weights);
934  };
935  /*auto &weights =*/ this->make_new_op( {1}, lambda_weights, {} );
936 
938 
939  /*auto &normal_vec =*/ this->make_new_op( {spacedim}, &side_reinit::ptop_normal_vec<dim>, {SideOps::opElInvJac} );
940  }
941  };
942 
943 } // closing namespace FeSide
944 
945 
946 
947 #endif /* PATCH_POINT_VALUES_HH_ */
Class represents element or FE operations.
std::string format_shape() const
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 for result_row_.
const std::vector< uint > & shape() const
Getter for shape_.
const std::vector< uint > & input_ops() const
Getter for input_ops_.
uint result_row_
First row to scalar, vector or matrix result.
uint n_dofs() const
Getter for n_dofs_.
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)
uint n_dofs_
Number of DOFs of FE operations (or 1 in case of element operations)
uint dim() const
Getter for 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
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.
void print_data_tables(ostream &stream, bool points, bool ints) const
ElOp< spacedim > & make_expansion(ElOp< spacedim > &el_op, std::initializer_list< uint > shape, ReinitFunction reinit_f)
void initialize(uint int_cols)
void resize_tables(uint n_points)
Resize data tables. Method is called before reinit of patch.
uint n_rows() const
Getter for 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)
uint n_rows_
Number of columns of point_vals table.
uint register_element(arma::mat coords, uint element_patch_idx)
Tensor tensor_val(uint result_row, uint point_idx) const
uint register_side(arma::mat elm_coords, arma::mat side_coords)
std::vector< ElOp< spacedim > > operations_
Vector of all defined operations.
Quadrature * get_quadrature() const
Getter for quadrature.
uint n_elems() const
Getter for n_elems_.
void reset()
Reset number of columns (points and elements)
uint n_points() const
Getter for n_points_.
std::vector< uint > elements_map_
Map of element patch indices to el_vals_ table.
uint dim() const
Getter for dim_.
ElOp< spacedim > & make_fe_op(std::initializer_list< uint > shape, ReinitFunction reinit_f, std::vector< uint > input_ops_vec, uint n_dofs)
Vector vector_val(uint result_row, uint point_idx) const
uint n_points_
Number of points in patch.
void print_operations(ostream &stream, uint bulk_side) const
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)
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)
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
Array< double > array
Definition: armor.hh:938
@ 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.
@ 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 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 expd_sd_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, 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)