Flow123d  DF_patch_fe_data_tables-3c41206
patch_point_values.hh
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file patch_point_values.hh
15  * @brief Store finite element data on the actual patch
16  * such as shape function values, gradients, Jacobian
17  * of the mapping from the reference cell etc.
18  * @author David Flanderka
19  */
20 
21 #ifndef PATCH_POINT_VALUES_HH_
22 #define PATCH_POINT_VALUES_HH_
23 
24 #include <Eigen/Dense>
25 
26 #include "fem/eigen_tools.hh"
27 #include "fem/dh_cell_accessor.hh"
28 #include "fem/element_values.hh"
30 
31 
32 template<unsigned int spacedim> class PatchFEValues;
33 template<unsigned int spacedim> class ElOp;
34 template<unsigned int dim> class BulkValues;
35 template<unsigned int dim> class SideValues;
36 using Scalar = double;
39 
40 
41 
42 
43 /// Type for conciseness
44 using ReinitFunction = std::function<void(std::vector<ElOp<3>> &, TableDbl &, TableInt &)>;
45 
46 
47 namespace FeBulk {
48  /// enum of bulk operation
49  enum BulkOps
50  {
61  opJxW
62  };
63 }
64 
65 
66 namespace FeSide {
67  /// enum of side operation
68  enum SideOps
69  {
86  };
87 }
88 
89 
90 /**
91  * @brief Class for storing FE data of quadrature points.
92  *
93  * Store data of bulk or side quadrature points of one dimension.
94  */
95 template<unsigned int spacedim = 3>
97 {
98 public:
99  /**
100  * Constructor
101  *
102  * Set dimension
103  */
105  : dim_(dim), n_rows_(0), elements_map_(300, 0), points_map_(300, 0) {}
106 
107  /**
108  * Initialize object, set number of columns (quantities) in tables.
109  *
110  * Number of columns of int_vals_ table is passed by argument, number of columns
111  * of other tables is given by n_rows_ value.
112  */
113  void initialize(uint int_cols) {
114  this->reset();
115 
116  point_vals_.resize(n_rows_);
117  int_vals_.resize(int_cols);
118  }
119 
120  /// Reset number of rows (points and elements)
121  inline void reset() {
122  n_points_ = 0;
123  n_elems_ = 0;
124  }
125 
126  /// Getter of dim_
127  inline uint dim() const {
128  return dim_;
129  }
130 
131  /// Getter of n_rows_
132  inline uint n_rows() const {
133  return n_rows_;
134  }
135 
136  /// Getter of n_elems_
137  inline uint n_elems() const {
138  return n_elems_;
139  }
140 
141  /// Getter of n_points_
142  inline uint n_points() const {
143  return n_points_;
144  }
145 
146  /// Return quadrature
148  return quad_;
149  }
150 
151  /// Adds the number of rows equal to n_added, returns index of first of them
152  /// Temporary method allow acces to old structure PatchFEValues::DimPatchFEValues
153  inline uint add_rows(uint n_added) {
154  uint old_size = n_rows_;
155  n_rows_ += n_added;
156  return old_size;
157  }
158 
159  /// Resize data tables
163  }
164 
165  /// Register element, add to point_vals_ table
166  uint register_element(arma::mat coords, uint element_patch_idx) {
167  uint res_column = operations_[FeBulk::BulkOps::opElCoords].result_row();
168  for (uint i_col=0; i_col<coords.n_cols; ++i_col)
169  for (uint i_row=0; i_row<coords.n_rows; ++i_row) {
170  point_vals_(res_column)(n_elems_) = coords(i_row, i_col);
171  ++res_column;
172  }
173 
174  elements_map_[element_patch_idx] = n_elems_;
175  return n_elems_++;
176  }
177 
178  /// Register side, add to point_vals_ table
179  uint register_side(arma::mat elm_coords, arma::mat side_coords) {
180  uint res_column = operations_[FeSide::SideOps::opElCoords].result_row();
181  for (uint i_col=0; i_col<elm_coords.n_cols; ++i_col)
182  for (uint i_row=0; i_row<elm_coords.n_rows; ++i_row) {
183  point_vals_(res_column)(n_elems_) = elm_coords(i_row, i_col);
184  ++res_column;
185  }
186 
187  res_column = operations_[FeSide::SideOps::opSdCoords].result_row();
188  for (uint i_col=0; i_col<side_coords.n_cols; ++i_col)
189  for (uint i_row=0; i_row<side_coords.n_rows; ++i_row) {
190  point_vals_(res_column)(n_elems_) = side_coords(i_row, i_col);
191  ++res_column;
192  }
193 
194  return n_elems_++;
195  }
196 
197  /// Register bulk point, add to int_vals_ table
198  uint register_bulk_point(uint elem_table_row, uint value_patch_idx, uint elem_idx) {
199  int_vals_(0)(n_points_) = value_patch_idx;
200  int_vals_(1)(n_points_) = elem_table_row;
201  int_vals_(2)(n_points_) = elem_idx;
202 
203  points_map_[value_patch_idx] = n_points_;
204  return n_points_++;
205  }
206 
207  /// Register side point, add to int_vals_ table
208  uint register_side_point(uint elem_table_row, uint value_patch_idx, uint elem_idx, uint side_idx) {
209  int_vals_(0)(n_points_) = value_patch_idx;
210  int_vals_(1)(n_points_) = elem_table_row;
211  int_vals_(2)(n_points_) = elem_idx;
212  int_vals_(3)(n_points_) = side_idx;
213 
214  points_map_[value_patch_idx] = n_points_;
215  return n_points_++;
216  }
217 
218  /// Add accessor to operations_ vector
219  ElOp<spacedim> &make_new_op(std::initializer_list<uint> shape, ReinitFunction reinit_f, std::vector<uint> input_ops_vec) {
220  ElOp<spacedim> op_accessor(this->dim_, shape, this->n_rows_, reinit_f, input_ops_vec);
221  this->n_rows_ += op_accessor.n_comp();
222  operations_.push_back(op_accessor);
223  return operations_[operations_.size()-1];
224  }
225 
226  /// Add accessor to operations_ vector
227  ElOp<spacedim> &make_expansion(ElOp<spacedim> &el_op, std::initializer_list<uint> shape, ReinitFunction reinit_f) {
228  ElOp<spacedim> op_accessor(this->dim_, shape, el_op.result_row(), reinit_f);
229  // shape passed from el_op throws:
230  // C++ exception with description "std::bad_alloc" thrown in the test body.
231  operations_.push_back(op_accessor);
232  return operations_[operations_.size()-1];
233  }
234 
235  /// Add accessor to operations_ vector
236  ElOp<spacedim> &make_fe_op(std::initializer_list<uint> shape, ReinitFunction reinit_f, std::vector<uint> input_ops_vec, uint n_dofs) {
237  ElOp<spacedim> op_accessor(this->dim_, shape, this->n_rows_, reinit_f, input_ops_vec, n_dofs);
238  this->n_rows_ += op_accessor.n_comp() * n_dofs;
239  operations_.push_back(op_accessor);
240  return operations_[operations_.size()-1];
241  }
242 
243 
244  /// Reinit data.
245  void reinit_patch() {
246  if (n_elems_ == 0) return; // skip if tables are empty
247  // precompute data on point_vals_ table
248  for (uint i=0; i<operations_.size(); ++i)
249  operations_[i].reinit_function(operations_, point_vals_, int_vals_);
250 
251 // // copy data from el_vals_ to point_vals_
252 // std::vector<uint> copied; // list of columns that will be copied
253 // for (uint i_op=0; i_op<operations_.size(); ++i_op)
254 // if (operations_[i_op].copy_vals()) {
255 // uint res_col = operations_[i_op].result_row();
256 // for (uint i_col=res_col; i_col<res_col+operations_[i_op].n_comp(); ++i_col)
257 // copied.push_back(i_col);
258 // }
259 // for (uint i_pt=0; i_pt<n_points_; ++i_pt) {
260 // uint el_table_idx = int_vals_(1)(i_pt);
261 // for (uint i_q=0; i_q<copied.size(); ++i_q)
262 // point_vals_(copied[i_q])(i_pt) = el_vals_(copied[i_q])(el_table_idx);
263 // }
264 //
265 // // precompute data on point_vals_ table
266 // for (uint i=0; i<operations_.size(); ++i)
267 // operations_[i].reinit_points(operations_, point_vals_);
268  }
269 
270  inline Scalar scalar_val(uint result_row, uint point_idx) const {
271  return point_vals_(result_row)(elements_map_[point_idx]);
272  }
273 
274  inline Vector vector_val(uint result_row, uint point_idx) const {
275  Vector val;
276  for (uint i=0; i<3; ++i)
277  val(i) = point_vals_(result_row+i)(elements_map_[point_idx]);
278  return val;
279  }
280 
281  inline Tensor tensor_val(uint result_row, uint point_idx) const {
282  Tensor val;
283  for (uint i=0; i<3; ++i)
284  for (uint j=0; j<3; ++j)
285  val(i,j) = point_vals_(result_row+3*i+j)(elements_map_[point_idx]);
286  return val;
287  }
288 
289  /// Temporary development method
290  void print(bool points, bool ints) const {
291  std::cout << "** Dimension: " << dim_ << std::endl;
292  if (points) {
293  std::cout << "Point vals: " << point_vals_.rows() << " - " << point_vals_.cols() << std::endl;
294  for (uint i_row=0; i_row<n_points_; ++i_row) {
295  for (uint i_col=0; i_col<n_rows_; ++i_col)
296  std::cout << point_vals_(i_col)(i_row) << " ";
297  std::cout << std::endl;
298  }
299  std::cout << std::endl;
300  }
301  if (ints) {
302  std::cout << "Int vals: " << int_vals_.rows() << " - " << int_vals_.cols() << std::endl;
303  for (uint i_row=0; i_row<n_points_; ++i_row) {
304  for (uint i_col=0; i_col<3; ++i_col)
305  std::cout << int_vals_(i_col)(i_row) << " ";
306  std::cout << std::endl;
307  }
308  std::cout << std::endl;
309  }
310  std::cout << "*****************" << std::endl;
311  }
312 
313 protected:
314  /**
315  * Store data of bulk or side quadrature points of one dimension
316  *
317  * Number of columns is given by n_rows_, number of used rows by n_points_.
318  */
320  /**
321  * Hold integer values of quadrature points of previous table.
322  *
323  * Table contains following columns:
324  * 0: Index of quadrature point on patch
325  * 1: Row of element/side in point_vals_ table in registration step (before expansion)
326  * 2: Element idx in Mesh
327  * 3: Index of side in element (column is allocated only for side point table)
328  * Number of used rows is given by n_points_.
329  */
331 
332  /// Vector of all defined operations
334 
335  uint dim_; ///< Dimension
336  uint n_rows_; ///< Number of columns of \p point_vals table
337  uint n_points_; ///< Number of points in patch
338  uint n_elems_; ///< Number of elements in patch
339  Quadrature *quad_; ///< Quadrature of given dimension and order passed in constructor.
340 
341  std::vector<uint> elements_map_; ///< Map of element patch indices to el_vals_ table
342  std::vector<uint> points_map_; ///< Map of point patch indices to point_vals_ and int_vals_ tables
343 
344  friend class PatchFEValues<spacedim>;
345  friend class ElOp<spacedim>;
346  template<unsigned int dim>
347  friend class BulkValues;
348  template<unsigned int dim>
349  friend class SideValues;
350  template<unsigned int dim>
351  friend class JoinValues;
352 };
353 
354 /**
355  * @brief Class represents FE operations.
356  */
357 template<unsigned int spacedim = 3>
358 class ElOp {
359 public:
360  /// Constructor
361  ElOp(uint dim, std::initializer_list<uint> shape, uint result_row, ReinitFunction reinit_f, std::vector<uint> input_ops = {}, uint n_dofs = 1)
363  {}
364 
365  /// Number of components computed from shape_ vector
366  inline uint n_comp() const {
367  if (shape_.size() == 1) return shape_[0];
368  else return shape_[0] * shape_[1];
369  }
370 
371  /// Getter of dimension
372  inline uint dim() const {
373  return dim_;
374  }
375 
376  /// Getter of result_row_
377  inline uint result_row() const {
378  return result_row_;
379  }
380 
381  /// Getter of input_ops_
382  inline const std::vector<uint> &input_ops() const {
383  return input_ops_;
384  }
385 
386  /// Getter of shape_
387  inline const std::vector<uint> &shape() const {
388  return shape_;
389  }
390 
391  /// Call reinit function on element table if function is defined
392  inline void reinit_function(std::vector<ElOp<spacedim>> &operations, TableDbl &data_table, TableInt &int_table) {
393  reinit_func(operations, data_table, int_table);
394  }
395 
396  /// Return map referenced Eigen::Matrix of given dimension
397  template<unsigned int dim1, unsigned int dim2>
398  Eigen::Map<Eigen::Matrix<ArrayDbl, dim1, dim2>> value(TableDbl &op_results, uint i_dof = 0) const {
399  return Eigen::Map<Eigen::Matrix<ArrayDbl, dim1, dim2>>(op_results.data() + result_row_ + i_dof * n_comp(), dim1, dim2);
400  }
401 
402 
403 protected:
404  uint dim_; ///< Dimension
405  std::vector<uint> shape_; ///< Shape of stored data (size of vector or number of rows and cols of matrix)
406  uint result_row_; ///< First row to scalar, vector or matrix result
407  std::vector<uint> input_ops_; ///< Indices of operations in PatchPointValues::operations_ vector on which ElOp is depended
408  uint n_dofs_; ///< Number of DOFs of FE operations (or 1 in case of element operations)
409 
410  ReinitFunction reinit_func; ///< Pointer to patch reinit function of element data table specialized by operation
411 };
412 
413 
414 /// Defines common functionality of reinit operations.
416  // empty base operation
417  static inline void op_base(FMT_UNUSED std::vector<ElOp<3>> &operations, FMT_UNUSED TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
418  // empty
419  }
420 
421  // expansion operation
422  static inline void expand_data(ElOp<3> &op, TableDbl &op_results, TableInt &el_table) {
423  uint row_begin = op.result_row();
424  uint row_end = row_begin + op.n_comp();
425  uint size = op_results(row_begin).rows();
426  for (int i_pt=size-1; i_pt>=0; --i_pt) {
427  uint el_table_idx = el_table(1)(i_pt);
428  for (uint i_q=row_begin; i_q<row_end; ++i_q)
429  op_results(i_q)(i_pt) = op_results(i_q)(el_table_idx);
430  }
431  }
432 };
433 
434 /// Defines reinit operations on bulk points.
435 struct bulk_reinit {
436  // element operations
437  template<unsigned int dim>
438  static inline void elop_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
439  // result matrix(spacedim, dim), input matrix(spacedim, dim+1)
440  auto &op = operations[FeBulk::BulkOps::opJac];
441  auto jac_value = op.value<3, dim>(op_results);
442  auto coords_value = operations[ op.input_ops()[0] ].value<3, dim+1>(op_results);
443  jac_value = eigen_tools::jacobian<3,dim>(coords_value);
444  }
445  template<unsigned int dim>
446  static inline void elop_inv_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
447  // result matrix(spacedim, dim), input matrix(spacedim, dim+1)
448  auto &op = operations[FeBulk::BulkOps::opInvJac];
449  auto inv_jac_value = op.value<dim, 3>(op_results);
450  auto jac_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
451  inv_jac_value = eigen_tools::inverse<3, dim>(jac_value);
452  }
453  template<unsigned int dim>
454  static inline void elop_jac_det(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
455  // result double, input matrix(spacedim, dim)
456  auto &op = operations[FeBulk::BulkOps::opJacDet];
457  auto &jac_det_value = op_results(op.result_row());
458  auto jac_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
459  jac_det_value = eigen_tools::determinant<Eigen::Matrix<ArrayDbl, 3, dim>>(jac_value);
460  }
461 
462  // expansion operations
463  static inline void expd_coords(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
464  auto &op = operations[FeBulk::BulkOps::opExpdCoords];
465  common_reinit::expand_data(op, op_results, el_table);
466  }
467  static inline void expd_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
468  auto &op = operations[FeBulk::BulkOps::opExpdJac];
469  common_reinit::expand_data(op, op_results, el_table);
470  }
471  static inline void expd_inv_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
472  auto &op = operations[FeBulk::BulkOps::opExpdInvJac];
473  common_reinit::expand_data(op, op_results, el_table);
474  }
475  static inline void expd_jac_det(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
476  auto &op = operations[FeBulk::BulkOps::opExpdJacDet];
477  common_reinit::expand_data(op, op_results, el_table);
478  }
479 
480  // point operations
481  static inline void ptop_coords(FMT_UNUSED std::vector<ElOp<3>> &operations, FMT_UNUSED TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
482  // Implement
483  }
484  static inline void ptop_weights(std::vector<ElOp<3>> &operations, TableDbl &op_results, std::vector<double> point_weights) {
485  auto &op = operations[FeBulk::BulkOps::opWeights];
486  ArrayDbl &result_row = op_results( op.result_row() );
487  auto n_points = point_weights.size();
488  for (uint i=0; i<result_row.rows(); ++i)
489  result_row(i) = point_weights[i%n_points];
490  }
491  static inline void ptop_JxW(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
492  auto &op = operations[FeBulk::BulkOps::opJxW];
493  ArrayDbl &weights_row = op_results( operations[op.input_ops()[0]].result_row() );
494  ArrayDbl &jac_det_row = op_results( operations[op.input_ops()[1]].result_row() );
495  ArrayDbl &result_row = op_results( op.result_row() );
496  result_row = jac_det_row * weights_row;
497  }
498  static inline void ptop_scalar_shape(std::vector<ElOp<3>> &operations, TableDbl &op_results,
499  std::vector< std::vector<double> > shape_values, uint scalar_shape_op_idx) {
500  auto &op = operations[scalar_shape_op_idx];
501  uint n_points = shape_values.size();
502 
503  for (uint i_row=0; i_row<shape_values[0].size(); ++i_row) {
504  ArrayDbl &result_row = op_results( op.result_row()+i_row );
505  for (uint i_pt=0; i_pt<result_row.rows(); ++i_pt)
506  result_row(i_pt) = shape_values[i_pt % n_points][i_row];
507  }
508  }
509  template<unsigned int dim>
510  static inline void ptop_scalar_shape_grads(std::vector<ElOp<3>> &operations, TableDbl &op_results,
511  std::vector< std::vector<arma::mat> > ref_shape_grads, uint scalar_shape_grads_op_idx) {
512  auto &op = operations[scalar_shape_grads_op_idx];
513  uint n_points = ref_shape_grads.size();
514  uint n_dofs = ref_shape_grads[0].size();
515 
516  Eigen::Vector<ArrayDbl, Eigen::Dynamic> ref_shape_grads_expd;
517  ref_shape_grads_expd.resize(dim*n_dofs);
518  for (uint i=0; i<ref_shape_grads_expd.rows(); ++i)
519  ref_shape_grads_expd(i).resize(op_results(0).rows());
520 
521  for (uint i_dof=0; i_dof<n_dofs; ++i_dof)
522  for (uint i_c=0; i_c<dim; ++i_c) {
523  ArrayDbl &shape_grad_row = ref_shape_grads_expd(i_dof*dim+i_c);
524  for (uint i_pt=0; i_pt<shape_grad_row.rows(); ++i_pt)
525  shape_grad_row(i_pt) = ref_shape_grads[i_pt % n_points][i_dof][i_c];
526  }
527 
528  auto inv_jac_value = operations[ op.input_ops()[0] ].value<dim, 3>(op_results);
529  for (uint i_dof=0; i_dof<n_dofs; ++i_dof) {
530  auto shape_grad_value = op.value<3, 1>(op_results, i_dof);
531  Eigen::Map<Eigen::Matrix<ArrayDbl, dim, 1>> ref_shape_grads_dof_value(ref_shape_grads_expd.data() + dim*i_dof, dim, 1);
532  shape_grad_value = inv_jac_value.transpose() * ref_shape_grads_dof_value;
533  }
534  }
535 };
536 
537 
538 
539 /// Defines reinit operations on side points.
540 struct side_reinit {
541  // element operations
542  template<unsigned int dim>
543  static inline void elop_el_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
544  // result matrix(spacedim, dim), input matrix(spacedim, dim+1)
545  auto &op = operations[FeSide::SideOps::opElJac];
546  auto jac_value = op.value<3, dim>(op_results);
547  auto coords_value = operations[ op.input_ops()[0] ].value<3, dim+1>(op_results);
548  jac_value = eigen_tools::jacobian<3,dim>(coords_value);
549  }
550  template<unsigned int dim>
551  static inline void elop_el_inv_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
552  auto &op = operations[FeSide::SideOps::opElInvJac];
553  auto inv_jac_value = op.value<dim, 3>(op_results);
554  auto jac_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
555  inv_jac_value = eigen_tools::inverse<3,dim>(jac_value);
556  }
557  template<unsigned int dim>
558  static inline void elop_sd_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
559  // result matrix(spacedim, dim), input matrix(spacedim, dim+1)
560  auto &op = operations[FeSide::SideOps::opSdJac];
561  auto jac_value = op.value<3, dim-1>(op_results);
562  auto coords_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
563  jac_value = eigen_tools::jacobian<3, dim-1>(coords_value);
564  }
565  template<unsigned int dim>
566  static inline void elop_sd_jac_det(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
567  // result double, input matrix(spacedim, dim)
568  auto &op = operations[FeSide::SideOps::opSdJacDet];
569  ArrayDbl &det_value = op_results( op.result_row() );
570  auto jac_value = operations[ op.input_ops()[0] ].value<3, dim-1>(op_results);
571  det_value = eigen_tools::determinant<Eigen::Matrix<ArrayDbl, 3, dim-1>>(jac_value);
572  }
573 
574  // expansion operations
575  static inline void expd_el_coords(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
576  auto &op = operations[FeSide::SideOps::opExpdElCoords];
577  common_reinit::expand_data(op, op_results, el_table);
578  }
579  static inline void expd_el_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
580  auto &op = operations[FeSide::SideOps::opExpdElJac];
581  common_reinit::expand_data(op, op_results, el_table);
582  }
583  static inline void expd_el_inv_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
584  auto &op = operations[FeSide::SideOps::opExpdElInvJac];
585  common_reinit::expand_data(op, op_results, el_table);
586  }
587 
588  static inline void expd_sd_coords(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
589  auto &op = operations[FeSide::SideOps::opExpdSdCoords];
590  common_reinit::expand_data(op, op_results, el_table);
591  }
592  static inline void expd_sd_jac(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
593  auto &op = operations[FeSide::SideOps::opExpdSdJac];
594  common_reinit::expand_data(op, op_results, el_table);
595  }
596  static inline void expd_sd_jac_det(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
597  auto &op = operations[FeSide::SideOps::opExpdSdJacDet];
598  common_reinit::expand_data(op, op_results, el_table);
599  }
600 
601  // Point operations
602  static inline void ptop_coords(FMT_UNUSED std::vector<ElOp<3>> &operations, FMT_UNUSED TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
603  // Implement
604  }
605  static inline void ptop_weights(std::vector<ElOp<3>> &operations, TableDbl &op_results, std::vector<double> point_weights) {
606  auto &op = operations[FeSide::SideOps::opWeights];
607  ArrayDbl &result_row = op_results( op.result_row() );
608  auto n_points = point_weights.size();
609  for (uint i=0; i<result_row.rows(); ++i)
610  result_row(i) = point_weights[i%n_points];
611  }
612  static inline void ptop_JxW(std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
613  auto &op = operations[FeSide::SideOps::opJxW];
614  ArrayDbl &weights_row = op_results( operations[op.input_ops()[0]].result_row() );
615  ArrayDbl &jac_det_row = op_results( operations[op.input_ops()[1]].result_row() );
616  ArrayDbl &result_row = op_results( op.result_row() );
617  result_row = jac_det_row * weights_row;
618  }
619  template<unsigned int dim>
620  static inline void ptop_normal_vec(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table) {
621  auto &op = operations[FeSide::SideOps::opNormalVec];
622  auto normal_value = op.value<3, 1>(op_results);
623  auto inv_jac_mat_value = operations[ op.input_ops()[0] ].value<dim, 3>(op_results);
624  normal_value = inv_jac_mat_value.transpose() * RefElement<dim>::normal_vector_array( el_table(3) );
625  }
626  static inline void ptop_scalar_shape(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table,
627  std::vector< std::vector< std::vector<double> > > shape_values, uint scalar_shape_op_idx) {
628  auto &op = operations[scalar_shape_op_idx];
629  uint n_points = shape_values[0].size();
630 
631  for (uint i_row=0; i_row<shape_values[0][0].size(); ++i_row) {
632  ArrayDbl &result_row = op_results( op.result_row()+i_row );
633  for (uint i_pt=0; i_pt<result_row.rows(); ++i_pt)
634  result_row(i_pt) = shape_values[el_table(3)(i_pt)][i_pt % n_points][i_row];
635  }
636  }
637  template<unsigned int dim>
638  static inline void ptop_scalar_shape_grads(std::vector<ElOp<3>> &operations, TableDbl &op_results, TableInt &el_table,
639  std::vector< std::vector< std::vector<arma::mat> > > ref_shape_grads, uint scalar_shape_grads_op_idx) {
640  auto &op = operations[scalar_shape_grads_op_idx];
641  uint n_points = ref_shape_grads[0].size();
642  uint n_dofs = ref_shape_grads[0][0].size();
643 
644  Eigen::Vector<ArrayDbl, Eigen::Dynamic> ref_shape_grads_expd;
645  ref_shape_grads_expd.resize(dim*n_dofs);
646  for (uint i=0; i<ref_shape_grads_expd.rows(); ++i)
647  ref_shape_grads_expd(i).resize(op_results(0).rows());
648 
649  for (uint i_dof=0; i_dof<n_dofs; ++i_dof)
650  for (uint i_c=0; i_c<dim; ++i_c) {
651  ArrayDbl &shape_grad_row = ref_shape_grads_expd(i_dof*dim+i_c);
652  for (uint i_pt=0; i_pt<shape_grad_row.rows(); ++i_pt)
653  shape_grad_row(i_pt) = ref_shape_grads[el_table(3)(i_pt)][i_pt % n_points][i_dof][i_c];
654  }
655 
656  auto inv_jac_value = operations[ op.input_ops()[0] ].value<dim, 3>(op_results);
657  for (uint i_dof=0; i_dof<n_dofs; ++i_dof) {
658  auto shape_grad_value = op.value<3, 1>(op_results, i_dof);
659  Eigen::Map<Eigen::Matrix<ArrayDbl, dim, 1>> ref_shape_grads_dof_value(ref_shape_grads_expd.data() + dim*i_dof, dim, 1);
660  shape_grad_value = inv_jac_value.transpose() * ref_shape_grads_dof_value;
661  }
662  }
663 };
664 
665 
666 
667 namespace FeBulk {
668 
669  /// Bulk data specialization, order of item in operations_ vector corresponds to the BulkOps enum
670  template<unsigned int spacedim = 3>
671  class PatchPointValues : public ::PatchPointValues<spacedim> {
672  public:
673  /// Constructor
675  : ::PatchPointValues<spacedim>(dim) {
676  this->quad_ = new QGauss(dim, 2*quad_order);
677  switch (dim) {
678  case 1:
679  init<1>();
680  break;
681  case 2:
682  init<2>();
683  break;
684  case 3:
685  init<3>();
686  break;
687  }
688  }
689 
690  private:
691  /// Initialize operations vector
692  template<unsigned int dim>
693  void init() {
694  // First step: adds element values operations
695  auto &el_coords = this->make_new_op( {spacedim, this->dim_+1}, &common_reinit::op_base, {} );
696 
697  auto &el_jac = this->make_new_op( {spacedim, this->dim_}, &bulk_reinit::elop_jac<dim>, {BulkOps::opElCoords} );
698 
699  auto &el_inv_jac = this->make_new_op( {this->dim_, spacedim}, &bulk_reinit::elop_inv_jac<dim>, {BulkOps::opJac} );
700 
701  auto &el_jac_det = this->make_new_op( {1}, &bulk_reinit::elop_jac_det<dim>, {BulkOps::opJac} );
702 
703  // Second step: adds expand operations (element values to point values)
704  this->make_expansion( el_coords, {spacedim, this->dim_+1}, &bulk_reinit::expd_coords );
705 
706  this->make_expansion( el_jac, {spacedim, this->dim_}, &bulk_reinit::expd_jac );
707 
708  this->make_expansion( el_inv_jac, {this->dim_, spacedim}, &bulk_reinit::expd_inv_jac );
709 
710  this->make_expansion( el_jac_det, {1}, &bulk_reinit::expd_jac_det );
711 
712  // Third step: adds point values operations
713  /*auto &pt_coords =*/ this->make_new_op( {spacedim}, &bulk_reinit::ptop_coords, {} );
714 
715  // use lambda reinit function
716  std::vector<double> point_weights = this->quad_->get_weights();
717  auto lambda_weights = [point_weights](std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
718  bulk_reinit::ptop_weights(operations, op_results, point_weights);
719  };
720  /*auto &weights =*/ this->make_new_op( {1}, lambda_weights, {} );
721 
722  /*auto &JxW =*/ this->make_new_op( {1}, &bulk_reinit::ptop_JxW, {BulkOps::opWeights, BulkOps::opJacDet} );
723  }
724  };
725 
726 } // closing namespace FeBulk
727 
728 
729 
730 namespace FeSide {
731 
732 /// Bulk Side specialization, order of item in operations_ vector corresponds to the SideOps enum
733  template<unsigned int spacedim = 3>
734  class PatchPointValues : public ::PatchPointValues<spacedim> {
735  public:
736  /// Constructor
738  : ::PatchPointValues<spacedim>(dim) {
739  this->quad_ = new QGauss(dim-1, 2*quad_order);
740  switch (dim) {
741  case 1:
742  init<1>();
743  break;
744  case 2:
745  init<2>();
746  break;
747  case 3:
748  init<3>();
749  break;
750  }
751  }
752 
753  private:
754  /// Initialize operations vector
755  template<unsigned int dim>
756  void init() {
757  // First step: adds element values operations
758  auto &el_coords = this->make_new_op( {spacedim, this->dim_+1}, &common_reinit::op_base, {} );
759 
760  auto &el_jac = this->make_new_op( {spacedim, this->dim_}, &side_reinit::elop_el_jac<dim>, {SideOps::opElCoords} );
761 
762  auto &el_inv_jac = this->make_new_op( {this->dim_, spacedim}, &side_reinit::elop_el_inv_jac<dim>, {SideOps::opElJac} );
763 
764  auto &sd_coords = this->make_new_op( {spacedim, this->dim_}, &common_reinit::op_base, {} );
765 
766  auto &sd_jac = this->make_new_op( {spacedim, this->dim_-1}, &side_reinit::elop_sd_jac<dim>, {SideOps::opSdCoords} );
767 
768  auto &sd_jac_det = this->make_new_op( {1}, &side_reinit::elop_sd_jac_det<dim>, {SideOps::opSdJac} );
769 
770  // Second step: adds expand operations (element values to point values)
771  this->make_expansion( el_coords, {spacedim, this->dim_+1}, &side_reinit::expd_el_coords );
772 
773  this->make_expansion( el_jac, {spacedim, this->dim_}, &side_reinit::expd_el_jac );
774 
775  this->make_expansion( el_inv_jac, {this->dim_, spacedim}, &side_reinit::expd_el_inv_jac );
776 
777  this->make_expansion( sd_coords, {spacedim, this->dim_}, &side_reinit::expd_sd_coords );
778 
779  this->make_expansion( sd_jac, {spacedim, this->dim_-1}, &side_reinit::expd_sd_jac );
780 
781  this->make_expansion( sd_jac_det, {1}, &side_reinit::expd_sd_jac_det );
782 
783  // Third step: adds point values operations
784  /*auto &coords =*/ this->make_new_op( {spacedim}, &side_reinit::ptop_coords, {} );
785 
786  // use lambda reinit function
787  std::vector<double> point_weights = this->quad_->get_weights();
788  auto lambda_weights = [point_weights](std::vector<ElOp<3>> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table) {
789  side_reinit::ptop_weights(operations, op_results, point_weights);
790  };
791  /*auto &weights =*/ this->make_new_op( {1}, lambda_weights, {} );
792 
793  /*auto &JxW =*/ this->make_new_op( {1}, &side_reinit::ptop_JxW, {SideOps::opWeights, SideOps::opSdJacDet} );
794 
795  /*auto &normal_vec =*/ this->make_new_op( {spacedim}, &side_reinit::ptop_normal_vec<dim>, {SideOps::opElInvJac} );
796  }
797  };
798 
799 } // closing namespace FeSide
800 
801 
802 
803 #endif /* PATCH_POINT_VALUES_HH_ */
Class represents FE operations.
std::vector< uint > input_ops_
Indices of operations in PatchPointValues::operations_ vector on which ElOp is depended.
uint dim_
Dimension.
std::vector< uint > shape_
Shape of stored data (size of vector or number of rows and cols of matrix)
uint result_row() const
Getter of result_row_.
const std::vector< uint > & shape() const
Getter of shape_.
const std::vector< uint > & input_ops() const
Getter of input_ops_.
uint result_row_
First row to scalar, vector or matrix result.
ReinitFunction reinit_func
Pointer to patch reinit function of element data table specialized by operation.
ElOp(uint dim, std::initializer_list< uint > shape, uint result_row, ReinitFunction reinit_f, std::vector< uint > input_ops={}, uint n_dofs=1)
Constructor.
uint n_dofs_
Number of DOFs of FE operations (or 1 in case of element operations)
uint dim() const
Getter of dimension.
Eigen::Map< Eigen::Matrix< ArrayDbl, dim1, dim2 > > value(TableDbl &op_results, uint i_dof=0) const
Return map referenced Eigen::Matrix of given dimension.
uint n_comp() const
Number of components computed from shape_ vector.
void reinit_function(std::vector< ElOp< spacedim >> &operations, TableDbl &data_table, TableInt &int_table)
Call reinit function on element table if function is defined.
Bulk data specialization, order of item in operations_ vector corresponds to the BulkOps enum.
PatchPointValues(uint dim, uint quad_order)
Constructor.
void init()
Initialize operations vector.
Bulk Side specialization, order of item in operations_ vector corresponds to the SideOps enum.
PatchPointValues(uint dim, uint quad_order)
Constructor.
void init()
Initialize operations vector.
Class for storing FE data of quadrature points.
ElOp< spacedim > & make_expansion(ElOp< spacedim > &el_op, std::initializer_list< uint > shape, ReinitFunction reinit_f)
Add accessor to operations_ vector.
uint add_rows(uint n_added)
Temporary method allow acces to old structure PatchFEValues::DimPatchFEValues.
void initialize(uint int_cols)
void resize_tables(uint n_points)
Resize data tables.
uint n_rows() const
Getter of 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)
Register side point, add to int_vals_ table.
uint n_rows_
Number of columns of point_vals table.
uint register_element(arma::mat coords, uint element_patch_idx)
Register element, add to point_vals_ table.
Tensor tensor_val(uint result_row, uint point_idx) const
uint register_side(arma::mat elm_coords, arma::mat side_coords)
Register side, add to point_vals_ table.
std::vector< ElOp< spacedim > > operations_
Vector of all defined operations.
Quadrature * get_quadrature() const
Return quadrature.
uint n_elems() const
Getter of n_elems_.
void reset()
Reset number of rows (points and elements)
uint n_points() const
Getter of n_points_.
std::vector< uint > elements_map_
Map of element patch indices to el_vals_ table.
uint dim() const
Getter of dim_.
void print(bool points, bool ints) const
Temporary development method.
ElOp< spacedim > & make_fe_op(std::initializer_list< uint > shape, ReinitFunction reinit_f, std::vector< uint > input_ops_vec, uint n_dofs)
Add accessor to operations_ vector.
Vector vector_val(uint result_row, uint point_idx) const
void reinit_patch()
Reinit data.
uint n_points_
Number of points in patch.
Quadrature * quad_
Quadrature of given dimension and order passed in constructor.
Scalar scalar_val(uint result_row, uint point_idx) const
uint register_bulk_point(uint elem_table_row, uint value_patch_idx, uint elem_idx)
Register bulk point, add to int_vals_ table.
uint n_elems_
Number of elements in patch.
ElOp< spacedim > & make_new_op(std::initializer_list< uint > shape, ReinitFunction reinit_f, std::vector< uint > input_ops_vec)
Add accessor to operations_ vector.
Symmetric Gauss-Legendre quadrature formulae on simplices.
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
const std::vector< double > & get_weights() const
Return a reference to the whole array of weights.
Definition: quadrature.hh:111
static Eigen::Vector< ArrayDbl, dim > normal_vector_array(Eigen::Array< uint, Eigen::Dynamic, 1 > loc_side_idx_array)
Definition: ref_element.cc:279
Store finite element data on the actual patch such as shape function values, gradients,...
Eigen::Vector< ArrayDbl, Eigen::Dynamic > TableDbl
Definition: eigen_tools.hh:33
Eigen::Array< double, Eigen::Dynamic, 1 > ArrayDbl
Definitions of Eigen structures.
Definition: eigen_tools.hh:31
Eigen::Vector< ArrayInt, Eigen::Dynamic > TableInt
Definition: eigen_tools.hh:34
Class ElementValues calculates data related to transformation of reference cell to actual cell (Jacob...
unsigned int uint
ArmaMat< double, N, M > mat
Definition: armor.hh:936
BulkOps
enum of bulk operation
SideOps
enum of side operation
Eigen::Matrix< ArrayDbl, spacedim, dim > jacobian(const Eigen::Matrix< ArrayDbl, spacedim, dim+1 > &coords)
Definition: eigen_tools.hh:228
void resize_table(typename Eigen::Vector< ET, Eigen::Dynamic > &table, uint size)
Resize vector of Eigen::Array to given size.
Definition: eigen_tools.hh:41
ArrayDbl determinant(const T &M)
Calculates determinant of a rectangular matrix.
std::function< void(std::vector< ElOp< 3 > > &, TableDbl &, TableInt &)> ReinitFunction
Type for conciseness.
double Scalar
#define FMT_UNUSED
Definition: posix.h:75
Definitions of particular quadrature rules on simplices.
Defines reinit operations on bulk points.
static void expd_jac_det(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
static void ptop_scalar_shape_grads(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, std::vector< std::vector< arma::mat > > ref_shape_grads, uint scalar_shape_grads_op_idx)
static void ptop_coords(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void ptop_scalar_shape(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, std::vector< std::vector< double > > shape_values, uint scalar_shape_op_idx)
static void elop_jac_det(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void elop_inv_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void elop_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void ptop_weights(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, std::vector< double > point_weights)
static void expd_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
static void ptop_JxW(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void expd_coords(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
static void expd_inv_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
Defines common functionality of reinit operations.
static void expand_data(ElOp< 3 > &op, TableDbl &op_results, TableInt &el_table)
static void op_base(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED TableDbl &op_results, FMT_UNUSED TableInt &el_table)
Defines reinit operations on side points.
static void ptop_normal_vec(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
static void elop_el_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void ptop_weights(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, std::vector< double > point_weights)
static void expd_sd_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
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 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)