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