Flow123d  DF_patch_fe_darcy_complete-579fe1e
patch_op_impl.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_op_impl.hh
15  * @brief Store finite element reinit functions.
16  * @author David Flanderka
17  */
18 
19 #ifndef OP_FUNCTION_IMPL_HH_
20 #define OP_FUNCTION_IMPL_HH_
21 
22 #include <Eigen/Dense>
23 
24 #include "fem/patch_op.hh"
25 #include "fem/patch_fe_values.hh"
26 
27 
28 template<>
29 template<>
30 inline Scalar PatchOp<3>::elem_value<Scalar>(uint point_idx) const {
31  PatchPointValues<3> &ppv = patch_fe_->patch_point_vals_[domain_][dim_-1];
32  return result_(0)( ppv.int_table_(domain_on_quads)(ppv.points_map_[point_idx]) );
33 }
34 
35 template<>
36 template<>
37 inline Vector PatchOp<3>::elem_value<Vector>(uint point_idx) const {
38  Vector val;
39  PatchPointValues<3> &ppv = patch_fe_->patch_point_vals_[domain_][dim_-1];
40  uint op_matrix_idx = ppv.int_table_(domain_on_quads)(ppv.points_map_[point_idx]);
41  for (uint i=0; i<3; ++i)
42  val(i) = result_(i)(op_matrix_idx);
43  return val;
44 }
45 
46 template<>
47 template<>
48 inline Tensor PatchOp<3>::elem_value<Tensor>(uint point_idx) const {
49  Tensor val;
50  PatchPointValues<3> &ppv = patch_fe_->patch_point_vals_[domain_][dim_-1];
51  uint op_matrix_idx = ppv.int_table_(domain_on_quads)(ppv.points_map_[point_idx]);
52  for (uint i=0; i<3; ++i)
53  for (uint j=0; j<3; ++j)
54  val(i,j) = result_(i+j*3)(op_matrix_idx);
55  return val;
56 }
57 
58 template<>
59 template<>
60 inline Scalar PatchOp<3>::point_value<Scalar>(uint point_idx, uint i_dof) const {
61  PatchPointValues<3> &ppv = patch_fe_->patch_point_vals_[domain_][dim_-1];
62  return result_(i_dof)(ppv.points_map_[point_idx]);
63 }
64 
65 template<>
66 template<>
67 inline Vector PatchOp<3>::point_value<Vector>(uint point_idx, uint i_dof) const {
68  Vector val;
69  PatchPointValues<3> &ppv = patch_fe_->patch_point_vals_[domain_][dim_-1];
70  uint op_matrix_idx = ppv.points_map_[point_idx];
71  for (uint i=0; i<3; ++i)
72  val(i) = result_(i + 3*i_dof)(op_matrix_idx);
73  return val;
74 }
75 
76 template<>
77 template<>
78 inline Tensor PatchOp<3>::point_value<Tensor>(uint point_idx, uint i_dof) const {
79  Tensor val;
80  PatchPointValues<3> &ppv = patch_fe_->patch_point_vals_[domain_][dim_-1];
81  uint op_matrix_idx = ppv.points_map_[point_idx];
82  for (uint i=0; i<9; ++i)
83  val(i) = result_(i+9*i_dof)(op_matrix_idx);
84  return val;
85 }
86 
87 
88 
89 
90 // explicit instantiation of template classes
91 
92 template class PatchOp<3>;
93 
94 
95 #endif /* OP_FUNCTION_IMPL_HH_ */
Class represents element or FE operations.
Definition: patch_op.hh:45
Eigen::Vector< ArenaVec< double >, Eigen::Dynamic > result_
Result matrix of operation.
Definition: patch_op.hh:179
uint dim_
Dimension.
Definition: patch_op.hh:176
ElemDomain domain_
Flag: BulkOp = 0, SideOp = 1.
Definition: patch_op.hh:177
PatchPointValues< spacedim > & ppv() const
Return reference of PatchPointValues.
Definition: patch_op.hh:150
PatchFEValues< spacedim > * patch_fe_
Pointer to PatchFEValues object.
Definition: patch_op.hh:182
unsigned int uint
double Scalar
Definition: op_accessors.hh:25
Class FEValues calculates finite element data on the actual cells such as shape function values,...
Base class of FE operations.
@ domain_on_quads
Index of bulk element or side for each quadrature point in patch.