MAST
sub_cell_fe.cpp
Go to the documentation of this file.
1 /*
2  * MAST: Multidisciplinary-design Adaptation and Sensitivity Toolkit
3  * Copyright (C) 2013-2019 Manav Bhatia
4  *
5  * This library is free software; you can redistribute it and/or
6  * modify it under the terms of the GNU Lesser General Public
7  * License as published by the Free Software Foundation; either
8  * version 2.1 of the License, or (at your option) any later version.
9  *
10  * This library is distributed in the hope that it will be useful,
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13  * Lesser General Public License for more details.
14  *
15  * You should have received a copy of the GNU Lesser General Public
16  * License along with this library; if not, write to the Free Software
17  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
18  */
19 
20 
21 // MAST includes
22 #include "level_set/sub_cell_fe.h"
25 #include "base/nonlinear_system.h"
26 #include "mesh/geom_elem.h"
27 
28 
29 // libMesh includes
30 #include "libmesh/elem.h"
31 
32 
34  const MAST::LevelSetIntersection& intersection):
35 MAST::FEBase (sys),
36 _intersection (intersection) {
37 
38 }
39 
40 
41 
43 
44 }
45 
46 
47 void
49  bool init_grads,
50  const std::vector<libMesh::Point>* pts) {
51 
52  // if there was no intersection, then move back to the parent class's
53  // method
55  MAST::FEBase::init(elem, init_grads, pts);
56  return;
57  }
58 
59  libmesh_assert(!_initialized);
60 
61  // this method does not allow quadrature points to be spcified.
62  libmesh_assert(!pts);
63 
64  _elem = &elem;
65 
66  // the quadrature element is the subcell on which the quadrature is to be
67  // performed and the reference element is the element inside which the
68  // subcell is defined
69  const libMesh::Elem
70  *ref_elem = nullptr,
71  *q_elem = nullptr;
72 
73  if (elem.use_local_elem()) {
74  ref_elem = &elem.get_reference_local_elem();
75  q_elem = &elem.get_quadrature_local_elem();
76  }
77  else {
78  ref_elem = &elem.get_reference_elem();
79  q_elem = &elem.get_quadrature_elem();
80  }
81 
82  const unsigned int
83  nv = _sys.n_vars();
84  libMesh::FEType
85  fe_type = _sys.fetype(0); // all variables are assumed to be of same type
86 
87  for (unsigned int i=1; i != nv; ++i)
88  libmesh_assert(fe_type == _sys.fetype(i));
89 
90  // we first initialize the sub-cell in the nondimensional coordinate system
91  // and then use the non-dimensional location of the quadrature point to
92  // initialize the shape functions of the parent element. After finding the
93  // locations of these points, we reinitialize the sub-cell in the physical
94  // coordinates to get the JxW quantities
95 
96  // initialize the quadrature rule since we will need this later and it
97  // it does not change
98  // the quadrature rule on the subcell that will be used to identify the
99  // quadrature point weight in the
100  std::unique_ptr<libMesh::QBase>
101  subcell_qrule(fe_type.default_quadrature_rule
102  (q_elem->dim(),
103  _sys.system().extra_quadrature_order+_extra_quadrature_order).release());
104  // subcell and local FE are always going to be initialized for linear
105  // shape functions. This is because we create linear sub-elements
106  // and since libmesh does not allow high-order FE on linear elements.
107  // We will initialize the quadrature point locations on these sub elements
108  // and then map it back to the high-order element
109  libMesh::FEType
110  linear_fetype(libMesh::FIRST, libMesh::LAGRANGE);
111  std::unique_ptr<libMesh::FEBase>
112  local_fe (libMesh::FEBase::build(q_elem->dim(), linear_fetype).release()),
113  subcell_fe(libMesh::FEBase::build(q_elem->dim(), linear_fetype).release());
114 
115  // quadrature points in the non-dimensional coordinate system of
116  // reference element are obtained from the local_fe
117  local_fe->get_xyz();
118 
119  // the JxW values are obtained from the quadrature element, and its
120  // sum should be equal to the area of the quadrature element.
121  // initialize the sub-cell FE to get the JxW coordinates
122  subcell_fe->get_JxW();
123 
124  // create the FE object and tell what quantities are needed from it.
125  _fe = libMesh::FEBase::build(q_elem->dim(), fe_type).release();
126  _fe->get_phi();
127  _fe->get_xyz();
128  _fe->get_JxW();
129  if (init_grads) {
130  _fe->get_dphi();
131  _fe->get_dphidxi();
132  _fe->get_dphideta();
133  _fe->get_dphidzeta();
134  }
135  if (_init_second_order_derivatives) _fe->get_d2phi();
136 
137 
139  // create the nondimensional coordinate nodes
141  std::vector<libMesh::Node*>
142  local_nodes(q_elem->n_nodes(), nullptr);
143 
144  std::unique_ptr<libMesh::Elem>
145  local_coord_elem(libMesh::Elem::build(q_elem->type()).release());
146  local_coord_elem->set_id(q_elem->id());
147 
148  // note that the quadrature element is the sub-element that was created
149  // inside the reference element. We intend to use the dofs and shape
150  // functions in the reference element.
151  //
152  // this local element is in the xi-eta coordinate system of the reference
153  // element. Then, if we initialize the quadrature on this local element
154  // then the physical location of these quadrature points will be the
155  // locations in the xi-eta space of the reference element, which we then
156  // use to initialize the FE on reference element.
157  for (unsigned int i=0; i<q_elem->n_nodes(); i++) {
158  const libMesh::Node
159  *n = q_elem->node_ptr(i);
160  const libMesh::Point
162  local_nodes[i] = new libMesh::Node(p);
163  // set the node id since libMesh does not allow invalid ids. This
164  // should not influence the operations here
165  local_nodes[i]->set_id(n->id());
166  local_coord_elem->set_node(i) = local_nodes[i];
167  }
168 
169  // in the nondimensional coordinate
170  local_fe->attach_quadrature_rule(subcell_qrule.get());
171  local_fe->reinit(local_coord_elem.get());
172 
173  // we use the xyz points of the local_elem , which should be the location
174  // of the quadrature points in the reference element nondimensional c/s.
175  // The shape functions and their derivatives are going to come on
176  // the parent elmeent, since the computations use the solution vector
177  // for local computations on that element.
178  _fe->reinit(ref_elem, &local_fe->get_xyz());
179 
180  // We use the _subcell_fe to get the
181  // JxW, since the area needs to come from the element on which
182  // the integration is being performed.
183  subcell_fe->attach_quadrature_rule(subcell_qrule.get());
184  subcell_fe->reinit(q_elem);
185  _subcell_JxW = subcell_fe->get_JxW();
186  _qpoints = local_fe->get_xyz();
187 
188  // transform the normal and
189  if (elem.use_local_elem()) {
190 
191  // now initialize the global xyz locations and normals
192  const std::vector<libMesh::Point>
193  &local_xyz = _fe->get_xyz();
194 
195  unsigned int
196  n = (unsigned int) local_xyz.size();
197  _global_xyz.resize(n);
198 
199  for (unsigned int i=0; i<n; i++)
201  }
202  else
203  _global_xyz = _fe->get_xyz();
204 
205  // now delete the node pointers
206  for (unsigned int i=0; i<q_elem->n_nodes(); i++)
207  delete local_nodes[i];
208 
209  _initialized = true;
210 }
211 
212 
213 
214 void
216  unsigned int s,
217  bool if_calculate_dphi) {
218 
220 
221  // if there was no intersection, then move back to the parent class's
222  // method.
223 
224  MAST::FEBase::init_for_side(elem, s, if_calculate_dphi);
225  return;
226  }
227 
228  // Otherwise, we follow the same procedure as the initialization
229  // over the domain, where the quadrature points from the subcell
230  // are mapped back to the original element and then the original
231  // element is initialized.
232  //
233  // This has to be carefully done. The JxW are needed on the sub-elem,
234  // however, the shape functions (and derivatives) are needed on the
235  // parent element. Additionally, the surface normal are also
236  // evaluated at the sub-elem.
237  //
238 
239 
240  libmesh_assert(!_initialized);
241 
242  _elem = &elem;
243 
244  // the quadrature element is the subcell on which the quadrature is to be
245  // performed and the reference element is the element inside which the
246  // subcell is defined
247  const libMesh::Elem
248  *ref_elem = nullptr,
249  *q_elem = nullptr;
250 
251  if (elem.use_local_elem()) {
252  ref_elem = &elem.get_reference_local_elem();
253  q_elem = &elem.get_quadrature_local_elem();
254  }
255  else {
256  ref_elem = &elem.get_reference_elem();
257  q_elem = &elem.get_quadrature_elem();
258  }
259 
260  const unsigned int
261  nv = _sys.n_vars();
262  libMesh::FEType
263  fe_type = _sys.fetype(0); // all variables are assumed to be of same type
264 
265  for (unsigned int i=1; i != nv; ++i)
266  libmesh_assert(fe_type == _sys.fetype(i));
267 
268  // initialize the quadrature rule since we will need this later and it
269  // it does not change
270  std::unique_ptr<libMesh::QBase>
271  subcell_qrule(fe_type.default_quadrature_rule
272  (q_elem->dim()-1,
273  _sys.system().extra_quadrature_order+_extra_quadrature_order).release());
274 
275  // subcell and local FE are always going to be initialized for linear
276  // shape functions. This is because we create linear sub-elements
277  // and since libmesh does not allow high-order FE on linear elements.
278  // We will initialize the quadrature point locations on these sub elements
279  // and then map it back to the high-order element
280  libMesh::FEType
281  linear_fetype(libMesh::FIRST, libMesh::LAGRANGE);
282 
283  std::unique_ptr<libMesh::FEBase>
284  local_fe (libMesh::FEBase::build(q_elem->dim(), linear_fetype).release()),
285  subcell_fe(libMesh::FEBase::build(q_elem->dim(), linear_fetype).release());
286 
287  // quadrature points in the non-dimensional coordinate system of
288  // reference element are obtained from the local_fe
289  local_fe->get_xyz();
290 
291  // the JxW values are obtained from the quadrature element, and its
292  // sum should be equal to the area of the quadrature element.
293  // initialize the sub-cell FE to get the JxW coordinates
294  subcell_fe->get_JxW();
295  subcell_fe->get_normals();
296 
297  _fe = libMesh::FEBase::build(q_elem->dim(), fe_type).release();
298  _fe->get_phi();
299  _fe->get_xyz();
300  _fe->get_JxW();
301  if (if_calculate_dphi)
302  _fe->get_dphi();
303 
304 
306  // create the nondimensional coordinate nodes
308  std::vector<libMesh::Node*>
309  local_nodes(q_elem->n_nodes(), nullptr);
310 
311  std::unique_ptr<libMesh::Elem>
312  local_coord_elem(libMesh::Elem::build(q_elem->type()).release());
313  local_coord_elem->set_id(q_elem->id());
314 
315  for (unsigned int i=0; i<q_elem->n_nodes(); i++) {
316 
317  const libMesh::Node
318  *n = q_elem->node_ptr(i);
319  const libMesh::Point
321  local_nodes[i] = new libMesh::Node(p);
322  // set the node id since libMesh does not allow invalid ids. This
323  // should not influence the operations here
324  local_nodes[i]->set_id(n->id());
325  local_coord_elem->set_node(i) = local_nodes[i];
326  }
327 
328  // in the nondimensional coordinate, initialize for the side
329  local_fe->attach_quadrature_rule(subcell_qrule.get());
330  local_fe->reinit(local_coord_elem.get(), s);
331 
332  // initialize parent FE using location of these points
333  // The shape functions and their derivatives are going to come on
334  // the parent elmeent, since the computations use the solution vector
335  // for local computations on that element.
336  _fe->reinit(ref_elem, &local_fe->get_xyz());
337 
338  // We use the _subcell_fe to get the
339  // JxW and surface normals, since the area and normals need to come
340  // from the element on which the integration is being performed.
341  subcell_fe->attach_quadrature_rule(subcell_qrule.get());
342  try {
343  subcell_fe->reinit(q_elem, s);
344  _subcell_JxW = subcell_fe->get_JxW();
345  _qpoints = local_fe->get_xyz();
346  // normals in the coordinate system attached to the reference element
347  _local_normals = subcell_fe->get_normals();
348  }
349  catch (...) {
350  unsigned int npts = subcell_qrule->n_points();
351  _subcell_JxW.resize(npts, 0.);
352  _qpoints.resize(npts);
353  _local_normals.resize(npts);
354  }
355 
356  // transform the normal and
357  if (elem.use_local_elem()) {
358 
359  // now initialize the global xyz locations and normals
360  const std::vector<libMesh::Point>
361  &local_xyz = _fe->get_xyz();
362 
363  unsigned int
364  n = (unsigned int) local_xyz.size();
365  _global_xyz.resize(n);
366  _global_normals.resize(n);
367 
368  for (unsigned int i=0; i<n; i++) {
369 
372  }
373  }
374  else {
376  _global_xyz = _fe->get_xyz();
377  }
378 
379  // now delete the node pointers
380  for (unsigned int i=0; i<q_elem->n_nodes(); i++)
381  delete local_nodes[i];
382 
383 
384  _initialized = true;
385 }
386 
387 
388 
389 const
390 std::vector<Real>&
392  libmesh_assert(_initialized);
393 
394  // if there was no intersection, then move back to the parent class's
395  // method
396  if (_subcell_JxW.size())
397  return _subcell_JxW;
398  else
399  return MAST::FEBase::get_JxW();
400 }
401 
402 
403 const std::vector<libMesh::Point>&
405 
406  libmesh_assert(_initialized);
407 
408  if (_subcell_JxW.size())
409  return _qpoints;
410  else
411  return MAST::FEBase::get_qpoints();
412 }
413 
414 
bool use_local_elem() const
Vector and matrix quantities defined on one- and two-dimensional elements that are oriented in two or...
Definition: geom_elem.cpp:314
MAST::NonlinearSystem & system()
std::vector< libMesh::Point > _global_normals
Definition: fe_base.h:189
std::vector< Real > _subcell_JxW
Definition: sub_cell_fe.h:77
virtual const std::vector< Real > & get_JxW() const
Definition: fe_base.cpp:220
const MAST::LevelSetIntersection & _intersection
Definition: sub_cell_fe.h:75
bool _init_second_order_derivatives
Definition: fe_base.h:181
SubCellFE(const MAST::SystemInitialization &sys, const MAST::LevelSetIntersection &intersection)
Definition: sub_cell_fe.cpp:33
virtual const std::vector< libMesh::Point > & get_qpoints() const
Definition: fe_base.cpp:407
const libMesh::Point & get_nondimensional_coordinate_for_node(const libMesh::Node &n) const
std::vector< libMesh::Point > _local_normals
Definition: fe_base.h:188
std::vector< libMesh::Point > _qpoints
Definition: fe_base.h:186
bool _initialized
Definition: fe_base.h:182
virtual const libMesh::Elem & get_reference_elem() const
Definition: geom_elem.cpp:54
virtual void init_for_side(const MAST::GeomElem &elem, unsigned int s, bool if_calculate_dphi)
Initializes the quadrature and finite element for element side integration.
Definition: fe_base.cpp:137
const libMesh::FEType & fetype(unsigned int i) const
const MAST::SystemInitialization & _sys
Definition: fe_base.h:179
virtual const libMesh::Elem & get_reference_local_elem() const
Definition: geom_elem.cpp:63
virtual const std::vector< Real > & get_JxW() const
virtual ~SubCellFE()
Definition: sub_cell_fe.cpp:42
libMesh::FEBase * _fe
Definition: fe_base.h:184
virtual const std::vector< libMesh::Point > & get_qpoints() const
std::vector< libMesh::Point > _global_xyz
Definition: fe_base.h:187
void transform_point_to_global_coordinate(const libMesh::Point &local_pt, libMesh::Point &global_pt) const
Definition: geom_elem.cpp:238
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
Definition: geom_elem.h:59
virtual void init(const MAST::GeomElem &elem, bool init_grads, const std::vector< libMesh::Point > *pts=nullptr)
Initializes the quadrature and finite element for element volume integration.
Definition: fe_base.cpp:64
virtual void init(const MAST::GeomElem &elem, bool init_grads, const std::vector< libMesh::Point > *pts=nullptr)
This assumes that elem is the sub-cell and that the original element will be available as the parent ...
Definition: sub_cell_fe.cpp:48
virtual const libMesh::Elem & get_quadrature_elem() const
Definition: geom_elem.cpp:72
virtual const libMesh::Elem & get_quadrature_local_elem() const
Definition: geom_elem.cpp:81
void transform_vector_to_global_coordinate(const libMesh::Point &local_vec, libMesh::Point &global_vec) const
Definition: geom_elem.cpp:258
const MAST::GeomElem * _elem
Definition: fe_base.h:183
virtual void init_for_side(const MAST::GeomElem &elem, unsigned int s, bool if_calculate_dphi)
This assumes that elem is the sub-cell and that the original element will be available as the parent ...
unsigned int _extra_quadrature_order
Definition: fe_base.h:180