Flow123d  DF_patch_fe_mechanics-b866f5c
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>> &, IntTableArena &)>;
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  /// operations evaluated on quadrature points
91  opCoords, ///< coordinations of quadrature point
92  opWeights, ///< weight of quadrature point
93  opJxW, ///< JxW value of quadrature point
94  opNormalVec ///< normal vector of quadrature point
95  };
96 }
97 
98 
99 /**
100  * Enumeration of FE operations.
101  *
102  * PatchFEValues::fe_ops_indices_ vector holds indices of FE operations
103  * in operations_ vector or value of invalid operation if operations is
104  * not defined.
105  */
106 enum FEOps
107 {
108  opScalarShape, ///< Scalar shape operation
109  opVectorShape, ///< Vector shape operation
110  opGradScalarShape, ///< Scalar shape gradient
111  opGradVectorShape, ///< Vector shape gradient
112  opVectorSymGrad, ///< Vector symmetric gradient
113  opVectorDivergence, ///< Vector divergence
114  OpNItems ///< Holds number of valid FE operations and value of invalid FE operation
115 };
116 
117 
118 /// Distinguishes operations by type and size of output rows
120 {
121  elemOp, ///< operation is evaluated on elements or sides
122  pointOp, ///< operation is evaluated on quadrature points
123  fixedSizeOp ///< operation has fixed size and it is filled during initialization
124 };
125 
126 
127 
128 /**
129  * v Class for storing FE data of quadrature points on one patch.
130  *
131  * Store data of bulk or side quadrature points of one dimension.
132  */
133 template<unsigned int spacedim = 3>
135 {
136 public:
137  /**
138  * Constructor
139  *
140  * @param dim Set dimension
141  */
144  points_map_(300, 0), asm_arena_(asm_arena) {
145  reset();
146  }
147 
148  /**
149  * Initialize object, set number of columns (quantities) in tables.
150  */
151  void initialize() {
152  this->reset();
153  int_table_.resize(int_sizes_.size());
154  }
155 
156  inline void init_finalize(PatchArena *patch_arena) {
157  patch_arena_ = patch_arena;
158  }
159 
160  /// Reset number of columns (points and elements)
161  inline void reset() {
162  n_points_ = 0;
163  n_elems_ = 0;
164  i_elem_ = 0;
165  }
166 
167  /// Getter for dim_
168  inline uint dim() const {
169  return dim_;
170  }
171 
172  /// Getter for n_elems_
173  inline uint n_elems() const {
174  return n_elems_;
175  }
176 
177  /// Getter for n_points_
178  inline uint n_points() const {
179  return n_points_;
180  }
181 
182  /// Getter for quadrature
184  return quad_;
185  }
186 
187  /// Resize data tables. Method is called before reinit of patch.
189  n_elems_ = n_elems;
192  for (uint i=0; i<int_table_.rows(); ++i) {
193  int_table_(i) = ArenaVec<uint>(sizes[ int_sizes_[i] ], *patch_arena_);
194  }
195  for (auto &elOp : operations_)
196  if (elOp.size_type() != fixedSizeOp) {
197  elOp.allocate_result(sizes[elOp.size_type()], *patch_arena_);
198  }
199  std::fill(elements_map_.begin(), elements_map_.end(), (uint)-1);
200  }
201 
202  /**
203  * Register element, add to coords operation
204  *
205  * @param coords Coordinates of element nodes.
206  * @param element_patch_idx Index of element on patch.
207  */
208  uint register_element(arma::mat coords, uint element_patch_idx) {
209  if (elements_map_[element_patch_idx] != (uint)-1) {
210  // Return index of element on patch if it is registered repeatedly
211  return elements_map_[element_patch_idx];
212  }
213 
215  auto &coords_mat = op.result_matrix();
216  std::size_t i_elem = i_elem_;
217  for (uint i_col=0; i_col<coords.n_cols; ++i_col)
218  for (uint i_row=0; i_row<coords.n_rows; ++i_row) {
219  coords_mat(i_row, i_col)(i_elem) = coords(i_row, i_col);
220  }
221 
222  elements_map_[element_patch_idx] = i_elem_;
223  return i_elem_++;
224  }
225 
226  /**
227  * Register side, add to coords operations
228  *
229  * @param coords Coordinates of element nodes.
230  * @param side_coords Coordinates of side nodes.
231  * @param side_idx Index of side on element.
232  */
233  uint register_side(arma::mat elm_coords, arma::mat side_coords, uint side_idx) {
234  {
236  auto &coords_mat = op.result_matrix();
237  std::size_t i_elem = i_elem_;
238  for (uint i_col=0; i_col<elm_coords.n_cols; ++i_col)
239  for (uint i_row=0; i_row<elm_coords.n_rows; ++i_row) {
240  coords_mat(i_row, i_col)(i_elem) = elm_coords(i_row, i_col);
241  }
242  }
243 
244  {
246  auto &coords_mat = op.result_matrix();
247  std::size_t i_elem = i_elem_;
248  for (uint i_col=0; i_col<side_coords.n_cols; ++i_col)
249  for (uint i_row=0; i_row<side_coords.n_rows; ++i_row) {
250  coords_mat(i_row, i_col)(i_elem) = side_coords(i_row, i_col);
251  }
252  }
253 
254  int_table_(3)(i_elem_) = side_idx;
255 
256  return i_elem_++;
257  }
258 
259  /**
260  * Register bulk point, add to int_table_
261  *
262  * @param elem_table_row Index of element in temporary element table.
263  * @param value_patch_idx Index of point in ElementCacheMap.
264  * @param elem_idx Index of element in Mesh.
265  * @param i_point_on_elem Index of point on element
266  */
267  uint register_bulk_point(uint elem_table_row, uint value_patch_idx, uint elem_idx, uint i_point_on_elem) {
268  uint point_pos = i_point_on_elem * n_elems_ + elem_table_row; // index of bulk point on patch
269  int_table_(0)(point_pos) = value_patch_idx;
270  int_table_(1)(point_pos) = elem_table_row;
271  int_table_(2)(point_pos) = elem_idx;
272 
273  points_map_[value_patch_idx] = point_pos;
274  return point_pos;
275  }
276 
277  /**
278  * Register side point, add to int_table_
279  *
280  * @param elem_table_row Index of side in temporary element table.
281  * @param value_patch_idx Index of point in ElementCacheMap.
282  * @param elem_idx Index of element in Mesh.
283  * @param side_idx Index of side on element.
284  * @param i_point_on_side Index of point on side
285  */
286  uint register_side_point(uint elem_table_row, uint value_patch_idx, uint elem_idx, uint side_idx, uint i_point_on_side) {
287  uint point_pos = i_point_on_side * n_elems_ + elem_table_row; // index of side point on patch
288  int_table_(0)(point_pos) = value_patch_idx;
289  int_table_(1)(point_pos) = elem_table_row;
290  int_table_(2)(point_pos) = elem_idx;
291  int_table_(4)(point_pos) = side_idx;
292 
293  points_map_[value_patch_idx] = point_pos;
294  return point_pos;
295  }
296 
297  /**
298  * Adds accessor of new operation to operations_ vector
299  *
300  * @param shape Shape of function output
301  * @param reinit_f Reinitialize function
302  * @param input_ops_vec Indices of input operations in operations_ vector.
303  * @param size_type Type of operation by size of rows
304  */
305  ElOp<spacedim> &make_new_op(std::initializer_list<uint> shape, ReinitFunction reinit_f, std::vector<uint> input_ops_vec, OpSizeType size_type = pointOp) {
306  ElOp<spacedim> op_accessor(this->dim_, shape, reinit_f, size_type, input_ops_vec);
307  row_sizes_.insert(row_sizes_.end(), op_accessor.n_comp(), size_type);
308  operations_.push_back(op_accessor);
309  return operations_[operations_.size()-1];
310  }
311 
312  /**
313  * Adds accessor of new operation with fixed data size (ref data) to operations_ vector
314  *
315  * @param shape Shape of function output
316  * @param reinit_f Reinitialize function
317  */
318  ElOp<spacedim> &make_fixed_op(std::initializer_list<uint> shape, ReinitFunction reinit_f) {
319  return make_new_op(shape, reinit_f, {}, fixedSizeOp);
320  }
321 
322  /**
323  * Adds accessor of FE operation and adds operation dynamically to operations_ vector
324  *
325  * @param shape Shape of function output
326  * @param reinit_f Reinitialize function
327  * @param input_ops_vec Indices of input operations in operations_ vector.
328  * @param n_dofs Number of DOFs
329  * @param size_type Type of operation by size of rows
330  */
331  ElOp<spacedim> &make_fe_op(std::initializer_list<uint> shape, ReinitFunction reinit_f, std::vector<uint> input_ops_vec, uint n_dofs,
332  OpSizeType size_type = pointOp) {
333  ElOp<spacedim> op_accessor(this->dim_, shape, reinit_f, size_type, input_ops_vec, n_dofs);
334  row_sizes_.insert(row_sizes_.end(), op_accessor.n_comp() * n_dofs, size_type);
335  operations_.push_back(op_accessor);
336  return operations_[operations_.size()-1];
337  }
338 
339 
340  /**
341  * Reinitializes patch data.
342  *
343  * Calls reinit functions defined on each operations.
344  */
345  void reinit_patch() {
346  if (n_elems_ == 0) return; // skip if tables are empty
347  for (uint i=0; i<operations_.size(); ++i)
348  operations_[i].reinit_function(operations_, int_table_);
349  }
350 
351  /**
352  * Returns scalar output value.
353  *
354  * @param op_idx Index of operation in operations vector
355  * @param point_idx Index of quadrature point in ElementCacheMap
356  * @param i_dof Index of DOF
357  */
358  inline Scalar scalar_value(uint op_idx, uint point_idx, uint i_dof=0) const {
359  return operations_[op_idx].result_matrix()(0)(points_map_[point_idx] + i_dof*n_points_);
360  }
361 
362  /**
363  * Returns vector output value.
364  *
365  * @param op_idx Index of operation in operations vector
366  * @param point_idx Index of quadrature point in ElementCacheMap
367  * @param i_dof Index of DOF
368  */
369  inline Vector vector_value(uint op_idx, uint point_idx, uint i_dof=0) const {
370  Vector val;
371  const auto &op_matrix = operations_[op_idx].result_matrix();
372  uint op_matrix_idx = points_map_[point_idx] + i_dof*n_points_;
373  for (uint i=0; i<3; ++i)
374  val(i) = op_matrix(i)(op_matrix_idx);
375  return val;
376  }
377 
378  /**
379  * Returns tensor output value.
380  *
381  * @param op_idx Index of operation in operations vector
382  * @param point_idx Index of quadrature point in ElementCacheMap
383  * @param i_dof Index of DOF
384  */
385  inline Tensor tensor_value(uint op_idx, uint point_idx, uint i_dof=0) const {
386  Tensor val;
387  const auto &op_matrix = operations_[op_idx].result_matrix();
388  uint op_matrix_idx = points_map_[point_idx] + i_dof*n_points_;
389  for (uint i=0; i<3; ++i)
390  for (uint j=0; j<3; ++j)
391  val(i,j) = op_matrix(i,j)(op_matrix_idx);
392  return val;
393  }
394 
395  /// Set index of FE operation in operations_ vector.
396  inline void set_fe_op(FEOps fe_op, unsigned int op_vec_idx) {
397  fe_ops_indices_[fe_op] = op_vec_idx;
398  }
399 
400  /// Return index of FE operation in operations_ vector.
401  inline unsigned int get_fe_op(FEOps fe_op) const {
402  return fe_ops_indices_[fe_op];
403  }
404 
405  /**
406  * Performs output of data tables to stream.
407  *
408  * Development method.
409  * @param points Allows switched off output of point table,
410  * @param ints Allows switched off output of int (connectivity to elements) table,
411  */
412  void print_data_tables(ostream& stream, bool points, bool ints) const {
413  if (points) {
414  stream << "Point vals: " << std::endl;
415  for (auto &op : operations_) {
416  const auto &mat = op.result_matrix();
417  for (uint i_mat=0; i_mat<mat.rows()*mat.cols(); ++i_mat) {
418  if (mat(i_mat).data_size()==0) stream << "<empty>";
419  else {
420  const double *vals = mat(i_mat).data_ptr();
421  for (size_t i_val=0; i_val<mat(i_mat).data_size(); ++i_val)
422  stream << vals[i_val] << " ";
423  }
424  stream << std::endl;
425  }
426  stream << " --- end of operation ---" << std::endl;
427  }
428  }
429  if (ints) {
430  stream << "Int vals: " << int_table_.rows() << " - " << int_table_.cols() << std::endl;
431  for (uint i_row=0; i_row<int_table_.rows(); ++i_row) {
432  if (int_table_(i_row).data_size()==0) stream << "<empty>";
433  else {
434  const uint *vals = int_table_(i_row).data_ptr();
435  for (size_t i_val=0; i_val<int_table_(i_row).data_size(); ++i_val)
436  stream << vals[i_val] << " ";
437  }
438  stream << std::endl;
439  }
440  stream << std::endl;
441  }
442  }
443 
444  /**
445  * Performs table of fixed operations to stream.
446  *
447  * Development method.
448  * @param bulk_side Needs set 0 (bulk) or 1 (side) for correct output of operation names.
449  */
450  void print_operations(ostream& stream, uint bulk_side) const {
452  {
453  { "el_coords", "jacobian", "inv_jac", "jac_det", "pt_coords", "weights", "JxW", "", "", "", "", "" },
454  { "el_coords", "el_jac", "el_inv_jac", "side_coords", "side_jac", "side_jac_det", "exp_el_coords", "exp_el_jac", "exp_el_inv_jac",
455  "exp_side_coords", "exp_side_jac", "exp_side_jac_det", "pt_coords", "weights", "JxW", "normal_vec", "", "", "", "", "" }
456  };
457  stream << std::setfill(' ') << " Operation" << setw(12) << "" << "Shape" << setw(2) << ""
458  << "n DOFs" << setw(2) << "" << "Input operations" << endl;
459  for (uint i=0; i<operations_.size(); ++i) {
460  stream << " " << std::left << setw(20) << op_names[bulk_side][i] << "" << " " << setw(6) << operations_[i].format_shape() << "" << " "
461  << setw(7) << operations_[i].n_dofs() << "" << " ";
462  auto &input_ops = operations_[i].input_ops();
463  for (auto i_o : input_ops) stream << op_names[bulk_side][i_o] << " ";
464  stream << std::endl;
465  }
466  }
467 
468 protected:
469  /**
470  * Hold integer values of quadrature points of defined operations.
471  *
472  * Table contains following rows:
473  * 0: Index of quadrature point on patch
474  * 1: Row of element/side in ElOp::result_ table in registration step (before expansion)
475  * 2: Element idx in Mesh
476  * - last two rows are allocated only for side point table
477  * 3: Index of side in element - short vector, size of column = number of sides
478  * 4: Index of side in element - long vector, size of column = number of points
479  * Number of used rows is given by n_points_.
480  */
482 
483  /// Set size and type of rows of int_table_, value is set implicitly in constructor of descendants
485 
486 
487  /// Vector of all defined operations
489 
490  /// Indices of FE operations in operations_ vector
492 
493  uint dim_; ///< Dimension
494  uint n_points_; ///< Number of points in patch
495  uint n_elems_; ///< Number of elements in patch
496  uint i_elem_; ///< Index of registered element in table, helper value used during patch creating.
497  Quadrature *quad_; ///< Quadrature of given dimension and order passed in constructor.
498 
499  std::vector<uint> elements_map_; ///< Map of element patch indices to ElOp::result_ and int_table_ tables
500  std::vector<uint> points_map_; ///< Map of point patch indices to ElOp::result_ and int_table_ tables
501  std::vector<OpSizeType> row_sizes_; ///< hold sizes of rows by type of operation
502  AssemblyArena &asm_arena_; ///< Reference to global assembly arena of PatchFeValues
503  PatchArena *patch_arena_; ///< Pointer to global patch arena of PatchFeValues
504 
505  friend class PatchFEValues<spacedim>;
506  friend class ElOp<spacedim>;
507  template<unsigned int dim>
508  friend class BulkValues;
509  template<unsigned int dim>
510  friend class SideValues;
511  template<unsigned int dim>
512  friend class JoinValues;
513 };
514 
515 /**
516  * @brief Class represents element or FE operations.
517  */
518 template<unsigned int spacedim = 3>
519 class ElOp {
520 public:
521  /**
522  * Constructor
523  *
524  * Set all data members.
525  */
526  ElOp(uint dim, std::initializer_list<uint> shape, ReinitFunction reinit_f, OpSizeType size_type, std::vector<uint> input_ops = {}, uint n_dofs = 1)
528  {}
529 
530  /// Aligns shape_vec to 2 items (equal to matrix number of dimensions)
531  std::vector<uint> set_shape_vec(std::initializer_list<uint> shape) const {
532  std::vector<uint> shape_vec(shape);
533  if (shape_vec.size() == 1) shape_vec.push_back(1);
534  ASSERT_EQ(shape_vec.size(), 2);
535  return shape_vec;
536  }
537 
538  /**
539  * Return number of operation components
540  *
541  * Value is computed from shape_ vector
542  */
543  inline uint n_comp() const {
544  return shape_[0] * shape_[1];
545  }
546 
547  /// Getter for dimension
548  inline uint dim() const {
549  return dim_;
550  }
551 
552  /// Getter for size_type_
554  return size_type_;
555  }
556 
557  /// Getter for n_dofs_
558  inline uint n_dofs() const {
559  return n_dofs_;
560  }
561 
562  /// Getter for input_ops_
563  inline const std::vector<uint> &input_ops() const {
564  return input_ops_;
565  }
566 
567  /// Getter for shape_
568  inline const std::vector<uint> &shape() const {
569  return shape_;
570  }
571 
572  /**
573  * Format shape to string
574  *
575  * Method is used in output development method.
576  */
577  inline std::string format_shape() const {
578  stringstream ss;
579  ss << shape_[0] << "x" << shape_[1];
580  return ss.str();
581  }
582 
583  /// Call reinit function on element table if function is defined
584  inline void reinit_function(std::vector<ElOp<spacedim>> &operations, IntTableArena &int_table) {
585  reinit_func(operations, int_table);
586  }
587 
588  inline void allocate_result(size_t data_size, PatchArena &arena) {
589  result_ = Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic>(shape_[0], shape_[1]);
590  for (uint i=0; i<n_comp(); ++i)
591  result_(i) = ArenaVec<double>(data_size, arena);
592  }
593 
594  /// Return map referenced result as Eigen::Matrix
595  Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> &result_matrix() {
596  return result_;
597  }
598 
599  /// Same as previous but return const reference
600  const Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> &result_matrix() const {
601  return result_;
602  }
603 
604 
605 protected:
606  uint dim_; ///< Dimension
607  std::vector<uint> shape_; ///< Shape of stored data (size of vector or number of rows and cols of matrix)
608  Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> result_; ///< Result matrix of operation
609  OpSizeType size_type_; ///< Type of operation by size of vector (element, point or fixed size)
610  std::vector<uint> input_ops_; ///< Indices of operations in PatchPointValues::operations_ vector on which ElOp is depended
611  uint n_dofs_; ///< Number of DOFs of FE operations (or 1 in case of element operations)
612 
613  ReinitFunction reinit_func; ///< Pointer to patch reinit function of element data table specialized by operation
614 };
615 
616 
617 /// Defines common functionality of reinit operations.
619  // empty base operation
620  static inline void op_base(FMT_UNUSED std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
621  // empty
622  }
623 
624  /// Common reinit function of vector symmetric gradient on bulk and side points
625  static inline void ptop_vector_sym_grad(std::vector<ElOp<3>> &operations, uint vector_sym_grad_op_idx) {
626  auto &op = operations[vector_sym_grad_op_idx];
627  auto &grad_vector_value = operations[ op.input_ops()[0] ].result_matrix();
628  auto &sym_grad_value = op.result_matrix();
629  sym_grad_value = 0.5 * (grad_vector_value.transpose() + grad_vector_value);
630  }
631  /// Common reinit function of vector divergence on bulk and side points
632  static inline void ptop_vector_divergence(std::vector<ElOp<3>> &operations, uint vector_divergence_op_idx) {
633  auto &op = operations[vector_divergence_op_idx];
634  auto &grad_vector_value = operations[ op.input_ops()[0] ].result_matrix();
635  auto &divergence_value = op.result_matrix();
636  divergence_value(0,0) = grad_vector_value(0,0) + grad_vector_value(1,1) + grad_vector_value(2,2);
637  }
638 };
639 
640 /// Defines reinit operations on bulk points.
641 struct bulk_reinit {
642  // element operations
643  template<unsigned int dim>
644  static inline void elop_jac(std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
645  // result matrix(spacedim, dim), input matrix(spacedim, dim+1)
646  auto &op = operations[FeBulk::BulkOps::opJac];
647  auto &jac_value = op.result_matrix();
648  const auto &coords_value = operations[ op.input_ops()[0] ].result_matrix();
649  for (unsigned int i=0; i<3; i++)
650  for (unsigned int j=0; j<dim; j++)
651  jac_value(i,j) = coords_value(i,j+1) - coords_value(i,0);
652  }
653  template<unsigned int dim>
654  static inline void elop_inv_jac(std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
655  // result matrix(spacedim, dim), input matrix(spacedim, dim+1)
656  auto &op = operations[FeBulk::BulkOps::opInvJac];
657  auto &inv_jac_value = op.result_matrix();
658  const auto &jac_value = operations[ op.input_ops()[0] ].result_matrix();
659  inv_jac_value = eigen_arena_tools::inverse<3, dim>(jac_value);
660  }
661  template<unsigned int dim>
662  static inline void elop_jac_det(std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
663  // result double, input matrix(spacedim, dim)
664  auto &op = operations[FeBulk::BulkOps::opJacDet];
665  auto &jac_det_value = op.result_matrix();
666  const auto &jac_value = operations[ op.input_ops()[0] ].result_matrix();
667  jac_det_value(0,0) = eigen_arena_tools::determinant<3, dim>(jac_value).abs();
668  }
669 
670  // point operations
671  static inline void ptop_coords(FMT_UNUSED std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
672  // Implement
673  }
674  static inline void ptop_weights(std::vector<ElOp<3>> &operations, PatchArena *arena, const std::vector<double> &point_weights) {
675  auto &op = operations[FeBulk::BulkOps::opWeights];
676  op.allocate_result(point_weights.size(), *arena);
677  auto &weights_value = op.result_matrix();
678  for (uint i=0; i<point_weights.size(); ++i)
679  weights_value(0,0)(i) = point_weights[i];
680  }
681  static inline void ptop_JxW(std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
682  auto &op = operations[FeBulk::BulkOps::opJxW];
683  auto &weights_value = operations[ op.input_ops()[0] ].result_matrix();
684  auto &jac_det_value = operations[ op.input_ops()[1] ].result_matrix();
685  ArenaOVec<double> weights_ovec( weights_value(0,0) );
686  ArenaOVec<double> jac_det_ovec( jac_det_value(0,0) );
687  ArenaOVec<double> jxw_ovec = jac_det_ovec * weights_ovec;
688  auto &jxw_value = op.result_matrix();
689  jxw_value(0,0) = jxw_ovec.get_vec();
690  }
691  static inline void ptop_scalar_shape(std::vector<ElOp<3>> &operations,
692  std::vector< std::vector<double> > shape_values, uint scalar_shape_op_idx) {
693  auto &op = operations[scalar_shape_op_idx];
694  uint n_dofs = shape_values.size();
695  uint n_points = shape_values[0].size();
696  uint n_elem = op.result_matrix()(0).data_size() / n_points;
697 
698  auto &shape_matrix = op.result_matrix();
699  ArenaVec<double> ref_vec(n_points * n_dofs, op.result_matrix()(0).arena());
700  ArenaVec<double> elem_vec(n_elem, op.result_matrix()(0).arena());
701  for (uint i=0; i<n_elem; ++i) {
702  elem_vec(i) = 1.0;
703  }
704  ArenaOVec<double> ref_ovec(ref_vec);
705  ArenaOVec<double> elem_ovec(elem_vec);
706  for (uint i_dof=0; i_dof<n_dofs; ++i_dof)
707  for (uint i_p=0; i_p<n_points; ++i_p)
708  ref_vec(i_dof * n_points + i_p) = shape_values[i_dof][i_p];
709  ArenaOVec<double> shape_ovec = elem_ovec * ref_ovec;
710  shape_matrix(0) = shape_ovec.get_vec();
711  }
712  static inline void ptop_vector_shape(std::vector<ElOp<3>> &operations, std::vector< std::vector<arma::vec3> > shape_values, uint vector_shape_op_idx) {
713  auto &op = operations[vector_shape_op_idx];
714  uint n_dofs = shape_values.size();
715  uint n_points = shape_values[0].size();
716  uint n_elem = op.result_matrix()(0).data_size() / n_points;
717 
718  auto &shape_matrix = op.result_matrix(); // result: shape 3x1
719  Eigen::Matrix<ArenaVec<double>, 3, 1> ref_shape_vec;
720  Eigen::Matrix<ArenaOVec<double>, 3, 1> ref_shape_ovec;
721  ArenaVec<double> elem_vec(n_elem, op.result_matrix()(0).arena());
722  for (uint i=0; i<n_elem; ++i) {
723  elem_vec(i) = 1.0;
724  }
725  for (uint c=0; c<3; ++c) {
726  ref_shape_vec(c) = ArenaVec<double>(n_points * n_dofs, op.result_matrix()(0).arena());
727  ref_shape_ovec(c) = ArenaOVec(ref_shape_vec(c));
728  }
729  ArenaOVec<double> elem_ovec(elem_vec);
730  for (uint i_dof=0; i_dof<n_dofs; ++i_dof)
731  for (uint i_p=0; i_p<n_points; ++i_p)
732  for (uint c=0; c<3; ++c)
733  ref_shape_vec(c,0)(i_dof * n_points + i_p) = shape_values[i_dof][i_p](c);
734  Eigen::Matrix<ArenaOVec<double>, 1, 3> shape_omatrix = elem_ovec * ref_shape_ovec.transpose();
735  for (uint c=0; c<3; ++c)
736  shape_matrix(c) = shape_omatrix(0,c).get_vec();
737  }
739  FMT_UNUSED std::vector< std::vector<arma::vec3> > shape_values, FMT_UNUSED uint vector_shape_op_idx) {
740 // auto &op = operations[vector_sym_grad_op_idx];
741 // auto &grad_vector_value = operations[ op.input_ops()[0] ].result_matrix();
742 // auto &sym_grad_value = op.result_matrix();
743 // sym_grad_value = 0.5 * (grad_vector_value.transpose() + grad_vector_value);
744  }
745  static inline void ptop_vector_piola_shape(FMT_UNUSED std::vector<ElOp<3>> &operations,
746  FMT_UNUSED std::vector< std::vector<arma::vec3> > shape_values, FMT_UNUSED uint vector_shape_op_idx) {
747 // auto &op = operations[vector_sym_grad_op_idx];
748 // auto &grad_vector_value = operations[ op.input_ops()[0] ].result_matrix();
749 // auto &sym_grad_value = op.result_matrix();
750 // sym_grad_value = 0.5 * (grad_vector_value.transpose() + grad_vector_value);
751  }
752  template<unsigned int dim>
753  static inline void ptop_scalar_shape_grads(std::vector<ElOp<3>> &operations,
754  std::vector< std::vector<arma::mat> > ref_shape_grads, uint scalar_shape_grads_op_idx) {
755  auto &op = operations[scalar_shape_grads_op_idx];
756  auto &inv_jac_vec = operations[ op.input_ops()[0] ].result_matrix();
757  uint n_points = ref_shape_grads.size();
758  uint n_dofs = ref_shape_grads[0].size();
759 
760  Eigen::Matrix<ArenaVec<double>, dim, 1> ref_grads_vec;
761  Eigen::Matrix<ArenaOVec<double>, dim, 1> ref_grads_ovec;
762  for (uint i=0; i<ref_grads_vec.rows(); ++i) {
763  ref_grads_vec(i) = ArenaVec<double>(n_points * n_dofs, op.result_matrix()(0).arena());
764  for (uint i_dof=0; i_dof<n_dofs; ++i_dof)
765  for (uint i_p=0; i_p<n_points; ++i_p)
766  ref_grads_vec(i)(i_dof * n_points + i_p) = ref_shape_grads[i_p][i_dof](i);
767  ref_grads_ovec(i) = ArenaOVec(ref_grads_vec(i));
768  }
769 
770  Eigen::Matrix<ArenaOVec<double>, dim, 3> inv_jac_ovec;
771  for (uint i=0; i<dim*3; ++i) {
772  inv_jac_ovec(i) = ArenaOVec(inv_jac_vec(i));
773  }
774 
775  auto &result_vec = op.result_matrix();
776  Eigen::Matrix<ArenaOVec<double>, 3, 1> result_ovec = inv_jac_ovec.transpose() * ref_grads_ovec;
777  for (uint i=0; i<3; ++i) {
778  result_vec(i) = result_ovec(i).get_vec();
779  }
780  }
781  template<unsigned int dim>
782  static inline void ptop_vector_shape_grads(std::vector<ElOp<3>> &operations,
783  std::vector< std::vector<arma::mat> > ref_shape_grads, uint vector_shape_grads_op_idx) {
784  auto &op = operations[vector_shape_grads_op_idx];
785  auto &inv_jac_vec = operations[ op.input_ops()[0] ].result_matrix();
786  uint n_points = ref_shape_grads.size();
787  uint n_dofs = ref_shape_grads[0].size();
788 
789  Eigen::Matrix<ArenaVec<double>, dim, 3> ref_grads_vec;
790  Eigen::Matrix<ArenaOVec<double>, dim, 3> ref_grads_ovec;
791  for (uint i=0; i<ref_grads_vec.rows()*ref_grads_vec.cols(); ++i) {
792  ref_grads_vec(i) = ArenaVec<double>(n_points * n_dofs, op.result_matrix()(0).arena());
793  for (uint i_dof=0; i_dof<n_dofs; ++i_dof)
794  for (uint i_p=0; i_p<n_points; ++i_p)
795  ref_grads_vec(i)(i_dof * n_points + i_p) = ref_shape_grads[i_p][i_dof](i);
796  ref_grads_ovec(i) = ArenaOVec(ref_grads_vec(i));
797  }
798 
799  Eigen::Matrix<ArenaOVec<double>, dim, 3> inv_jac_ovec;
800  for (uint i=0; i<dim*3; ++i) {
801  inv_jac_ovec(i) = ArenaOVec(inv_jac_vec(i));
802  }
803 
804  auto &result_vec = op.result_matrix();
805  Eigen::Matrix<ArenaOVec<double>, 3, 3> result_ovec = inv_jac_ovec.transpose() * ref_grads_ovec;
806  for (uint i=0; i<9; ++i) {
807  result_vec(i) = result_ovec(i).get_vec();
808  }
809  }
810  template<unsigned int dim>
812  FMT_UNUSED std::vector< std::vector<arma::mat> > ref_shape_grads, FMT_UNUSED uint vector_shape_grads_op_idx) {
813 // auto &op = operations[scalar_shape_grads_op_idx];
814 // auto &inv_jac_vec = operations[ op.input_ops()[0] ].result_matrix();
815 // uint n_points = ref_shape_grads.size();
816 // uint n_dofs = ref_shape_grads[0].size();
817 //
818 // Eigen::Matrix<ArenaVec<double>, dim, 1> ref_grads_vec;
819 // Eigen::Matrix<ArenaOVec<double>, dim, 1> ref_grads_ovec;
820 // for (uint i=0; i<ref_grads_vec.rows(); ++i) {
821 // ref_grads_vec(i) = ArenaVec<double>(n_points * n_dofs, op.result_matrix()(0).arena());
822 // for (uint i_dof=0; i_dof<n_dofs; ++i_dof)
823 // for (uint i_p=0; i_p<n_points; ++i_p)
824 // ref_grads_vec(i)(i_dof * n_points + i_p) = ref_shape_grads[i_p][i_dof](i);
825 // ref_grads_ovec(i) = ArenaOVec(ref_grads_vec(i));
826 // }
827 //
828 // Eigen::Matrix<ArenaOVec<double>, dim, 3> inv_jac_ovec;
829 // for (uint i=0; i<dim*3; ++i) {
830 // inv_jac_ovec(i) = ArenaOVec(inv_jac_vec(i));
831 // }
832 //
833 // auto &result_vec = op.result_matrix();
834 // Eigen::Matrix<ArenaOVec<double>, 3, 1> result_ovec = inv_jac_ovec.transpose() * ref_grads_ovec;
835 // for (uint i=0; i<3; ++i) {
836 // result_vec(i) = result_ovec(i).get_vec();
837 // }
838  }
839  template<unsigned int dim>
841  FMT_UNUSED std::vector< std::vector<arma::mat> > ref_shape_grads, FMT_UNUSED uint vector_shape_grads_op_idx) {
842 // auto &op = operations[scalar_shape_grads_op_idx];
843 // auto &inv_jac_vec = operations[ op.input_ops()[0] ].result_matrix();
844 // uint n_points = ref_shape_grads.size();
845 // uint n_dofs = ref_shape_grads[0].size();
846 //
847 // Eigen::Matrix<ArenaVec<double>, dim, 1> ref_grads_vec;
848 // Eigen::Matrix<ArenaOVec<double>, dim, 1> ref_grads_ovec;
849 // for (uint i=0; i<ref_grads_vec.rows(); ++i) {
850 // ref_grads_vec(i) = ArenaVec<double>(n_points * n_dofs, op.result_matrix()(0).arena());
851 // for (uint i_dof=0; i_dof<n_dofs; ++i_dof)
852 // for (uint i_p=0; i_p<n_points; ++i_p)
853 // ref_grads_vec(i)(i_dof * n_points + i_p) = ref_shape_grads[i_p][i_dof](i);
854 // ref_grads_ovec(i) = ArenaOVec(ref_grads_vec(i));
855 // }
856 //
857 // Eigen::Matrix<ArenaOVec<double>, dim, 3> inv_jac_ovec;
858 // for (uint i=0; i<dim*3; ++i) {
859 // inv_jac_ovec(i) = ArenaOVec(inv_jac_vec(i));
860 // }
861 //
862 // auto &result_vec = op.result_matrix();
863 // Eigen::Matrix<ArenaOVec<double>, 3, 1> result_ovec = inv_jac_ovec.transpose() * ref_grads_ovec;
864 // for (uint i=0; i<3; ++i) {
865 // result_vec(i) = result_ovec(i).get_vec();
866 // }
867  }
868 };
869 
870 
871 
872 /// Defines reinit operations on side points.
873 struct side_reinit {
874  // element operations
875  template<unsigned int dim>
876  static inline void elop_el_jac(std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
877  // result matrix(spacedim, dim), input matrix(spacedim, dim+1)
878  auto &op = operations[FeSide::SideOps::opElJac];
879  auto &jac_value = op.result_matrix();
880  const auto &coords_value = operations[ op.input_ops()[0] ].result_matrix();
881  for (unsigned int i=0; i<3; i++)
882  for (unsigned int j=0; j<dim; j++)
883  jac_value(i,j) = coords_value(i,j+1) - coords_value(i,0);
884  }
885  template<unsigned int dim>
886  static inline void elop_el_inv_jac(std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
887  auto &op = operations[FeSide::SideOps::opElInvJac];
888  auto &inv_jac_value = op.result_matrix();
889  const auto &jac_value = operations[ op.input_ops()[0] ].result_matrix();
890  inv_jac_value = eigen_arena_tools::inverse<3, dim>(jac_value);
891  }
892  template<unsigned int dim>
893  static inline void elop_sd_jac(std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
894  // result matrix(spacedim, dim), input matrix(spacedim, dim+1)
895  auto &op = operations[FeSide::SideOps::opSideJac];
896  auto &jac_value = op.result_matrix();
897  const auto &coords_value = operations[ op.input_ops()[0] ].result_matrix();
898  for (unsigned int i=0; i<3; i++)
899  for (unsigned int j=0; j<dim-1; j++)
900  jac_value(i,j) = coords_value(i,j+1) - coords_value(i,0);
901  }
902  template<unsigned int dim>
903  static inline void elop_sd_jac_det(std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
904  // result double, input matrix(spacedim, dim)
905  auto &op = operations[FeSide::SideOps::opSideJacDet];
906  auto &jac_det_value = op.result_matrix();
907  const auto &jac_value = operations[ op.input_ops()[0] ].result_matrix();
908  jac_det_value(0,0) = eigen_arena_tools::determinant<3, dim-1>(jac_value).abs();
909  }
910 
911  // Point operations
912  static inline void ptop_coords(FMT_UNUSED std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
913  // Implement
914  }
915  static inline void ptop_weights(std::vector<ElOp<3>> &operations, PatchArena *arena, const std::vector<double> &point_weights) {
916  auto &op = operations[FeSide::SideOps::opWeights];
917  op.allocate_result(point_weights.size(), *arena);
918  auto &weights_value = op.result_matrix();
919  for (uint i=0; i<point_weights.size(); ++i)
920  weights_value(0,0)(i) = point_weights[i];
921  }
922  static inline void ptop_JxW(std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
923  auto &op = operations[FeSide::SideOps::opJxW];
924  auto &weights_value = operations[ op.input_ops()[0] ].result_matrix();
925  auto &jac_det_value = operations[ op.input_ops()[1] ].result_matrix();
926  ArenaOVec<double> weights_ovec( weights_value(0,0) );
927  ArenaOVec<double> jac_det_ovec( jac_det_value(0,0) );
928  ArenaOVec<double> jxw_ovec = jac_det_ovec * weights_ovec;
929  auto &jxw_value = op.result_matrix();
930  jxw_value(0,0) = jxw_ovec.get_vec();
931  }
932  template<unsigned int dim>
933  static inline void ptop_normal_vec(std::vector<ElOp<3>> &operations, IntTableArena &el_table) {
934  auto &op = operations[FeSide::SideOps::opNormalVec];
935  auto &normal_value = op.result_matrix();
936  auto &inv_jac_mat_value = operations[ op.input_ops()[0] ].result_matrix();
937  normal_value = inv_jac_mat_value.transpose() * RefElement<dim>::normal_vector_array( el_table(3) );
938 
939  ArenaVec<double> norm_vec( normal_value(0).data_size(), normal_value(0).arena() );
940  Eigen::VectorXd A(3);
941  for (uint i=0; i<normal_value(0).data_size(); ++i) {
942  A(0) = normal_value(0)(i);
943  A(1) = normal_value(1)(i);
944  A(2) = normal_value(2)(i);
945  norm_vec(i) = A.norm();
946  }
947 
948  size_t points_per_side = el_table(4).data_size() / el_table(3).data_size();
949  size_t n_points = el_table(3).data_size();
950  for (uint i=0; i<3; ++i) {
951  normal_value(i) = normal_value(i) / norm_vec;
952  ArenaVec<double> expand_vec( normal_value(i).data_size() * points_per_side, normal_value(i).arena() );
953  for (uint j=0; j<expand_vec.data_size(); ++j) {
954  expand_vec(j) = normal_value(i)(j % n_points);
955  }
956  normal_value(i) = expand_vec;
957  }
958  }
959  static inline void ptop_scalar_shape(std::vector<ElOp<3>> &operations, IntTableArena &el_table,
960  std::vector< std::vector< std::vector<double> > > shape_values, uint scalar_shape_op_idx) {
961  uint n_dofs = shape_values[0][0].size();
962  uint n_sides = el_table(3).data_size();
963  uint n_patch_points = el_table(4).data_size();
964 
965  auto &op = operations[scalar_shape_op_idx];
966  auto &scalar_shape_value = op.result_matrix();
967  scalar_shape_value(0) = ArenaVec<double>(n_dofs*n_patch_points, scalar_shape_value(0).arena());
968 
969  for (uint i_dof=0; i_dof<n_dofs; ++i_dof) {
970  uint dof_shift = i_dof * n_patch_points;
971  for (uint i_pt=0; i_pt<n_patch_points; ++i_pt)
972  scalar_shape_value(0)(i_pt + dof_shift) = shape_values[el_table(4)(i_pt)][i_pt / n_sides][i_dof];
973  }
974  }
975  static inline void ptop_vector_shape(std::vector<ElOp<3>> &operations, IntTableArena &el_table,
976  std::vector< std::vector< std::vector<arma::vec3> > > shape_values, uint vector_shape_op_idx) {
977  uint n_dofs = shape_values[0][0].size();
978  uint n_sides = el_table(3).data_size();
979  uint n_patch_points = el_table(4).data_size();
980 
981  auto &op = operations[vector_shape_op_idx];
982  auto &vector_shape_value = op.result_matrix();
983  for (uint c=0; c<3; c++)
984  vector_shape_value(c) = ArenaVec<double>(n_dofs*n_patch_points, vector_shape_value(c).arena());
985 
986  for (uint i_dof=0; i_dof<n_dofs; ++i_dof) {
987  uint dof_shift = i_dof * n_patch_points;
988  for (uint i_pt=0; i_pt<n_patch_points; ++i_pt)
989  for (uint c=0; c<3; c++)
990  vector_shape_value(c)(i_pt + dof_shift) = shape_values[el_table(4)(i_pt)][i_pt / n_sides][i_dof](c);
991  }
992  }
994  FMT_UNUSED std::vector< std::vector< std::vector<arma::vec3> > > shape_values, FMT_UNUSED uint vector_shape_op_idx) {}
995  static inline void ptop_vector_piola_shape(FMT_UNUSED std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table,
996  FMT_UNUSED std::vector< std::vector< std::vector<arma::vec3> > > shape_values, FMT_UNUSED uint vector_shape_op_idx) {}
997  template<unsigned int dim>
998  static inline void ptop_scalar_shape_grads(std::vector<ElOp<3>> &operations, IntTableArena &el_table,
999  std::vector< std::vector< std::vector<arma::mat> > > ref_shape_grads, uint scalar_shape_grads_op_idx) {
1000  uint n_points = ref_shape_grads[0].size();
1001  uint n_dofs = ref_shape_grads[0][0].size();
1002  uint n_sides = el_table(3).data_size();
1003  uint n_patch_points = el_table(4).data_size();
1004 
1005  // Result vector
1006  auto &op = operations[scalar_shape_grads_op_idx];
1007  auto &grad_scalar_shape_value = op.result_matrix();
1008 
1009  // Expands inverse jacobian to inv_jac_expd_value
1010  auto &inv_jac_value = operations[ op.input_ops()[0] ].result_matrix();
1011  Eigen::Matrix<ArenaVec<double>, dim, 3> inv_jac_expd_value;
1012  for (uint i=0; i<dim*3; ++i) {
1013  inv_jac_expd_value(i) = ArenaVec<double>( n_dofs*n_patch_points, inv_jac_value(i).arena() );
1014  for (uint j=0; j<n_dofs*n_patch_points; ++j)
1015  inv_jac_expd_value(i)(j) = inv_jac_value(i)(j%n_sides);
1016  }
1017 
1018  // Fill ref shape gradients by q_point. DOF and side_idx
1019  Eigen::Matrix<ArenaVec<double>, dim, 1> ref_shape_grads_expd;
1020  for (uint i=0; i<dim; ++i) {
1021  ref_shape_grads_expd(i) = ArenaVec<double>( n_dofs*n_patch_points, inv_jac_value(0).arena() );
1022  }
1023  for (uint i_dof=0; i_dof<n_dofs; ++i_dof) {
1024  for (uint i_pt=0; i_pt<n_points; ++i_pt) {
1025  uint i_begin = (i_dof * n_points + i_pt) * n_sides;
1026  for (uint i_sd=0; i_sd<n_sides; ++i_sd) {
1027  for (uint i_c=0; i_c<dim; ++i_c) {
1028  ref_shape_grads_expd(i_c)(i_begin + i_sd) = ref_shape_grads[el_table(3)(i_sd)][i_pt][i_dof][i_c];
1029  }
1030  }
1031  }
1032  }
1033 
1034  // computes operation result
1035  grad_scalar_shape_value = inv_jac_expd_value.transpose() * ref_shape_grads_expd;
1036  }
1037  template<unsigned int dim>
1038  static inline void ptop_vector_shape_grads(std::vector<ElOp<3>> &operations, IntTableArena &el_table,
1039  std::vector< std::vector< std::vector<arma::mat> > > ref_shape_grads, uint vector_shape_grads_op_idx) {
1040  uint n_points = ref_shape_grads[0].size();
1041  uint n_dofs = ref_shape_grads[0][0].size();
1042  uint n_sides = el_table(3).data_size();
1043  uint n_patch_points = el_table(4).data_size();
1044 
1045  // Result vector
1046  auto &op = operations[vector_shape_grads_op_idx];
1047  auto &grad_scalar_shape_value = op.result_matrix();
1048 
1049  // Expands inverse jacobian to inv_jac_expd_value
1050  auto &inv_jac_value = operations[ op.input_ops()[0] ].result_matrix();
1051  Eigen::Matrix<ArenaVec<double>, dim, 3> inv_jac_expd_value;
1052  for (uint i=0; i<dim*3; ++i) {
1053  inv_jac_expd_value(i) = ArenaVec<double>( n_dofs*n_patch_points, inv_jac_value(i).arena() );
1054  for (uint j=0; j<n_dofs*n_patch_points; ++j)
1055  inv_jac_expd_value(i)(j) = inv_jac_value(i)(j%n_sides);
1056  }
1057 
1058  // Fill ref shape gradients by q_point. DOF and side_idx
1059  Eigen::Matrix<ArenaVec<double>, dim, 3> ref_shape_grads_expd;
1060  for (uint i=0; i<3*dim; ++i) {
1061  ref_shape_grads_expd(i) = ArenaVec<double>( n_dofs*n_patch_points, inv_jac_value(0).arena() );
1062  }
1063  for (uint i_dof=0; i_dof<n_dofs; ++i_dof) {
1064  for (uint i_pt=0; i_pt<n_points; ++i_pt) {
1065  uint i_begin = (i_dof * n_points + i_pt) * n_sides;
1066  for (uint i_sd=0; i_sd<n_sides; ++i_sd) {
1067  for (uint i_c=0; i_c<3*dim; ++i_c) {
1068  ref_shape_grads_expd(i_c)(i_begin + i_sd) = ref_shape_grads[el_table(3)(i_sd)][i_pt][i_dof][i_c];
1069  }
1070  }
1071  }
1072  }
1073 
1074  // computes operation result
1075  grad_scalar_shape_value = inv_jac_expd_value.transpose() * ref_shape_grads_expd;
1076 
1077 // ---/ OLD CODE
1078 // Eigen::Matrix<ArenaOVec<double>, dim, 3> inv_jac_ovec;
1079 // for (uint i=0; i<dim*3; ++i) {
1080 // inv_jac_ovec(i) = ArenaOVec(inv_jac_vec(i));
1081 // }
1082 //
1083 // auto &result_vec = op.result_matrix();
1084 // Eigen::Matrix<ArenaOVec<double>, 3, 3> result_ovec = inv_jac_ovec.transpose() * ref_grads_ovec;
1085 // for (uint i=0; i<9; ++i) {
1086 // result_vec(i) = result_ovec(i).get_vec();
1087 // }
1088  }
1089  template<unsigned int dim>
1091  FMT_UNUSED std::vector< std::vector< std::vector<arma::mat> > > ref_shape_grads, FMT_UNUSED uint vector_shape_grads_op_idx) {}
1092  template<unsigned int dim>
1094  FMT_UNUSED std::vector< std::vector< std::vector<arma::mat> > > ref_shape_grads, FMT_UNUSED uint vector_shape_grads_op_idx) {}
1095 };
1096 
1097 
1098 // template specialization
1099 template<>
1100 inline void side_reinit::elop_sd_jac<1>(FMT_UNUSED std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
1101 }
1102 
1103 template<>
1104 inline void side_reinit::elop_sd_jac_det<1>(std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
1105  auto &op = operations[FeSide::SideOps::opSideJacDet];
1106  auto &result_vec = op.result_matrix();
1107  for (uint i=0;i<result_vec(0,0).data_size(); ++i) {
1108  result_vec(0,0)(i) = 1.0;
1109  }
1110 }
1111 
1112 
1113 
1114 namespace FeBulk {
1115 
1116  /// Bulk data specialization, order of item in operations_ vector corresponds to the BulkOps enum
1117  template<unsigned int spacedim = 3>
1118  class PatchPointValues : public ::PatchPointValues<spacedim> {
1119  public:
1120  /// Constructor
1121  PatchPointValues(uint dim, uint quad_order, AssemblyArena &asm_arena)
1122  : ::PatchPointValues<spacedim>(dim, asm_arena) {
1123  this->quad_ = new QGauss(dim, 2*quad_order);
1124  this->int_sizes_ = {pointOp, pointOp, pointOp};
1125  switch (dim) {
1126  case 1:
1127  init<1>();
1128  break;
1129  case 2:
1130  init<2>();
1131  break;
1132  case 3:
1133  init<3>();
1134  break;
1135  }
1136  }
1137 
1138  private:
1139  /// Initialize operations vector
1140  template<unsigned int dim>
1141  void init() {
1142  // First step: adds element values operations
1143  /*auto &el_coords =*/ this->make_new_op( {spacedim, this->dim_+1}, &common_reinit::op_base, {}, OpSizeType::elemOp );
1144 
1145  /*auto &el_jac =*/ this->make_new_op( {spacedim, this->dim_}, &bulk_reinit::elop_jac<dim>, {BulkOps::opElCoords}, OpSizeType::elemOp );
1146 
1147  /*auto &el_inv_jac =*/ this->make_new_op( {this->dim_, spacedim}, &bulk_reinit::elop_inv_jac<dim>, {BulkOps::opJac}, OpSizeType::elemOp );
1148 
1149  /*auto &el_jac_det =*/ this->make_new_op( {1}, &bulk_reinit::elop_jac_det<dim>, {BulkOps::opJac}, OpSizeType::elemOp );
1150 
1151  // Second step: adds point values operations
1152  /*auto &pt_coords =*/ this->make_new_op( {spacedim}, &bulk_reinit::ptop_coords, {} );
1153 
1154  // use lambda reinit function
1155  const std::vector<double> &point_weights_vec = this->quad_->get_weights();
1156  auto lambda_weights = [this, point_weights_vec](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
1157  bulk_reinit::ptop_weights(operations, this->patch_arena_, point_weights_vec);
1158  };
1159  /*auto &weights =*/ this->make_fixed_op( {1}, lambda_weights );
1160 
1161  /*auto &JxW =*/ this->make_new_op( {1}, &bulk_reinit::ptop_JxW, {BulkOps::opWeights, BulkOps::opJacDet} );
1162  }
1163  };
1164 
1165 } // closing namespace FeBulk
1166 
1167 
1168 
1169 namespace FeSide {
1170 
1171 /// Bulk Side specialization, order of item in operations_ vector corresponds to the SideOps enum
1172  template<unsigned int spacedim = 3>
1173  class PatchPointValues : public ::PatchPointValues<spacedim> {
1174  public:
1175  /// Constructor
1176  PatchPointValues(uint dim, uint quad_order, AssemblyArena &asm_arena)
1177  : ::PatchPointValues<spacedim>(dim, asm_arena) {
1178  this->quad_ = new QGauss(dim-1, 2*quad_order);
1179  this->int_sizes_ = {pointOp, pointOp, pointOp, elemOp, pointOp};
1180  switch (dim) {
1181  case 1:
1182  init<1>();
1183  break;
1184  case 2:
1185  init<2>();
1186  break;
1187  case 3:
1188  init<3>();
1189  break;
1190  }
1191  }
1192 
1193  private:
1194  /// Initialize operations vector
1195  template<unsigned int dim>
1196  void init() {
1197  // First step: adds element values operations
1198  /*auto &el_coords =*/ this->make_new_op( {spacedim, this->dim_+1}, &common_reinit::op_base, {}, OpSizeType::elemOp );
1199 
1200  /*auto &el_jac =*/ this->make_new_op( {spacedim, this->dim_}, &side_reinit::elop_el_jac<dim>, {SideOps::opElCoords}, OpSizeType::elemOp );
1201 
1202  /*auto &el_inv_jac =*/ this->make_new_op( {this->dim_, spacedim}, &side_reinit::elop_el_inv_jac<dim>, {SideOps::opElJac}, OpSizeType::elemOp );
1203 
1204  // Second step: adds side values operations
1205  /*auto &sd_coords =*/ this->make_new_op( {spacedim, this->dim_}, &common_reinit::op_base, {}, OpSizeType::elemOp );
1206 
1207  /*auto &sd_jac =*/ this->make_new_op( {spacedim, this->dim_-1}, &side_reinit::elop_sd_jac<dim>, {SideOps::opSideCoords}, OpSizeType::elemOp );
1208 
1209  /*auto &sd_jac_det =*/ this->make_new_op( {1}, &side_reinit::elop_sd_jac_det<dim>, {SideOps::opSideJac}, OpSizeType::elemOp );
1210 
1211  // Third step: adds point values operations
1212  /*auto &coords =*/ this->make_new_op( {spacedim}, &side_reinit::ptop_coords, {} );
1213 
1214  // use lambda reinit function
1215  const std::vector<double> &point_weights_vec = this->quad_->get_weights();
1216  auto lambda_weights = [this, point_weights_vec](std::vector<ElOp<3>> &operations, FMT_UNUSED IntTableArena &el_table) {
1217  side_reinit::ptop_weights(operations, this->patch_arena_, point_weights_vec);
1218  };
1219  /*auto &weights =*/ this->make_fixed_op( {1}, lambda_weights );
1220 
1222 
1223  /*auto &normal_vec =*/ this->make_new_op( {spacedim}, &side_reinit::ptop_normal_vec<dim>, {SideOps::opElInvJac} );
1224  }
1225  };
1226 
1227 } // closing namespace FeSide
1228 
1229 
1230 
1231 #endif /* PATCH_POINT_VALUES_HH_ */
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:333
ArenaVec< T > get_vec() const
Convert ArenaOVec to ArenaVec and its.
Definition: arena_vec.hh:328
size_t data_size() const
Getter for data_size_.
Definition: arena_vec.hh:104
Class represents element or FE operations.
ElOp(uint dim, std::initializer_list< uint > shape, ReinitFunction reinit_f, OpSizeType size_type, std::vector< uint > input_ops={}, uint n_dofs=1)
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)
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)
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.
void reinit_function(std::vector< ElOp< spacedim >> &operations, IntTableArena &int_table)
Call reinit function on element table if function is defined.
uint n_dofs_
Number of DOFs of FE operations (or 1 in case of element operations)
uint dim() const
Getter for dimension.
uint n_comp() const
OpSizeType size_type() const
Getter for size_type_.
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)
AssemblyArena & asm_arena_
Reference to global assembly arena of PatchFeValues.
std::vector< OpSizeType > row_sizes_
hold sizes of rows by type of operation
std::vector< uint > points_map_
Map of point patch indices to ElOp::result_ and int_table_ tables.
std::vector< unsigned int > fe_ops_indices_
Indices of FE operations in operations_ vector.
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.
unsigned int get_fe_op(FEOps fe_op) const
Return index of FE operation in operations_ vector.
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)
void set_fe_op(FEOps fe_op, unsigned int op_vec_idx)
Set index of FE operation in operations_ vector.
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)
std::vector< uint > elements_map_
Map of element patch indices to ElOp::result_ and int_table_ tables.
uint dim() const
Getter for dim_.
Tensor tensor_value(uint op_idx, uint point_idx, uint i_dof=0) const
Scalar scalar_value(uint op_idx, uint point_idx, uint i_dof=0) const
void init_finalize(PatchArena *patch_arena)
uint register_side(arma::mat elm_coords, arma::mat side_coords, uint side_idx)
uint register_side_point(uint elem_table_row, uint value_patch_idx, uint elem_idx, uint side_idx, uint i_point_on_side)
IntTableArena int_table_
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.
ElOp< spacedim > & make_new_op(std::initializer_list< uint > shape, ReinitFunction reinit_f, std::vector< uint > input_ops_vec, OpSizeType size_type=pointOp)
std::vector< OpSizeType > int_sizes_
Set size and type of rows of int_table_, value is set implicitly in constructor of descendants.
Vector vector_value(uint op_idx, uint point_idx, uint i_dof=0) 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::Matrix< ArenaVec< double >, dim, 1 > normal_vector_array(ArenaVec< uint > loc_side_idx_vec)
Definition: ref_element.cc:280
Store finite element data on the actual patch such as shape function values, gradients,...
Eigen::Vector< ArenaVec< uint >, Eigen::Dynamic > IntTableArena
Definition: eigen_tools.hh:39
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
@ 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.
@ opNormalVec
normal vector of quadrature point
@ opSideCoords
operations evaluated on sides
@ opJxW
JxW value of quadrature point.
@ opCoords
operations evaluated on quadrature points
@ opElCoords
operations evaluated on elements
@ opElJac
Jacobian of element.
@ opWeights
weight of quadrature point
@ opSideJac
Jacobian of element.
ArenaVec< double > determinant(const Eigen::Matrix< ArenaVec< double >, m, n > &A)
Calculates determinant of a rectangular matrix.
@ OpNItems
Holds number of valid FE operations and value of invalid FE operation.
@ opVectorShape
Vector shape operation.
@ opVectorDivergence
Vector divergence.
@ opVectorSymGrad
Vector symmetric gradient.
@ opGradScalarShape
Scalar shape gradient.
@ opGradVectorShape
Vector shape gradient.
@ opScalarShape
Scalar shape operation.
double Scalar
std::function< void(std::vector< ElOp< 3 > > &, IntTableArena &)> ReinitFunction
Type for conciseness.
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 ptop_vector_contravariant_shape_grads(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED std::vector< std::vector< arma::mat > > ref_shape_grads, FMT_UNUSED uint vector_shape_grads_op_idx)
static void elop_inv_jac(std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
static void ptop_vector_shape(std::vector< ElOp< 3 >> &operations, std::vector< std::vector< arma::vec3 > > shape_values, uint vector_shape_op_idx)
static void ptop_vector_contravariant_shape(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED std::vector< std::vector< arma::vec3 > > shape_values, FMT_UNUSED uint vector_shape_op_idx)
static void ptop_scalar_shape_grads(std::vector< ElOp< 3 >> &operations, std::vector< std::vector< arma::mat > > ref_shape_grads, uint scalar_shape_grads_op_idx)
static void ptop_vector_piola_shape(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED std::vector< std::vector< arma::vec3 > > shape_values, FMT_UNUSED uint vector_shape_op_idx)
static void elop_jac(std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
static void ptop_vector_piola_shape_grads(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED std::vector< std::vector< arma::mat > > ref_shape_grads, FMT_UNUSED uint vector_shape_grads_op_idx)
static void elop_jac_det(std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
static void ptop_JxW(std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
static void ptop_weights(std::vector< ElOp< 3 >> &operations, PatchArena *arena, const std::vector< double > &point_weights)
static void ptop_scalar_shape(std::vector< ElOp< 3 >> &operations, std::vector< std::vector< double > > shape_values, uint scalar_shape_op_idx)
static void ptop_vector_shape_grads(std::vector< ElOp< 3 >> &operations, std::vector< std::vector< arma::mat > > ref_shape_grads, uint vector_shape_grads_op_idx)
static void ptop_coords(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
Defines common functionality of reinit operations.
static void ptop_vector_divergence(std::vector< ElOp< 3 >> &operations, uint vector_divergence_op_idx)
Common reinit function of vector divergence on bulk and side points.
static void ptop_vector_sym_grad(std::vector< ElOp< 3 >> &operations, uint vector_sym_grad_op_idx)
Common reinit function of vector symmetric gradient on bulk and side points.
static void op_base(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
Defines reinit operations on side points.
static void ptop_vector_shape(std::vector< ElOp< 3 >> &operations, IntTableArena &el_table, std::vector< std::vector< std::vector< arma::vec3 > > > shape_values, uint vector_shape_op_idx)
static void ptop_vector_shape_grads(std::vector< ElOp< 3 >> &operations, IntTableArena &el_table, std::vector< std::vector< std::vector< arma::mat > > > ref_shape_grads, uint vector_shape_grads_op_idx)
static void elop_el_jac(std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
static void elop_sd_jac_det(std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
static void elop_sd_jac(std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
static void ptop_normal_vec(std::vector< ElOp< 3 >> &operations, IntTableArena &el_table)
static void ptop_weights(std::vector< ElOp< 3 >> &operations, PatchArena *arena, const std::vector< double > &point_weights)
static void ptop_scalar_shape_grads(std::vector< ElOp< 3 >> &operations, IntTableArena &el_table, std::vector< std::vector< std::vector< arma::mat > > > ref_shape_grads, uint scalar_shape_grads_op_idx)
static void ptop_scalar_shape(std::vector< ElOp< 3 >> &operations, IntTableArena &el_table, std::vector< std::vector< std::vector< double > > > shape_values, uint scalar_shape_op_idx)
static void ptop_vector_contravariant_shape_grads(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table, FMT_UNUSED std::vector< std::vector< std::vector< arma::mat > > > ref_shape_grads, FMT_UNUSED uint vector_shape_grads_op_idx)
static void ptop_vector_contravariant_shape(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table, FMT_UNUSED std::vector< std::vector< std::vector< arma::vec3 > > > shape_values, FMT_UNUSED uint vector_shape_op_idx)
static void elop_el_inv_jac(std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
static void ptop_vector_piola_shape(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table, FMT_UNUSED std::vector< std::vector< std::vector< arma::vec3 > > > shape_values, FMT_UNUSED uint vector_shape_op_idx)
static void ptop_JxW(std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
static void ptop_coords(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
static void ptop_vector_piola_shape_grads(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table, FMT_UNUSED std::vector< std::vector< std::vector< arma::mat > > > ref_shape_grads, FMT_UNUSED uint vector_shape_grads_op_idx)