8 #ifndef FIELD_FE_IMPL_HH_
9 #define FIELD_FE_IMPL_HH_
22 namespace it = Input::Type;
27 template <
int spacedim,
class Value>
40 template <
int spacedim,
class Value>
52 VecGetArray(*data_vec_, &data_);
54 unsigned int ndofs = max(dh_->fe<1>()->n_dofs(), max(dh_->fe<2>()->n_dofs(), dh_->fe<3>()->n_dofs()));
55 dof_indices =
new unsigned int[ndofs];
63 template <
int spacedim,
class Value>
71 arma::mat::fixed<1,3> im1 = pinv(m1);
79 dh_->get_dof_indices(cell, dof_indices);
81 if (dh_->fe<1>()->is_scalar()) {
83 for (
int i=0; i<dh_->fe<1>()->n_dofs(); i++)
84 value += data_[dof_indices[i]]*fe_values1.
shape_value(i, 0);
85 this->value_(0,0) = value;
90 for (
int i=0; i<dh_->fe<1>()->n_dofs(); i++)
91 value += data_[dof_indices[i]]*fe_values1.
shape_vector(i, 0);
92 for (
int i=0; i<3; i++)
93 this->value_(i,0) = value(i);
96 else if (elm.
dim() == 2) {
97 arma::mat::fixed<3,2> m2;
100 arma::mat::fixed<2,3> im2 = pinv(m2);
108 dh_->get_dof_indices(cell, dof_indices);
110 if (dh_->fe<2>()->is_scalar()) {
112 for (
int i=0; i<dh_->fe<2>()->n_dofs(); i++)
113 value += data_[dof_indices[i]]*fe_values2.
shape_value(i, 0);
114 this->value_(0,0) = value;
119 for (
int i=0; i<dh_->fe<2>()->n_dofs(); i++)
120 value += data_[dof_indices[i]]*fe_values2.
shape_vector(i, 0);
121 for (
int i=0; i<3; i++)
122 this->value_(i,0) = value(i);
137 dh_->get_dof_indices(cell, dof_indices);
139 if (dh_->fe<3>()->is_scalar()) {
141 for (
int i=0; i<dh_->fe<3>()->n_dofs(); i++)
142 value += data_[dof_indices[i]]*fe_values3.
shape_value(i, 0);
143 this->value_(0,0) = value;
148 for (
int i=0; i<dh_->fe<3>()->n_dofs(); i++)
149 value += data_[dof_indices[i]]*fe_values3.
shape_vector(i, 0);
150 for (
int i=0; i<3; i++)
151 this->value_(i,0) = value(i);
155 return this->r_value_;
163 template <
int spacedim,
class Value>
171 if (elm.
dim() == 1) {
173 arma::mat::fixed<1,3> im1 = pinv(m1);
175 dh_->get_dof_indices(cell, dof_indices);
177 for (
int k=0; k<point_list.size(); k++) {
185 Value envelope(value_list[k]);
187 if (dh_->fe<1>()->is_scalar()) {
189 for (
int i=0; i<dh_->fe<1>()->n_dofs(); i++)
190 value += data_[dof_indices[i]]*fe_values1.
shape_value(i, 0);
191 envelope(0,0) = value;
196 for (
int i=0; i<dh_->fe<1>()->n_dofs(); i++)
197 value += data_[dof_indices[i]]*fe_values1.
shape_vector(i, 0);
198 for (
int i=0; i<3; i++)
199 envelope(i,0) = value(i);
203 else if (elm.
dim() == 2) {
204 arma::mat::fixed<3,2> m2;
207 arma::mat::fixed<2,3> im2 = pinv(m2);
209 dh_->get_dof_indices(cell, dof_indices);
211 for (
int k=0; k<point_list.size(); k++) {
219 Value envelope(value_list[k]);
221 if (dh_->fe<2>()->is_scalar()) {
223 for (
int i=0; i<dh_->fe<2>()->n_dofs(); i++)
224 value += data_[dof_indices[i]]*fe_values2.
shape_value(i, 0);
225 envelope(0,0) = value;
230 for (
int i=0; i<dh_->fe<2>()->n_dofs(); i++)
231 value += data_[dof_indices[i]]*fe_values2.
shape_vector(i, 0);
232 for (
int i=0; i<3; i++)
233 envelope(i,0) = value(i);
244 dh_->get_dof_indices(cell, dof_indices);
246 for (
int k=0; k<point_list.size(); k++) {
254 Value envelope(value_list[k]);
256 if (dh_->fe<3>()->is_scalar()) {
258 for (
int i=0; i<dh_->fe<3>()->n_dofs(); i++)
259 value += data_[dof_indices[i]]*fe_values3.
shape_value(i, 0);
260 envelope(0,0) = value;
265 for (
int i=0; i<dh_->fe<3>()->n_dofs(); i++)
266 value += data_[dof_indices[i]]*fe_values3.
shape_vector(i, 0);
267 for (
int i=0; i<3; i++)
268 envelope(i,0) = value(i);
276 template <
int spacedim,
class Value>
279 if (dof_indices !=
nullptr)
delete[] dof_indices;
const Element * element() const
virtual void value_list(const std::vector< Point > &point_list, const ElementAccessor< spacedim > &elm, std::vector< typename Value::return_type > &value_list)
Class FEValues calculates finite element data on the actual cells such as shape function values...
void set_fe_data(const DOFHandlerMultiDim *dh, Mapping< 1, 3 > *map1, Mapping< 2, 3 > *map2, Mapping< 3, 3 > *map3, const Vec *data)
#define ASSERT_EQUAL(a, b)
Provides the numbering of the finite element degrees of freedom on the computational mesh...
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm)
Basic definitions of numerical quadrature rules.
Space< spacedim >::Point Point
const arma::vec::fixed< spacedim > shape_vector(const unsigned int function_no, const unsigned int point_no)
Return the value of the function_no-th shape function at the point_no-th quadrature point...
void set_point(const unsigned int i, const arma::vec::fixed< dim > &p)
Sets individual quadrature point coordinates.
FieldFE(unsigned int n_comp=0)
Implementation.
const double shape_value(const unsigned int function_no, const unsigned int point_no)
Return the value of the function_no-th shape function at the point_no-th quadrature point...
void reinit(ElementFullIter &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Calculates finite element data on the actual cell.
Abstract class for description of finite elements.