Reference documentation for deal.II version 9.0.0
fe_p1nc.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2015 - 2018 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #include <deal.II/fe/fe_p1nc.h>
18 #include <deal.II/base/std_cxx14/memory.h>
19 
20 DEAL_II_NAMESPACE_OPEN
21 
22 
24  :
25  FiniteElement<2,2>(FiniteElementData<2>(get_dpo_vector(),
26  1,
27  1,
28  FiniteElementData<2>::L2),
29  std::vector<bool> (1, false),
30  std::vector<ComponentMask>(4, ComponentMask(1,true)))
31 {
32  // face support points: 2 end vertices
33  unit_face_support_points.resize(2);
34  unit_face_support_points[0][0] = 0.0;
35  unit_face_support_points[1][0] = 1.0;
36 
37  // initialize constraints matrix
39 }
40 
41 
42 
43 std::string FE_P1NC::get_name () const
44 {
45  return "FE_P1NC";
46 }
47 
48 
49 
51 {
53 
54  if (flags & update_values)
56  if (flags & update_gradients)
57  out |= update_gradients;
58  if (flags & update_cell_normal_vectors)
60  if (flags & update_hessians)
61  out |= update_hessians;
62 
63  return out;
64 }
65 
66 
67 
68 std::unique_ptr<FiniteElement<2,2>>
69  FE_P1NC::clone () const
70 {
71  return std_cxx14::make_unique<FE_P1NC>(*this);
72 }
73 
74 
75 
76 std::vector<unsigned int>
78 {
79  std::vector<unsigned int> dpo(3);
80  dpo[0] = 1; // dofs per object: vertex
81  dpo[1] = 0; // line
82  dpo[2] = 0; // quad
83  return dpo;
84 }
85 
86 
87 
88 std::array<std::array<double,3>,4>
90 {
91  // edge midpoints
92  const Point<2> mpt[4] = { (cell->vertex(0) + cell->vertex(2)) / 2,
93  (cell->vertex(1) + cell->vertex(3)) / 2,
94  (cell->vertex(0) + cell->vertex(1)) / 2,
95  (cell->vertex(2) + cell->vertex(3)) / 2
96  };
97 
98  // center point
99  const Point<2> cpt = (cell->vertex(0) +
100  cell->vertex(1) +
101  cell->vertex(2) +
102  cell->vertex(3)) / 4;
103 
104  const double det = (mpt[0](0)-mpt[1](0))*(mpt[2](1)-mpt[3](1)) - (mpt[2](0)-mpt[3](0))*(mpt[0](1)-mpt[1](1));
105 
106  std::array<std::array<double,3>,4> coeffs;
107  coeffs[0][0] = ((mpt[2](1)-mpt[3](1))*(0.5) -(mpt[0](1)-mpt[1](1))*(0.5))/det;
108  coeffs[1][0] = ((mpt[2](1)-mpt[3](1))*(-0.5) -(mpt[0](1)-mpt[1](1))*(0.5))/det;
109  coeffs[2][0] = ((mpt[2](1)-mpt[3](1))*(0.5) -(mpt[0](1)-mpt[1](1))*(-0.5))/det;
110  coeffs[3][0] = ((mpt[2](1)-mpt[3](1))*(-0.5) -(mpt[0](1)-mpt[1](1))*(-0.5))/det;
111 
112  coeffs[0][1] = (-(mpt[2](0)-mpt[3](0))*(0.5) +(mpt[0](0)-mpt[1](0))*(0.5))/det;
113  coeffs[1][1] = (-(mpt[2](0)-mpt[3](0))*(-0.5) +(mpt[0](0)-mpt[1](0))*(0.5))/det;
114  coeffs[2][1] = (-(mpt[2](0)-mpt[3](0))*(0.5) +(mpt[0](0)-mpt[1](0))*(-0.5))/det;
115  coeffs[3][1] = (-(mpt[2](0)-mpt[3](0))*(-0.5) +(mpt[0](0)-mpt[1](0))*(-0.5))/det;
116 
117  coeffs[0][2] = 0.25 - cpt(0)*coeffs[0][0] - cpt(1)*coeffs[0][1];
118  coeffs[1][2] = 0.25 - cpt(0)*coeffs[1][0] - cpt(1)*coeffs[1][1];
119  coeffs[2][2] = 0.25 - cpt(0)*coeffs[2][0] - cpt(1)*coeffs[2][1];
120  coeffs[3][2] = 0.25 - cpt(0)*coeffs[3][0] - cpt(1)*coeffs[3][1];
121 
122  return coeffs;
123 }
124 
125 
126 
127 std::unique_ptr<FiniteElement<2,2>::InternalDataBase>
128 FE_P1NC::get_data (const UpdateFlags update_flags,
129  const Mapping<2,2> &,
130  const Quadrature<2> &quadrature,
132 {
133  auto data = std_cxx14::make_unique<FiniteElement<2,2>::InternalDataBase>();
134 
135  data->update_each = requires_update_flags(update_flags);
136 
137  const unsigned int n_q_points = quadrature.size();
138  output_data.initialize (n_q_points, FE_P1NC(), data->update_each);
139 
140  // this is a linear element, so its second derivatives are zero
141  if (data->update_each & update_hessians)
142  output_data.shape_hessians.fill(Tensor<2,2>());
143 
144  return data;
145 }
146 
147 
148 
149 std::unique_ptr<FiniteElement<2,2>::InternalDataBase>
150 FE_P1NC::get_face_data (const UpdateFlags update_flags,
151  const Mapping<2,2> &,
152  const Quadrature<1> &quadrature,
154 {
155  auto data = std_cxx14::make_unique<FiniteElement<2,2>::InternalDataBase>();
156 
157  data->update_each = requires_update_flags(update_flags);
158 
159  const unsigned int n_q_points = quadrature.size();
160  output_data.initialize (n_q_points, FE_P1NC(), data->update_each);
161 
162  // this is a linear element, so its second derivatives are zero
163  if (data->update_each & update_hessians)
164  output_data.shape_hessians.fill(Tensor<2,2>());
165 
166  return data;
167 }
168 
169 
170 
171 std::unique_ptr<FiniteElement<2,2>::InternalDataBase>
172 FE_P1NC::get_subface_data (const UpdateFlags update_flags,
173  const Mapping<2,2> &,
174  const Quadrature<1> &quadrature,
176 {
177  auto data = std_cxx14::make_unique<FiniteElement<2,2>::InternalDataBase>();
178 
179  data->update_each = requires_update_flags(update_flags);
180 
181  const unsigned int n_q_points = quadrature.size();
182  output_data.initialize (n_q_points, FE_P1NC(), data->update_each);
183 
184  // this is a linear element, so its second derivatives are zero
185  if (data->update_each & update_hessians)
186  output_data.shape_hessians.fill(Tensor<2,2>());
187 
188  return data;
189 }
190 
191 
192 
193 void
196  const Quadrature<2> &,
197  const Mapping<2,2> &,
200  const FiniteElement<2,2>::InternalDataBase &fe_internal,
202 {
203  const UpdateFlags flags(fe_internal.update_each);
204 
205  const unsigned int n_q_points = mapping_data.quadrature_points.size();
206 
207  // linear shape functions
208  std::array<std::array<double,3>,4> coeffs = get_linear_shape_coefficients (cell);
209 
210  // compute on the cell
211  if (flags & update_values)
212  for (unsigned int i=0; i<n_q_points; ++i)
213  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
214  output_data.shape_values[k][i] = (coeffs[k][0]*mapping_data.quadrature_points[i](0) +
215  coeffs[k][1]*mapping_data.quadrature_points[i](1) +
216  coeffs[k][2]);
217 
218  if (flags & update_gradients)
219  for (unsigned int i=0; i<n_q_points; ++i)
220  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
221  output_data.shape_gradients[k][i] = Point<2>(coeffs[k][0],
222  coeffs[k][1]);
223 }
224 
225 
226 
227 void
229  const unsigned int face_no,
230  const Quadrature<1> &quadrature,
231  const Mapping<2,2> &mapping,
233  const ::internal::FEValuesImplementation::MappingRelatedData<2,2> &,
234  const InternalDataBase &fe_internal,
236 {
237  const UpdateFlags flags(fe_internal.update_each);
238 
239  // linear shape functions
240  const std::array<std::array<double,3>,4> coeffs
242 
243  // compute on the face
244  const Quadrature<2> quadrature_on_face = QProjector<2>::project_to_face(quadrature, face_no);
245 
246  if (flags & update_values)
247  for (unsigned int i=0; i<quadrature_on_face.size(); ++i)
248  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
249  {
250  const Point<2> quadrature_point
251  = mapping.transform_unit_to_real_cell(cell, quadrature_on_face.point(i));
252 
253  output_data.shape_values[k][i] = (coeffs[k][0]*quadrature_point(0) +
254  coeffs[k][1]*quadrature_point(1) +
255  coeffs[k][2]);
256  }
257 
258  if (flags & update_gradients)
259  for (unsigned int i=0; i<quadrature_on_face.size(); ++i)
260  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
261  output_data.shape_gradients[k][i] = Point<2>(coeffs[k][0],
262  coeffs[k][1]);
263 }
264 
265 
266 
267 void
269  const unsigned int face_no,
270  const unsigned int sub_no,
271  const Quadrature<1> &quadrature,
272  const Mapping<2,2> &mapping,
274  const ::internal::FEValuesImplementation::MappingRelatedData<2,2> &,
275  const InternalDataBase &fe_internal,
277 {
278  const UpdateFlags flags(fe_internal.update_each);
279 
280  // linear shape functions
281  const std::array<std::array<double,3>,4> coeffs
283 
284  // compute on the subface
285  const Quadrature<2> quadrature_on_subface = QProjector<2>::project_to_subface(quadrature, face_no, sub_no);
286 
287  if (flags & update_values)
288  for (unsigned int i=0; i<quadrature_on_subface.size(); ++i)
289  {
290  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
291  {
292  const Point<2> quadrature_point
293  = mapping.transform_unit_to_real_cell(cell, quadrature_on_subface.point(i));
294 
295  output_data.shape_values[k][i] = (coeffs[k][0]*quadrature_point(0) +
296  coeffs[k][1]*quadrature_point(1) +
297  coeffs[k][2]);
298  }
299  }
300 
301  if (flags & update_gradients)
302  for (unsigned int i=0; i<quadrature_on_subface.size(); ++i)
303  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
304  output_data.shape_gradients[k][i] = Point<2>(coeffs[k][0],
305  coeffs[k][1]);
306 }
307 
308 
309 
311 {
312  // coefficient relation between children and mother
314  .TableBase<2,double>::reinit (interface_constraints_size());
315 
316  interface_constraints(0,0) = 0.5;
317  interface_constraints(0,1) = 0.5;
318 }
319 
320 
321 DEAL_II_NAMESPACE_CLOSE
Transformed quadrature weights.
Shape function values.
FE_P1NC()
Definition: fe_p1nc.cc:23
virtual std::unique_ptr< FiniteElement< 2, 2 > > clone() const
Definition: fe_p1nc.cc:69
FullMatrix< double > interface_constraints
Definition: fe.h:2401
virtual void fill_fe_values(const Triangulation< 2, 2 >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< 2 > &quadrature, const Mapping< 2, 2 > &mapping, const Mapping< 2, 2 >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< 2, 2 > &mapping_data, const FiniteElement< 2, 2 >::InternalDataBase &fe_internal, internal::FEValuesImplementation::FiniteElementRelatedData< 2, 2 > &output_data) const
Definition: fe_p1nc.cc:194
virtual void fill_fe_subface_values(const Triangulation< 2, 2 >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< 1 > &quadrature, const Mapping< 2, 2 > &mapping, const Mapping< 2, 2 >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< 2, 2 > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< 2, 2 > &output_data) const
Definition: fe_p1nc.cc:268
const Point< dim > & point(const unsigned int i) const
STL namespace.
Transformed quadrature points.
virtual std::unique_ptr< FiniteElement< 2, 2 >::InternalDataBase > get_data(const UpdateFlags update_flags, const Mapping< 2, 2 > &, const Quadrature< 2 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< 2, 2 > &output_data) const
Definition: fe_p1nc.cc:128
static std::vector< unsigned int > get_dpo_vector()
Definition: fe_p1nc.cc:77
No update.
static std::array< std::array< double, 3 >, 4 > get_linear_shape_coefficients(const Triangulation< 2, 2 >::cell_iterator &cell)
Definition: fe_p1nc.cc:89
UpdateFlags
static void project_to_subface(const SubQuadrature &quadrature, const unsigned int face_no, const unsigned int subface_no, std::vector< Point< dim > > &q_points, const RefinementCase< dim-1 > &ref_case=RefinementCase< dim-1 >::isotropic_refinement)
virtual UpdateFlags requires_update_flags(const UpdateFlags flags) const
Definition: fe_p1nc.cc:50
Abstract base class for mapping classes.
Definition: dof_tools.h:46
std::vector< Point< dim-1 > > unit_face_support_points
Definition: fe.h:2420
void initialize(const unsigned int n_quadrature_points, const FiniteElement< dim, spacedim > &fe, const UpdateFlags flags)
Definition: fe_values.cc:2662
Second derivatives of shape functions.
virtual std::string get_name() const
Definition: fe_p1nc.cc:43
unsigned int size() const
const unsigned int dofs_per_cell
Definition: fe_base.h:295
std::vector< Point< spacedim > > quadrature_points
Definition: mpi.h:53
static void project_to_face(const SubQuadrature &quadrature, const unsigned int face_no, std::vector< Point< dim > > &q_points)
Shape function gradients.
virtual void fill_fe_face_values(const Triangulation< 2, 2 >::cell_iterator &cell, const unsigned int face_no, const Quadrature< 1 > &quadrature, const Mapping< 2, 2 > &mapping, const Mapping< 2, 2 >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< 2, 2 > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< 2, 2 > &output_data) const
Definition: fe_p1nc.cc:228
void initialize_constraints()
Definition: fe_p1nc.cc:310
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const =0
TableIndices< 2 > interface_constraints_size() const
Definition: fe.cc:820
void fill(InputIterator entries, const bool C_style_indexing=true)