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