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