ug4
fe_geom_impl.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2013-2015: G-CSC, Goethe University Frankfurt
3  * Author: Andreas Vogel
4  *
5  * This file is part of UG4.
6  *
7  * UG4 is free software: you can redistribute it and/or modify it under the
8  * terms of the GNU Lesser General Public License version 3 (as published by the
9  * Free Software Foundation) with the following additional attribution
10  * requirements (according to LGPL/GPL v3 §7):
11  *
12  * (1) The following notice must be displayed in the Appropriate Legal Notices
13  * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14  *
15  * (2) The following notice must be displayed at a prominent place in the
16  * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17  *
18  * (3) The following bibliography is recommended for citation and must be
19  * preserved in all covered files:
20  * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21  * parallel geometric multigrid solver on hierarchically distributed grids.
22  * Computing and visualization in science 16, 4 (2013), 151-164"
23  * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24  * flexible software system for simulating pde based models on high performance
25  * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26  *
27  * This program is distributed in the hope that it will be useful,
28  * but WITHOUT ANY WARRANTY; without even the implied warranty of
29  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30  * GNU Lesser General Public License for more details.
31  */
32 
33 #include "fe_geom.h"
34 
35 namespace ug{
36 
38 // FEGeometry
40 
41 template < typename TElem, int TWorldDim,
42  typename TTrialSpace, typename TQuadratureRule>
44 FEGeometry()
45 : m_rQuadRule(Provider<quad_rule_type>::get()),
46  m_rTrialSpace(Provider<trial_space_type>::get())
47 {
48  // evaluate local shapes and gradients
49  for(size_t ip = 0; ip < nip; ++ip)
50  for(size_t sh = 0; sh < nsh; ++sh)
51  {
52  m_vvShape[ip][sh] = m_rTrialSpace.shape(sh, m_rQuadRule.point(ip));
53  m_rTrialSpace.grad(m_vvGradLocal[ip][sh], sh, m_rQuadRule.point(ip));
54  }
55 }
56 
57 template < typename TElem, int TWorldDim,
58  typename TTrialSpace, typename TQuadratureRule>
59 void
61 update_local(ReferenceObjectID roid, const LFEID& lfeID, size_t orderQuad)
62 {
64  UG_THROW("FEGeometry::update: Geometry only for "
66  <<roid<<" requested.");
67 
68  if(orderQuad > m_rQuadRule.order())
69  UG_THROW("FEGeometry::update: Geometry only for order "
70  << m_rQuadRule.order()<<", but order "<<
71  orderQuad<<" requested.");
72 }
73 
74 template < typename TElem, int TWorldDim,
75  typename TTrialSpace, typename TQuadratureRule>
76 void
78 update(GridObject* elem, const MathVector<worldDim>* vCorner,
79  const LFEID& lfeID, size_t orderQuad)
80 {
81 // check
82  UG_ASSERT(orderQuad <= m_rQuadRule.order(), "Wrong order requested.");
83 
84  UG_ASSERT(dynamic_cast<TElem*>(elem) != NULL, "Wrong element type.");
85  TElem* pElem = static_cast<TElem*>(elem);
86 
87 // check if element changed
88  if(pElem == m_pElem) return;
89  else m_pElem = pElem;
90 
91 // update the mapping for the new corners
92  m_mapping.update(vCorner);
93 
94 // compute global integration points
95  m_mapping.local_to_global(&m_vIPGlobal[0], local_ips(), nip);
96 
97 // evaluate global data
98  m_mapping.jacobian_transposed_inverse(&m_vJTInv[0], &m_vDetJ[0],
99  local_ips(), nip);
100 
101 // compute global gradients
102  for(size_t ip = 0; ip < nip; ++ip)
103  for(size_t sh = 0; sh < nsh; ++sh)
104  MatVecMult(m_vvGradGlobal[ip][sh],
105  m_vJTInv[ip], m_vvGradLocal[ip][sh]);
106 }
107 
108 } // end namespace ug
TTrialSpace trial_space_type
type of trial space
Definition: fe_geom.h:62
static const size_t nip
number of integration points
Definition: fe_geom.h:71
const quad_rule_type & m_rQuadRule
Quadrature rule.
Definition: fe_geom.h:152
TQuadratureRule quad_rule_type
type of quadrature rule
Definition: fe_geom.h:65
FEGeometry()
Constructor.
Definition: fe_geom_impl.h:44
number m_vvShape[nip][nsh]
shape functions evaluated at ip
Definition: fe_geom.h:164
static const size_t nsh
number of shape functions
Definition: fe_geom.h:68
const trial_space_type & m_rTrialSpace
Quadrature rule.
Definition: fe_geom.h:155
MathVector< dim > m_vvGradLocal[nip][nsh]
global gradient evaluated at ip
Definition: fe_geom.h:167
void update(GridObject *pElem, const MathVector< worldDim > *vCorner)
update Geometry for corners
Definition: fe_geom.h:134
void update_local(ReferenceObjectID roid, const LFEID &lfeID, size_t orderQuad)
update Geometry for roid
Definition: fe_geom_impl.h:61
The base class for all geometric objects, such as vertices, edges, faces, volumes,...
Definition: grid_base_objects.h:157
Identifier for Local Finite Elements.
Definition: local_finite_element_id.h:98
Provider, holding a single instance of an object.
Definition: provider.h:48
Definition: grid_base_object_traits.h:68
#define UG_ASSERT(expr, msg)
Definition: assert.h:70
#define UG_THROW(msg)
Definition: error.h:57
void MatVecMult(vector_t_out &vOut, const matrix_t &m, const vector_t_in &v)
Matrix - Vector Multiplication.
Definition: math_matrix_vector_functions_common_impl.hpp:49
SM_edge_weight_map< typename T::value_type, ug::BidirectionalMatrix< T > > get(edge_weight_t, ug::BidirectionalMatrix< T > const &g)
Definition: bidirectional_boost.h:157
the ug namespace
ReferenceObjectID
these ids are used to identify the shape of a geometric object.
Definition: grid_base_objects.h:74