MAST
geom_elem.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 "mesh/geom_elem.h"
23 #include "mesh/fe_base.h"
24 #include "base/nonlinear_system.h"
26 
27 // libMesh includes
28 #include "libmesh/mesh_base.h"
29 #include "libmesh/boundary_info.h"
30 
31 
33 _sys_init (nullptr),
34 _use_local_elem (false),
35 _ref_elem (nullptr),
36 _local_elem (nullptr) {
37 
38 }
39 
40 
42 
43  if (_local_elem) {
44  delete _local_elem;
45 
46  for (unsigned int i=0; i<_local_nodes.size(); i++)
47  delete _local_nodes[i];
48  }
49 }
50 
51 
52 
53 const libMesh::Elem&
55 
56  libmesh_assert(_ref_elem);
57 
58  return *_ref_elem;
59 }
60 
61 
62 const libMesh::Elem&
64 
65  libmesh_assert(_local_elem);
66 
67  return *_local_elem;
68 }
69 
70 
71 const libMesh::Elem&
73 
74  libmesh_assert(_ref_elem);
75 
76  return *_ref_elem;
77 }
78 
79 
80 const libMesh::Elem&
82 
83  libmesh_assert(_local_elem);
84 
85  return *_local_elem;
86 }
87 
88 
89 
90 unsigned int
92 
93  libmesh_assert(_ref_elem);
94 
95  return _ref_elem->dim();
96 }
97 
98 
99 
100 unsigned int
102 
103  libmesh_assert(_ref_elem);
104 
105  return _ref_elem->n_sides();
106 }
107 
108 
109 libMesh::FEType
110 MAST::GeomElem::get_fe_type(unsigned int i) const {
111 
112  libmesh_assert(_ref_elem);
113 
114  return _sys_init->fetype(i);
115 }
116 
117 
118 void
120 
121  libmesh_assert_equal_to(y_vec.size(), 3);
122 
123  _local_y = y_vec;
124 }
125 
126 
127 void
128 MAST::GeomElem::init(const libMesh::Elem& elem,
129  const MAST::SystemInitialization& sys_init) {
130 
131  libmesh_assert(!_ref_elem);
132 
133  _ref_elem = &elem;
134  _sys_init = &sys_init;
135 
136  // initialize the local element if needed.
138 }
139 
140 
141 std::unique_ptr<MAST::FEBase>
142 MAST::GeomElem::init_fe(bool init_grads,
143  bool init_second_order_derivative,
144  int extra_quadrature_order) const {
145 
146  libmesh_assert(_ref_elem);
147 
148  std::unique_ptr<MAST::FEBase> fe(new MAST::FEBase(*_sys_init));
149  fe->set_extra_quadrature_order(extra_quadrature_order);
150  fe->set_evaluate_second_order_derivatives(init_second_order_derivative);
151 
152  fe->init(*this, init_grads);
153 
154  return fe;
155 }
156 
157 
158 std::unique_ptr<MAST::FEBase>
160  bool init_grads,
161  int extra_quadrature_order) const {
162 
163  std::unique_ptr<MAST::FEBase> fe(new MAST::FEBase(*_sys_init));
164  fe->set_extra_quadrature_order(extra_quadrature_order);
165 
166  fe->init_for_side(*this, s, init_grads);
167 
168  return fe;
169 }
170 
171 
172 
173 void
175 (std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc,
176  std::map<unsigned int, std::vector<MAST::BoundaryConditionBase*>>& loads) const {
177 
178  // this assumes that the quadrature element is the same as the local element
179 
180  typedef std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*> maptype;
181  std::pair<maptype::const_iterator, maptype::const_iterator> range;
182 
183  loads.clear();
184 
185  const libMesh::BoundaryInfo&
186  binfo = *_sys_init->system().get_mesh().boundary_info;
187 
188  std::vector<libMesh::boundary_id_type> bids;
189 
190  for (unsigned int i=0; i<_ref_elem->n_sides(); i++) {
191 
192  bids.clear();
193  binfo.boundary_ids(_ref_elem, i, bids);
194 
195  for (unsigned int j=0; j<bids.size(); j++) {
196 
197  range = bc.equal_range(bids[j]);
198 
199  maptype::const_iterator it = range.first;
200  for ( ; it != range.second; it++) loads[i].push_back(it->second);
201  }
202  }
203 }
204 
205 
206 
207 void
209 (unsigned int s, std::vector<libMesh::boundary_id_type>& bc_ids) const {
210 
211  // this assumes that the quadrature element is the same as the local element
212  bc_ids.clear();
213 
214  _sys_init->system().get_mesh().boundary_info->boundary_ids(_ref_elem, s, bc_ids);
215 }
216 
217 
218 const RealVectorX&
220 
221  libmesh_assert(_ref_elem);
222 
223  return _domain_surface_normal;
224 }
225 
226 
227 const RealMatrixX&
229 
230  libmesh_assert(_local_elem);
231 
232  return _T_mat;
233 }
234 
235 
236 void
238 transform_point_to_global_coordinate(const libMesh::Point& local_pt,
239  libMesh::Point& global_pt) const {
240 
241  libmesh_assert(_local_elem);
242 
243  global_pt.zero();
244 
245  // now calculate the global coordinates with respect to the origin
246  for (unsigned int j=0; j<3; j++)
247  for (unsigned int k=0; k<3; k++)
248  global_pt(j) += _T_mat(j,k)*local_pt(k);
249 
250  // shift to the global coordinate
251  global_pt += (*_ref_elem->node_ptr(0));
252 }
253 
254 
255 
256 void
258 transform_vector_to_global_coordinate(const libMesh::Point& local_vec,
259  libMesh::Point& global_vec) const {
260 
261  libmesh_assert(_local_elem);
262 
263  global_vec.zero();
264 
265  // now calculate the global coordinates with respect to the origin
266  for (unsigned int j=0; j<3; j++)
267  for (unsigned int k=0; k<3; k++)
268  global_vec(j) += _T_mat(j,k)*local_vec(k);
269 }
270 
271 
272 
273 void
275 transform_vector_to_local_coordinate(const libMesh::Point& global_vec,
276  libMesh::Point& local_vec) const {
277 
278  libmesh_assert(_local_elem);
279 
280  local_vec.zero();
281 
282  // now calculate the global coordinates with respect to the origin
283  for (unsigned int j=0; j<3; j++)
284  for (unsigned int k=0; k<3; k++)
285  local_vec(j) += _T_mat(k,j) * global_vec(k);
286 }
287 
288 
289 
290 void
293  RealVectorX& global_vec) const {
294 
295  libmesh_assert(_local_elem);
296 
297  global_vec = _T_mat * local_vec;
298 }
299 
300 
301 void
304  RealVectorX& local_vec) const {
305 
306  libmesh_assert(_local_elem);
307 
308  local_vec = _T_mat.transpose() * global_vec;
309 }
310 
311 
312 
313 bool
315 
316  libmesh_assert(_ref_elem);
317  return _use_local_elem;
318 }
319 
320 
321 
322 void
324 
325  libmesh_assert(_ref_elem);
326  libmesh_assert(!_local_elem);
327 
328  switch (_ref_elem->dim()) {
329 
330  case 1: {
331 
332  // if the y and z cooedinates of the nodes are the same then
333  // the element is along the x-axis and should not need any
334  // local element.
335  const Real
336  y = _ref_elem->point(0)(1),
337  z = _ref_elem->point(0)(2);
338 
339  for (unsigned int i=1; i<_ref_elem->n_nodes(); i++) {
340  if (std::fabs(y-_ref_elem->point(i)(1)) != 0. ||
341  std::fabs(z-_ref_elem->point(i)(2)) != 0.)
342  _use_local_elem = true;
343  }
344 
345  if (_use_local_elem)
347  }
348  break;
349 
350  case 2: {
351 
352  // if the z cooedinate of the nodes are the same then
353  // the element is in the xy-plane and should not need any
354  // local element.
355  const Real
356  z = _ref_elem->point(0)(2);
357 
358  for (unsigned int i=1; i<_ref_elem->n_nodes(); i++) {
359  if (std::fabs(z-_ref_elem->point(i)(2)) != 0.)
360  _use_local_elem = true;
361  }
362 
363  if (_use_local_elem)
365  }
366  break;
367 
368  case 3:
369  _use_local_elem = false;
370  break;
371 
372  default:
373  libmesh_error(); // should not get here.
374  }
375 
376 }
377 
378 
379 void
381 
382  libmesh_assert(_ref_elem);
383  libmesh_assert(_use_local_elem);
384  libmesh_assert(!_local_elem);
385  libmesh_assert(_local_y.size());
386  libmesh_assert_equal_to(_ref_elem->dim(), 1);
387 
388 
389  // first node is the origin of the new cs
390  // calculate the coordinate system for the plane of the element
391  libMesh::Point v1, v2, v3, p;
392  v1 = *_ref_elem->node_ptr(1); v1 -= *_ref_elem->node_ptr(0); v1 /= v1.norm(); // local x
393  for (unsigned int i=0; i<3; i++) v2(i) = _local_y(i); // vector in local x-y plane
394  v3 = v1.cross(v2); // local z
395  libmesh_assert_greater(v3.norm(), 0.); // 0. implies x == y
396  v3 /= v3.norm();
397  v2 = v3.cross(v1); v2 /= v2.norm(); // local y
398 
399  _T_mat = RealMatrixX::Zero(3,3);
400 
401  _local_elem = libMesh::Elem::build(_ref_elem->type()).release();
402  _local_nodes.resize(_ref_elem->n_nodes());
403  for (unsigned int i=0; i<_ref_elem->n_nodes(); i++) {
404  _local_nodes[i] = new libMesh::Node;
405  _local_nodes[i]->set_id() = _ref_elem->node_ptr(i)->id();
406  _local_elem->set_node(i) = _local_nodes[i];
407  }
408 
409  // now the transformation matrix from old to new cs
410  // an_i vn_i = a_i v_i
411  // an_j = a_i v_i.vn_j = a_i t_ij = T^t a_i
412  // t_ij = v_i.vn_j
413 
414  for (unsigned int i=0; i<3; i++) {
415  _T_mat(i,0) = v1(i);
416  _T_mat(i,1) = v2(i);
417  _T_mat(i,2) = v3(i);
418  }
419 
420  // now calculate the new coordinates with respect to the origin
421  for (unsigned int i=0; i<_local_nodes.size(); i++) {
422  p = *_ref_elem->node_ptr(i);
423  p -= *_ref_elem->node_ptr(0); // local wrt origin
424  for (unsigned int j=0; j<3; j++)
425  for (unsigned int k=0; k<3; k++)
426  (*_local_nodes[i])(j) += _T_mat(k,j)*p(k);
427  }
428 }
429 
430 
431 void
433 
434  libmesh_assert(_ref_elem);
435  libmesh_assert(_use_local_elem);
436  libmesh_assert(!_local_elem);
437  libmesh_assert_equal_to(_ref_elem->dim(), 2);
438 
439  // first node is the origin of the new cs
440  // calculate the coordinate system for the plane of the element
441  libMesh::Point v1, v2, v3, p;
442  v1 = *_ref_elem->node_ptr(1); v1 -= *_ref_elem->node_ptr(0); v1 /= v1.norm(); // local x
443  v2 = *_ref_elem->node_ptr(2); v2 -= *_ref_elem->node_ptr(0); v2 /= v2.norm();
444  v3 = v1.cross(v2); v3 /= v3.norm(); // local z
445  v2 = v3.cross(v1); v2 /= v2.norm(); // local y
446 
447  // copy the surface normal
448  _domain_surface_normal = RealVectorX::Zero(3);
449  for (unsigned int i=0; i<3; i++)
450  _domain_surface_normal(i) = v3(i);
451 
452  _T_mat = RealMatrixX::Zero(3,3);
453 
454  _local_elem = libMesh::Elem::build(_ref_elem->type()).release();
455  _local_nodes.resize(_ref_elem->n_nodes());
456  for (unsigned int i=0; i<_ref_elem->n_nodes(); i++) {
457  _local_nodes[i] = new libMesh::Node;
458  _local_nodes[i]->set_id() = _ref_elem->node_ptr(i)->id();
459  _local_elem->set_node(i) = _local_nodes[i];
460  }
461 
462  // now the transformation matrix from old to new cs
463  // an_i vn_i = a_i v_i
464  // an_j = a_i v_i.vn_j = a_i t_ij = T^t a_i
465  // t_ij = v_i.vn_j
466 
467  for (unsigned int i=0; i<3; i++) {
468  _T_mat(i,0) = v1(i);
469  _T_mat(i,1) = v2(i);
470  _T_mat(i,2) = v3(i);
471  }
472 
473  // now calculate the new coordinates with respect to the origin
474  for (unsigned int i=0; i<_local_nodes.size(); i++) {
475  p = *_ref_elem->node_ptr(i);
476  p -= *_ref_elem->node_ptr(0); // local wrt origin
477  for (unsigned int j=0; j<3; j++)
478  for (unsigned int k=0; k<3; k++)
479  (*_local_nodes[i])(j) += _T_mat(k,j)*p(k);
480  }
481 
482 }
bool _use_local_elem
Definition: geom_elem.h:224
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()
const MAST::SystemInitialization * _sys_init
system initialization object for this element
Definition: geom_elem.h:222
void _init_local_elem_1d()
initializes the local element
Definition: geom_elem.cpp:380
virtual void clear() libmesh_override
Clear all the data structures associated with the system.
const libMesh::Elem * _ref_elem
reference element in the mesh for which the data structure is initialized
Definition: geom_elem.h:230
RealVectorX _local_y
the y-axis that is used to define the coordinate system for a 1-D element.
Definition: geom_elem.h:241
void _init_local_elem()
initializes the local element
Definition: geom_elem.cpp:323
RealMatrixX _T_mat
Transformation matrix defines T_ij = V_i^t .
Definition: geom_elem.h:260
void set_local_y_vector(const RealVectorX &y_vec)
for 1D elements the transformed coordinate system attached to the element defines the local x-axis al...
Definition: geom_elem.cpp:119
libMesh::Real Real
virtual const libMesh::Elem & get_reference_elem() const
Definition: geom_elem.cpp:54
libMesh::Elem * _local_elem
a local element is created if
Definition: geom_elem.h:235
const libMesh::FEType & fetype(unsigned int i) const
unsigned int n_sides_quadrature_elem() const
number of sides on quadrature element.
Definition: geom_elem.cpp:101
virtual const libMesh::Elem & get_reference_local_elem() const
Definition: geom_elem.cpp:63
libMesh::FEType get_fe_type(unsigned int i) const
Definition: geom_elem.cpp:110
Matrix< Real, Dynamic, Dynamic > RealMatrixX
const RealMatrixX & T_matrix() const
O.
Definition: geom_elem.cpp:228
RealVectorX _domain_surface_normal
surface normal of the 1D/2D element.
Definition: geom_elem.h:246
Matrix< Real, Dynamic, 1 > RealVectorX
std::vector< libMesh::Node * > _local_nodes
nodes for local element
Definition: geom_elem.h:251
virtual void get_boundary_ids_on_quadrature_elem_side(unsigned int s, std::vector< libMesh::boundary_id_type > &bc_ids) const
Definition: geom_elem.cpp:209
void transform_point_to_global_coordinate(const libMesh::Point &local_pt, libMesh::Point &global_pt) const
Definition: geom_elem.cpp:238
virtual std::unique_ptr< MAST::FEBase > init_side_fe(unsigned int s, bool init_grads, int extra_quadrature_order=0) const
initializes the finite element shape function and quadrature object for the side with the order of qu...
Definition: geom_elem.cpp:159
unsigned int dim() const
Definition: geom_elem.cpp:91
const RealVectorX & domain_surface_normal() const
Definition: geom_elem.cpp:219
virtual void external_side_loads_for_quadrature_elem(std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase * > &bc, std::map< unsigned int, std::vector< MAST::BoundaryConditionBase * >> &loads) const
From the given list of boundary loads, this identifies the sides of the quadrature element and the lo...
Definition: geom_elem.cpp:175
virtual const libMesh::Elem & get_quadrature_elem() const
Definition: geom_elem.cpp:72
void _init_local_elem_2d()
initializes the local element
Definition: geom_elem.cpp:432
void transform_vector_to_local_coordinate(const libMesh::Point &global_vec, libMesh::Point &local_vec) const
Definition: geom_elem.cpp:275
virtual void init(const libMesh::Elem &elem, const MAST::SystemInitialization &sys_init)
initialize the object for the specified reference elem.
Definition: geom_elem.cpp:128
virtual ~GeomElem()
Definition: geom_elem.cpp:41
virtual const libMesh::Elem & get_quadrature_local_elem() const
Definition: geom_elem.cpp:81
virtual std::unique_ptr< MAST::FEBase > init_fe(bool init_grads, bool init_second_order_derivative, int extra_quadrature_order=0) const
initializes the finite element shape function and quadrature object with the order of quadrature rule...
Definition: geom_elem.cpp:142
void transform_vector_to_global_coordinate(const libMesh::Point &local_vec, libMesh::Point &global_vec) const
Definition: geom_elem.cpp:258