Loading [MathJax]/extensions/tex2jax.js
ug4
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
hfvcr_geom.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2013-2015: G-CSC, Goethe University Frankfurt
3 * Author: Christian Wehner
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/*
34 * Node centered finite volume geometry for Crouzeix-Raviart-Elements
35 */
36
37#ifndef __H__UG__LIB_DISC__SPATIAL_DISC__DISC_HELPER__HANGING_CR_FINITE_VOLUME_GEOMETRY__
38#define __H__UG__LIB_DISC__SPATIAL_DISC__DISC_HELPER__HANGING_CR_FINITE_VOLUME_GEOMETRY__
39
40// extern libraries
41#include <cmath>
42#include <map>
43#include <vector>
44
45// other ug4 modules
46#include "common/common.h"
47
48// library intern includes
58
59#include "fv_geom_base.h"
60#include "fv_util.h"
61
62namespace ug{
63
64/* hcrfv traits */
65
66template <int TWorldDim,int nrfaceco> struct hcrfv_traits
67{
68 typedef void scv_type;
69 typedef void face_type;
70 static const size_t maxNumSCVF;
71 static const size_t maxNumSCV;
72 static const size_t maxNSH;
73 static const size_t maxNumCo;
74};
75
76template <> struct hcrfv_traits<2, 2>
77{
80 static const size_t maxNumSCVF = 8;
81 static const size_t maxNumSCV = 8;
82 static const size_t maxNSH = maxNumSCV;
83 static const size_t maxNumCo = 4;
84};
85
86template <> struct hcrfv_traits<2, 3>
87{
90 static const size_t maxNumSCVF = 8;
91 static const size_t maxNumSCV = 8;
92 static const size_t maxNSH = maxNumSCV;
93 static const size_t maxNumCo = 4;
94};
95
96template <> struct hcrfv_traits<3, 3>
97{
100 static const size_t maxNumSCVF = 40;
101 static const size_t maxNumSCV = 24;
102 static const size_t maxNSH = maxNumSCV;
103 static const size_t maxNumCo = 8;
104};
105
106template <> struct hcrfv_traits<3, 4>
107{
110 static const size_t maxNumSCVF = 40;
111 static const size_t maxNumSCV = 24;
112 static const size_t maxNSH = maxNumSCV;
113 static const size_t maxNumCo = 8;
114};
115
118// Dimension-indipendent Finite Volume Geometry
121
122template < typename TElem, int TWorldDim>
124{
125 public:
127 typedef TElem elem_type;
128
131
133 static const int dim = ref_elem_type::dim;
134
136 static const int worldDim = TWorldDim;
137
139 static const bool usesHangingNodes = true;
140
142 static const bool staticLocalData = true;
143
146
148 static const size_t nsh = local_shape_fct_set_type::nsh;
149
151 static const size_t numNaturalSCV = nsh;
152
154 static const size_t numNaturalSCVF = ref_elem_type::numEdges;
155
158
159 static const size_t maxNumSCV = traits::maxNumSCV;
160
161 static const size_t maxNumSCVF = traits::maxNumSCVF;
162
163 public:
165 static const int order = 1;
166
172
174 static const size_t nip = 1;
175
176 // local/global element barycenter
179
180 public:
182
191 class SCVF
192 {
193 public:
195 static const size_t numCo = dim;
196
197 public:
198 SCVF() {}
199
201 inline size_t from() const {return From;}
202
204 inline size_t to() const {return To;}
205
207 inline size_t num_ip() const {return nip;}
208
210 inline const MathVector<dim>& local_ip() const {return localIP;}
211
213 inline const MathVector<worldDim>& global_ip() const {return globalIP;}
214
216 inline const MathVector<worldDim>& normal() const {return Normal;}
217
219 inline const MathMatrix<worldDim,dim>& JTInv() const {return JtInv;}
220
222 inline number detJ() const {return detj;}
223
225 inline size_t num_sh() const {return nsh;}
226
228 inline number shape(size_t sh) const {return vShape[sh];}
229
231 inline const number* shape_vector() const {return vShape;}
232
234 inline const MathVector<dim>& local_grad(size_t sh) const
235 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vLocalGrad[sh];}
236
238 inline const MathVector<dim>* local_grad_vector() const {return vLocalGrad;}
239
241 inline const MathVector<worldDim>& global_grad(size_t sh) const
242 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vGlobalGrad[sh];}
243
246
248 inline size_t num_corners() const {return numCo;}
249
251 inline const MathVector<dim>& local_corner(size_t co) const
252 {UG_ASSERT(co < num_corners(), "Invalid index."); return vLocPos[co];}
253
255 inline const MathVector<worldDim>& global_corner(size_t co) const
256 {UG_ASSERT(co < num_corners(), "Invalid index."); return vGloPos[co];}
257
258 private:
259 // let outer class access private members
260 friend class HCRFVGeometry<TElem, TWorldDim>;
261
262 // This scvf separates the scv with the ids given in "from" and "to"
263 // The computed normal points in direction from->to
264 size_t From, To;
265
266 // The normal on the SCVF pointing (from -> to)
267 MathVector<worldDim> Normal; // normal (incl. area)
268
269 MathVector<dim> vLocPos[numCo]; // local corners of scvf
270 MathVector<worldDim> vGloPos[numCo]; // global corners of scvf
271
272 // scvf part
273 MathVector<dim> localIP; // local integration point
274 MathVector<worldDim> globalIP; // global integration point
275
276 // shapes and derivatives
277 number vShape[nsh]; // shapes at ip
278 MathVector<dim> vLocalGrad[nsh]; // local grad at ip
279 MathVector<worldDim> vGlobalGrad[nsh]; // global grad at ip
280 MathMatrix<worldDim,dim> JtInv; // Jacobian transposed at ip
281 number detj; // Jacobian det at ip
282 };
283
285 class SCV
286 {
287 public:
289 static const size_t maxNumCo = 5;
290
291 public:
293
295 inline number volume() const {return Vol;}
296
298 inline size_t num_corners() const {return numCorners;}
299
301 inline const MathVector<dim>& local_corner(size_t co) const
302 {UG_ASSERT(co < num_corners(), "Invalid index."); return vLocPos[co];}
303
305 inline const MathVector<worldDim>& global_corner(size_t co) const
306 {UG_ASSERT(co < num_corners(), "Invalid index."); return vGloPos[co];}
307
309 inline size_t node_id() const {return nodeID;}
310
312 inline size_t num_ip() const {return nip;}
313
315 inline const MathVector<dim>& local_ip() const {return vLocIP;}
316
318 inline const MathVector<worldDim>& global_ip() const {return vGlobIP;}
319
321 inline const MathVector<worldDim>& normal() const {return Normal;}
322
324 inline const MathMatrix<worldDim,dim>& JTInv() const {return JtInv;}
325
327 inline number detJ() const {return detj;}
328
330 inline size_t num_sh() const {return nsh;}
331
333 inline number shape(size_t sh) const {return vShape[sh];}
334
336 inline const number* shape_vector() const {return vShape;}
337
339 inline const MathVector<dim>& local_grad(size_t sh) const
340 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vLocalGrad[sh];}
341
343 inline const MathVector<dim>* local_grad_vector() const {return vLocalGrad;}
344
346 inline const MathVector<worldDim>& global_grad(size_t sh) const
347 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vGlobalGrad[sh];}
348
351
352 private:
353 // let outer class access private members
354 friend class HCRFVGeometry<TElem, TWorldDim>;
355
356 // The normal on the associated face pointing outward
357 MathVector<worldDim> Normal; // normal (incl. area)
358
359 // node id of associated node
360 size_t nodeID;
361
362 // volume of scv
364
365 // number of corners of this element
367
368 // local and global positions of this element
369 MathVector<dim> vLocPos[maxNumCo]; // local position of node
370 MathVector<worldDim> vGloPos[maxNumCo]; // global position of node
371
372 MathVector<dim> vLocIP; // local position of node
373 MathVector<worldDim> vGlobIP; // global position of node
374
375 // shapes and derivatives
376 number vShape[nsh]; // shapes at ip
377 MathVector<dim> vLocalGrad[nsh]; // local grad at ip
378 MathVector<worldDim> vGlobalGrad[nsh]; // global grad at ip
379 MathMatrix<worldDim,dim> JtInv; // Jacobian transposed at ip
380 number detj; // Jacobian det at ip
381 };
382
384 {
385 public:
386 static const size_t maxNumConstrainingDofs = 4;
387 inline size_t constraining_dofs_index(size_t i) const{
388 return cDofInd[i];
389 }
390 inline number constraining_dofs_weight(size_t i) const{
391 return cDofWeights[i];
392 }
393 inline size_t index() const{
394 return i;
395 }
396 inline size_t num_constraining_dofs() const{
397 return numConstrainingDofs;
398 }
399 private:
400 // let outer class access private members
401 friend class HCRFVGeometry<TElem, TWorldDim>;
402
403 // constraining dofs indices
405 // weights
407 // local index of dof in element
408 size_t i;
409 // nr of constraining dofs
411 };
412
414 {
415 public:
416 size_t index;
417 size_t associatedSCV[2];
418 size_t scvfIndex;
419 // indicates if the handled side is from or to
420 bool from;
421 };
422
423 public:
426
428 void update_local_data();
429
431 void update(GridObject* elem, const MathVector<worldDim>* vCornerCoords,
432 const ISubsetHandler* ish = NULL);
433
435 void print();
436
437 const MathVector<worldDim>* corners() const {return m_vCo;}
438
440 inline size_t num_scvf() const {return numSCVF;};
441
443 inline const SCVF& scvf(size_t i) const
444 {UG_ASSERT(i < maxNumSCVF, "Invalid Index."); return m_vSCVF[i];}
445
447 inline size_t num_scv() const {return numSCV;}
448
450 inline const SCV& scv(size_t i) const
451 {UG_ASSERT(i < maxNumSCV, "Invalid Index."); return m_vSCV[i];}
452
454 inline size_t num_constrained_dofs() const {return numConstrainedDofs;}
455
457 inline const CONSTRAINED_DOF& constrained_dof(size_t i) const
458 {UG_ASSERT(i < numConstrainedDofs, "Invalid Index."); return m_vCD[i];}
459
461 inline size_t num_sh() const {return nsh;};
462
463 public:
465 size_t num_scvf_ips() const {return numSCVF;}
466
469
472
474 size_t num_scv_ips() const {return numSCV;}
475
478
481
483 const MathVector<dim> local_bary() const {return localBary;}
484
487
488 protected:
489 // global and local ips on SCVF
492 // coord of location for unknowns in faces (edge/face barycenter)
495
496 static const size_t numMaxCo = 8;
497 // corner coordinates
499
500 private:
502
505
508
511
512 std::vector<HandledEdge> handledEdges;
513
515 TElem* m_pElem;
516
519
522
525
526 size_t numSCV;
527 size_t numSCVF;
529 // numDofs number of all dofs including constraining and constrained dofs
530 size_t numDofs;
531
533
534 static const size_t deleted = 117;
535
536};
537
538
539}
540
541#endif /* __H__UG__LIB_DISC__SPATIAL_DISC__DISC_HELPER__HANGING_CR_FINITE_VOLUME_GEOMETRY__ */
Lagrange Shape Function Set without virtual functions and fixed order.
Definition crouzeix_raviart.h:92
base class for all FVGeometries
Definition fv_geom_base.h:76
The base class for all geometric objects, such as vertices, edges, faces, volumes,...
Definition grid_base_objects.h:157
Definition hfvcr_geom.h:384
size_t index() const
Definition hfvcr_geom.h:393
size_t cDofInd[maxNumConstrainingDofs]
Definition hfvcr_geom.h:404
size_t i
Definition hfvcr_geom.h:408
size_t num_constraining_dofs() const
Definition hfvcr_geom.h:396
number cDofWeights[maxNumConstrainingDofs]
Definition hfvcr_geom.h:406
size_t constraining_dofs_index(size_t i) const
Definition hfvcr_geom.h:387
number constraining_dofs_weight(size_t i) const
Definition hfvcr_geom.h:390
size_t numConstrainingDofs
Definition hfvcr_geom.h:410
static const size_t maxNumConstrainingDofs
Definition hfvcr_geom.h:386
Definition hfvcr_geom.h:414
bool from
Definition hfvcr_geom.h:420
size_t scvfIndex
Definition hfvcr_geom.h:418
size_t index
Definition hfvcr_geom.h:416
size_t associatedSCV[2]
Definition hfvcr_geom.h:417
Sub-Control Volume Face structure.
Definition hfvcr_geom.h:192
size_t from() const
index of SubControlVolume on one side of the scvf
Definition hfvcr_geom.h:201
size_t num_corners() const
number of corners, that bound the scvf
Definition hfvcr_geom.h:248
MathVector< dim > vLocPos[numCo]
Definition hfvcr_geom.h:269
const MathMatrix< worldDim, dim > & JTInv() const
Transposed Inverse of Jacobian in integration point.
Definition hfvcr_geom.h:219
size_t From
Definition hfvcr_geom.h:264
const number * shape_vector() const
vector of shape functions in ip point
Definition hfvcr_geom.h:231
const MathVector< dim > & local_grad(size_t sh) const
value of local gradient of shape function i in integration point
Definition hfvcr_geom.h:234
const MathVector< worldDim > & global_grad(size_t sh) const
value of global gradient of shape function i in integration point
Definition hfvcr_geom.h:241
const MathVector< worldDim > & global_corner(size_t co) const
return glbal corner number i
Definition hfvcr_geom.h:255
number vShape[nsh]
Definition hfvcr_geom.h:277
size_t num_sh() const
number of shape functions
Definition hfvcr_geom.h:225
static const size_t numCo
Number of corners of scvf.
Definition hfvcr_geom.h:195
const MathVector< worldDim > * global_grad_vector() const
vector of global gradients in ip point
Definition hfvcr_geom.h:245
MathVector< worldDim > vGlobalGrad[nsh]
Definition hfvcr_geom.h:279
size_t to() const
index of SubControlVolume on one side of the scvf
Definition hfvcr_geom.h:204
const MathVector< worldDim > & normal() const
normal on scvf (points direction "from"->"to"). Norm is equal to area
Definition hfvcr_geom.h:216
const MathVector< dim > * local_grad_vector() const
vector of local gradients in ip point
Definition hfvcr_geom.h:238
MathVector< worldDim > Normal
Definition hfvcr_geom.h:267
const MathVector< dim > & local_ip() const
local integration point of scvf
Definition hfvcr_geom.h:210
SCVF()
Definition hfvcr_geom.h:198
MathVector< dim > vLocalGrad[nsh]
Definition hfvcr_geom.h:278
size_t To
Definition hfvcr_geom.h:264
number detj
Definition hfvcr_geom.h:281
MathVector< worldDim > vGloPos[numCo]
Definition hfvcr_geom.h:270
number shape(size_t sh) const
value of shape function i in integration point
Definition hfvcr_geom.h:228
number detJ() const
Determinant of Jacobian in integration point.
Definition hfvcr_geom.h:222
const MathVector< dim > & local_corner(size_t co) const
return local corner number i
Definition hfvcr_geom.h:251
size_t num_ip() const
number of integration points on scvf
Definition hfvcr_geom.h:207
const MathVector< worldDim > & global_ip() const
global integration point of scvf
Definition hfvcr_geom.h:213
MathVector< dim > localIP
Definition hfvcr_geom.h:273
MathVector< worldDim > globalIP
Definition hfvcr_geom.h:274
MathMatrix< worldDim, dim > JtInv
Definition hfvcr_geom.h:280
sub control volume structure
Definition hfvcr_geom.h:286
const MathMatrix< worldDim, dim > & JTInv() const
Transposed Inverse of Jacobian in integration point.
Definition hfvcr_geom.h:324
const MathVector< worldDim > & global_ip() const
global integration point
Definition hfvcr_geom.h:318
SCV()
Definition hfvcr_geom.h:292
const MathVector< dim > & local_corner(size_t co) const
return local corner number i
Definition hfvcr_geom.h:301
const number * shape_vector() const
vector of shape functions in ip point
Definition hfvcr_geom.h:336
const MathVector< worldDim > * global_grad_vector() const
vector of global gradients in ip point
Definition hfvcr_geom.h:350
number detj
Definition hfvcr_geom.h:380
number volume() const
volume of scv
Definition hfvcr_geom.h:295
const MathVector< dim > & local_grad(size_t sh) const
value of local gradient of shape function i in integration point
Definition hfvcr_geom.h:339
size_t num_ip() const
number of integration points
Definition hfvcr_geom.h:312
size_t node_id() const
node id that this scv is associated to
Definition hfvcr_geom.h:309
const MathVector< dim > & local_ip() const
local integration point of scv
Definition hfvcr_geom.h:315
int numCorners
Definition hfvcr_geom.h:366
const MathVector< dim > * local_grad_vector() const
vector of local gradients in ip point
Definition hfvcr_geom.h:343
MathVector< worldDim > vGlobalGrad[nsh]
Definition hfvcr_geom.h:378
size_t num_sh() const
number of shape functions
Definition hfvcr_geom.h:330
MathMatrix< worldDim, dim > JtInv
Definition hfvcr_geom.h:379
number shape(size_t sh) const
value of shape function i in integration point
Definition hfvcr_geom.h:333
const MathVector< worldDim > & normal() const
normal on scvf (points direction "from"->"to"). Norm is equal to area
Definition hfvcr_geom.h:321
static const size_t maxNumCo
Number of corners of scv.
Definition hfvcr_geom.h:289
number Vol
Definition hfvcr_geom.h:363
MathVector< worldDim > vGlobIP
Definition hfvcr_geom.h:373
size_t num_corners() const
number of corners, that bound the scvf
Definition hfvcr_geom.h:298
MathVector< worldDim > Normal
Definition hfvcr_geom.h:357
number vShape[nsh]
Definition hfvcr_geom.h:376
MathVector< worldDim > vGloPos[maxNumCo]
Definition hfvcr_geom.h:370
MathVector< dim > vLocalGrad[nsh]
Definition hfvcr_geom.h:377
MathVector< dim > vLocIP
Definition hfvcr_geom.h:372
number detJ() const
Determinant of Jacobian in integration point.
Definition hfvcr_geom.h:327
size_t nodeID
Definition hfvcr_geom.h:360
const MathVector< worldDim > & global_grad(size_t sh) const
value of global gradient of shape function i in integration point
Definition hfvcr_geom.h:346
MathVector< dim > vLocPos[maxNumCo]
Definition hfvcr_geom.h:369
const MathVector< worldDim > & global_corner(size_t co) const
return global corner number i
Definition hfvcr_geom.h:305
Definition hfvcr_geom.h:124
ReferenceMapping< ref_elem_type, worldDim > m_mapping
Reference Mapping.
Definition hfvcr_geom.h:518
SCV m_vSCV[maxNumSCV]
SubControlVolumes.
Definition hfvcr_geom.h:507
const MathVector< worldDim > * corners() const
Definition hfvcr_geom.h:437
const SCVF & scvf(size_t i) const
const access to SubControlVolumeFace number i
Definition hfvcr_geom.h:443
MathVector< dim > m_vLocSCVF_IP[maxNumSCVF]
Definition hfvcr_geom.h:491
reference_element_traits< TElem >::reference_element_type ref_elem_type
type of reference element
Definition hfvcr_geom.h:130
MathVector< worldDim > globalBary
Definition hfvcr_geom.h:178
hcrfv_traits< TWorldDim, TWorldDim+1 >::scv_type scv_type1
Definition hfvcr_geom.h:170
void update_local_data()
update local data
Definition hfvcr_geom.cpp:63
size_t num_scv() const
number of SubControlVolumes
Definition hfvcr_geom.h:447
static const size_t nip
number of integration points
Definition hfvcr_geom.h:174
size_t num_scvf() const
number of SubControlVolumeFaces
Definition hfvcr_geom.h:440
MathVector< dim > m_vLocUnkCoords[maxNumSCV]
Definition hfvcr_geom.h:494
MathVector< worldDim > m_vCo[numMaxCo]
Definition hfvcr_geom.h:498
const MathVector< worldDim > global_bary() const
returns global barycenter
Definition hfvcr_geom.h:486
std::vector< HandledEdge > handledEdges
Definition hfvcr_geom.h:512
void update(GridObject *elem, const MathVector< worldDim > *vCornerCoords, const ISubsetHandler *ish=NULL)
update data for given element
Definition hfvcr_geom.cpp:138
static const size_t numMaxCo
Definition hfvcr_geom.h:496
static const int worldDim
dimension of world
Definition hfvcr_geom.h:136
const CONSTRAINED_DOF & constrained_dof(size_t i) const
const access to constrained dof i
Definition hfvcr_geom.h:457
const SCV & scv(size_t i) const
const access to SubControlVolume number i
Definition hfvcr_geom.h:450
static const int order
order
Definition hfvcr_geom.h:165
static const size_t deleted
Definition hfvcr_geom.h:534
static const size_t nsh
number of shape functions
Definition hfvcr_geom.h:148
HCRFVGeometry()
construct object and initialize local values and sizes
Definition hfvcr_geom.cpp:54
size_t num_sh() const
number of shape functions
Definition hfvcr_geom.h:461
size_t numConstrainedDofs
Definition hfvcr_geom.h:528
MathVector< dim > m_ipCoord[maxNumSCVF]
Definition hfvcr_geom.h:501
SCVF m_vSCVF[maxNumSCVF]
SubControlVolumeFaces.
Definition hfvcr_geom.h:504
const MathVector< dim > * scvf_local_ips() const
returns all ips of scvf as they appear in scv loop
Definition hfvcr_geom.h:468
size_t num_scvf_ips() const
returns number of all scvf ips
Definition hfvcr_geom.h:465
size_t num_scv_ips() const
returns number of all scv ips
Definition hfvcr_geom.h:474
MathVector< worldDim > m_vGlobSCVF_IP[maxNumSCVF]
Definition hfvcr_geom.h:490
void print()
debug output
Definition hfvcr_geom.cpp:560
size_t numSCV
Definition hfvcr_geom.h:526
const MathVector< worldDim > * scvf_global_ips() const
returns all ips of scvf as they appear in scv loop
Definition hfvcr_geom.h:471
TElem * m_pElem
pointer to current element
Definition hfvcr_geom.h:515
hcrfv_traits< TWorldDim, TWorldDim >::face_type face_type0
Definition hfvcr_geom.h:169
const MathVector< dim > * scv_local_ips() const
returns all ips of scv as they appear in scv loop
Definition hfvcr_geom.h:477
hcrfv_traits< TWorldDim, TWorldDim >::scv_type scv_type0
traits
Definition hfvcr_geom.h:168
CONSTRAINED_DOF m_vCD[maxNumSCV]
constrained Dofs
Definition hfvcr_geom.h:510
static const size_t maxNumSCVF
Definition hfvcr_geom.h:161
static const size_t numNaturalSCVF
number of SubControlVolumeFaces
Definition hfvcr_geom.h:154
size_t num_constrained_dofs() const
number of constrained dofs
Definition hfvcr_geom.h:454
static const int dim
dimension of reference element
Definition hfvcr_geom.h:133
hcrfv_traits< dim, worldDim > traits
used traits
Definition hfvcr_geom.h:157
const local_shape_fct_set_type & m_rTrialSpace
Shape function set.
Definition hfvcr_geom.h:524
const MathVector< dim > local_bary() const
returns local barycenter
Definition hfvcr_geom.h:483
const MathVector< worldDim > * scv_global_ips() const
returns all ips of scv as they appear in scv loop
Definition hfvcr_geom.h:480
static const size_t numNaturalSCV
number of SubControlVolumes
Definition hfvcr_geom.h:151
size_t numDofs
Definition hfvcr_geom.h:530
size_t numSCVF
Definition hfvcr_geom.h:527
MathVector< worldDim > m_vGlobUnkCoords[maxNumSCV]
Definition hfvcr_geom.h:493
const ref_elem_type & m_rRefElem
Reference Element.
Definition hfvcr_geom.h:521
hcrfv_traits< TWorldDim, TWorldDim+1 >::face_type face_type1
Definition hfvcr_geom.h:171
static const bool usesHangingNodes
Hanging node flag: this geometry supports hanging nodes.
Definition hfvcr_geom.h:139
CrouzeixRaviartLSFS< ref_elem_type > local_shape_fct_set_type
type of Shape function used
Definition hfvcr_geom.h:145
bool localUpdateNecessary
Definition hfvcr_geom.h:532
TElem elem_type
type of element
Definition hfvcr_geom.h:127
MathVector< dim > localBary
Definition hfvcr_geom.h:177
static const size_t maxNumSCV
Definition hfvcr_geom.h:159
static const bool staticLocalData
flag indicating if local data may change
Definition hfvcr_geom.h:142
Definition subset_handler_interface.h:223
A class for fixed size, dense matrices.
Definition math_matrix.h:63
a mathematical Vector with N entries.
Definition math_vector.h:97
Definition reference_element.h:345
Definition reference_mapping.h:65
Definition reference_element.h:546
Definition reference_element.h:445
Definition reference_element.h:494
Definition reference_element.h:394
#define UG_ASSERT(expr, msg)
Definition assert.h:70
double number
Definition types.h:124
the ug namespace
ReferenceTriangle scv_type
Definition hfvcr_geom.h:78
ReferenceEdge face_type
Definition hfvcr_geom.h:79
ReferenceTriangle scv_type
Definition hfvcr_geom.h:88
ReferenceEdge face_type
Definition hfvcr_geom.h:89
ReferenceTriangle face_type
Definition hfvcr_geom.h:99
ReferenceTetrahedron scv_type
Definition hfvcr_geom.h:98
ReferenceQuadrilateral face_type
Definition hfvcr_geom.h:109
ReferencePyramid scv_type
Definition hfvcr_geom.h:108
Definition hfvcr_geom.h:67
static const size_t maxNumSCVF
Definition hfvcr_geom.h:70
static const size_t maxNumSCV
Definition hfvcr_geom.h:71
static const size_t maxNSH
Definition hfvcr_geom.h:72
void scv_type
Definition hfvcr_geom.h:68
static const size_t maxNumCo
Definition hfvcr_geom.h:73
void face_type
Definition hfvcr_geom.h:69
traits for reference elements
Definition reference_element_traits.h:48