18 #ifndef FE_VALUES_MAP_HH_
19 #define FE_VALUES_MAP_HH_
29 template<
class FV,
unsigned int spacedim = 3>
34 FMT_UNUSED const typename FV::FEInternalData &fe_data) {}
38 const typename FV::FEInternalData &fe_data) {
41 for (
unsigned int i = 0; i < fe_data.n_points; i++)
42 for (
unsigned int j = 0; j < fe_data.n_dofs; j++) {
43 fe_values.set_shape_value(i, j, fe_data.ref_shape_values[i][j][0] );
49 const typename FV::FEInternalData &fe_data) {
52 for (
unsigned int i = 0; i < fe_data.n_points; i++)
53 for (
unsigned int j = 0; j < fe_data.n_dofs; j++)
54 fe_values.set_shape_gradient(i, j, trans(elm_values.
inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j] );
62 template<
class FV,
unsigned int spacedim = 3>
67 FMT_UNUSED const typename FV::FEInternalData &fe_data) {}
71 const typename FV::FEInternalData &fe_data) {
74 for (
unsigned int i = 0; i < fe_data.n_points; i++)
75 for (
unsigned int j = 0; j < fe_data.n_dofs; j++)
78 for (
unsigned int c=0; c<spacedim; c++)
79 fe_values.set_shape_value(i, j*spacedim+c, fv_vec(c));
85 const typename FV::FEInternalData &fe_data) {
88 for (
unsigned int i = 0; i < fe_data.n_points; i++)
89 for (
unsigned int j = 0; j < fe_data.n_dofs; j++)
93 for (
unsigned int c=0; c<spacedim; c++)
94 fe_values.set_shape_gradient(i, j*spacedim+c, grads.col(c));
103 template<
class FV,
unsigned int spacedim = 3>
108 FMT_UNUSED const typename FV::FEInternalData &fe_data) {}
112 const typename FV::FEInternalData &fe_data) {
115 for (
unsigned int i = 0; i < fe_data.n_points; i++)
116 for (
unsigned int j = 0; j < fe_data.n_dofs; j++)
119 for (
unsigned int c=0; c<spacedim; c++)
120 fe_values.set_shape_value(i, j*spacedim+c, fv_vec[c]);
126 const typename FV::FEInternalData &fe_data) {
129 for (
unsigned int i = 0; i < fe_data.n_points; i++)
130 for (
unsigned int j = 0; j < fe_data.n_dofs; j++)
133 for (
unsigned int c=0; c<spacedim; c++)
134 fe_values.set_shape_gradient(i, j*spacedim+c, grads.col(c));
143 template<
class FV,
unsigned int spacedim = 3>
148 FMT_UNUSED const typename FV::FEInternalData &fe_data) {}
152 const typename FV::FEInternalData &fe_data) {
155 for (
unsigned int i = 0; i < fe_data.n_points; i++)
156 for (
unsigned int j = 0; j < fe_data.n_dofs; j++)
158 arma::vec fv_vec = fe_data.ref_shape_values[i][j];
159 for (
unsigned int c=0; c<spacedim; c++)
160 fe_values.set_shape_value(i, j*spacedim+c, fv_vec[c]);
166 const typename FV::FEInternalData &fe_data) {
169 for (
unsigned int i = 0; i < fe_data.n_points; i++)
170 for (
unsigned int j = 0; j < fe_data.n_dofs; j++)
173 for (
unsigned int c=0; c<spacedim; c++)
174 fe_values.set_shape_gradient(i, j*spacedim+c, grads.col(c));
183 template<
class FV,
unsigned int spacedim = 3>
188 FMT_UNUSED const typename FV::FEInternalData &fe_data) {}
192 const typename FV::FEInternalData &fe_data) {
195 for (
unsigned int i = 0; i < fe_data.n_points; i++)
196 for (
unsigned int j = 0; j < fe_data.n_dofs; j++)
198 arma::vec fv_vec = fe_data.ref_shape_values[i][j];
199 for (
unsigned int c=0; c<spacedim*spacedim; c++)
200 fe_values.set_shape_value(i, j*spacedim*spacedim+c, fv_vec[c]);
206 const typename FV::FEInternalData &fe_data) {
209 for (
unsigned int i = 0; i < fe_data.n_points; i++)
210 for (
unsigned int j = 0; j < fe_data.n_dofs; j++)
213 for (
unsigned int c=0; c<spacedim*spacedim; c++)
214 fe_values.set_shape_gradient(i, j*spacedim*spacedim+c, grads.col(c));
223 template<
class FV,
unsigned int spacedim = 3>
228 const typename FV::FEInternalData &fe_data) {
231 unsigned int comp_offset = 0;
232 for (
unsigned int f=0; f<fe_values.fe_sys_dofs_.size(); f++)
235 typename FV::FEInternalData vec_fe_data(fe_data, fe_values.fe_sys_dofs_[f], comp_offset, fe_values.fe_sys_n_components_[f]);
236 fe_values.fe_values_vec[f].fill_data(elm_values, vec_fe_data);
238 comp_offset += fe_values.fe_sys_n_components_[f];
244 const typename FV::FEInternalData &fe_data) {
248 unsigned int comp_offset = 0;
249 unsigned int shape_offset = 0;
250 for (
unsigned int f=0; f<fe_values.fe_sys_dofs_.size(); f++)
253 for (
unsigned int i=0; i<fe_data.n_points; i++)
254 for (
unsigned int n=0; n<fe_values.fe_sys_dofs_[f].size(); n++)
255 for (
unsigned int c=0; c<fe_values.fe_sys_n_space_components_[f]; c++)
256 fe_values.set_shape_value( i, shape_offset+fe_values.n_components_*n+comp_offset+c,
257 fe_values.fe_values_vec[f].shape_value(i, n*fe_values.fe_sys_n_space_components_[f]+c) );
259 comp_offset += fe_values.fe_sys_n_space_components_[f];
260 shape_offset += fe_values.fe_sys_dofs_[f].size()*fe_values.n_components_;
266 const typename FV::FEInternalData &fe_data) {
270 unsigned int comp_offset = 0;
271 unsigned int shape_offset = 0;
272 for (
unsigned int f=0; f<fe_values.fe_sys_dofs_.size(); f++)
275 for (
unsigned int i=0; i<fe_data.n_points; i++)
276 for (
unsigned int n=0; n<fe_values.fe_sys_dofs_[f].size(); n++)
277 for (
unsigned int c=0; c<fe_values.fe_sys_n_space_components_[f]; c++)
278 fe_values.set_shape_gradient(
280 shape_offset+fe_values.n_components_*n+comp_offset+c,
281 fe_values.fe_values_vec[f].shape_grad(i, n*fe_values.fe_sys_n_space_components_[f]+c) );
283 comp_offset += fe_values.fe_sys_n_space_components_[f];
284 shape_offset += fe_values.fe_sys_dofs_[f].size()*fe_values.n_components_;
Class for computation of data on cell and side.
arma::mat jacobian(const unsigned int point_no) const
Return Jacobian matrix at point point_no.
double determinant(const unsigned int point_no) const
Return the relative volume change of the cell (Jacobian determinant).
arma::mat inverse_jacobian(const unsigned int point_no) const
Return inverse Jacobian matrix at point point_no.
Helper class allows update values and gradients of FEValues of FEVectorContravariant type.
void update_values(FV &fe_values, const ElementValues< spacedim > &elm_values, const typename FV::FEInternalData &fe_data)
Update shape_values of given FEValues object.
void fill_values_vec(FMT_UNUSED FV &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, FMT_UNUSED const typename FV::FEInternalData &fe_data)
Empty method.
void update_gradients(FV &fe_values, const ElementValues< spacedim > &elm_values, const typename FV::FEInternalData &fe_data)
Update shape_gradients of given FEValues object.
Helper class allows update values and gradients of FEValues of FEVectorPiola type.
void update_gradients(FV &fe_values, const ElementValues< spacedim > &elm_values, const typename FV::FEInternalData &fe_data)
Update shape_gradients of given FEValues object.
void fill_values_vec(FMT_UNUSED FV &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, FMT_UNUSED const typename FV::FEInternalData &fe_data)
Empty method.
void update_values(FV &fe_values, const ElementValues< spacedim > &elm_values, const typename FV::FEInternalData &fe_data)
Update shape_values of given FEValues object.
Helper class allows update values and gradients of FEValues of FEScalar type.
void fill_values_vec(FMT_UNUSED FV &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, FMT_UNUSED const typename FV::FEInternalData &fe_data)
Empty method.
void update_values(FV &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, const typename FV::FEInternalData &fe_data)
Update shape_values of given FEValues object.
void update_gradients(FV &fe_values, const ElementValues< spacedim > &elm_values, const typename FV::FEInternalData &fe_data)
Update shape_gradients of given FEValues object.
Helper class allows update values and gradients of FEValues of FEMixedSystem type.
void update_values(FV &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, const typename FV::FEInternalData &fe_data)
Update shape_values of given FEValues object.
void fill_values_vec(FV &fe_values, const ElementValues< spacedim > &elm_values, const typename FV::FEInternalData &fe_data)
Fill fe_values_vec of components of mixed system FEValues object.
void update_gradients(FV &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, const typename FV::FEInternalData &fe_data)
Update shape_gradients of given FEValues object.
Helper class allows update values and gradients of FEValues of FETensor type.
void update_gradients(FV &fe_values, const ElementValues< spacedim > &elm_values, const typename FV::FEInternalData &fe_data)
Update shape_gradients of given FEValues object.
void fill_values_vec(FMT_UNUSED FV &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, FMT_UNUSED const typename FV::FEInternalData &fe_data)
Empty method.
void update_values(FV &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, const typename FV::FEInternalData &fe_data)
Update shape_values of given FEValues object.
Helper class allows update values and gradients of FEValues of FEVector type.
void fill_values_vec(FMT_UNUSED FV &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, FMT_UNUSED const typename FV::FEInternalData &fe_data)
Empty method.
void update_values(FV &fe_values, FMT_UNUSED const ElementValues< spacedim > &elm_values, const typename FV::FEInternalData &fe_data)
Update shape_values of given FEValues object.
void update_gradients(FV &fe_values, const ElementValues< spacedim > &elm_values, const typename FV::FEInternalData &fe_data)
Update shape_gradients of given FEValues object.
Class ElementValues calculates data related to transformation of reference cell to actual cell (Jacob...
Class FEValues calculates finite element data on the actual cells such as shape function values,...
ArmaMat< double, N, M > mat