Loading [MathJax]/extensions/tex2jax.js
ug4
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
fv1_geom.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2010-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#ifndef __H__UG__LIB_DISC__SPATIAL_DISC__DISC_HELPER__FV1_GEOMETRY__
34#define __H__UG__LIB_DISC__SPATIAL_DISC__DISC_HELPER__FV1_GEOMETRY__
35
36// extern libraries
37#include <cmath>
38#include <map>
39#include <vector>
40
41// other ug4 modules
42#include "common/common.h"
43
44// library intern includes
53#include "fv_util.h"
54#include "fv_geom_base.h"
55
56namespace ug{
57
59// FV1 Geometry for Reference Element Type
61
62// forward declaration
63template < typename TElem, int TWorldDim, bool TCondensed> class FV1Geometry_gen;
64
66
75template <typename TElem, int TWorldDim>
76class FV1Geometry : public FV1Geometry_gen<TElem, TWorldDim, false> {};
77
79
94template <typename TElem, int TWorldDim>
95class FV1CondensedGeometry : public FV1Geometry_gen<TElem, TWorldDim, true> {};
96
98
112template <typename TElem, int TWorldDim, bool TCondensed>
114{
115 public:
117 typedef TElem elem_type;
118
121
124
125 public:
127 static const int dim = ref_elem_type::dim;
128
130 static const int worldDim = TWorldDim;
131
133 static const bool usesHangingNodes = false;
134
136 static const bool staticLocalData = true;
137
139 static const bool condensed_scvf_ips = TCondensed;
140
141 public:
143 static const int order = 1;
144
146 static const size_t numSCV = traits::numSCV;
147
149 typedef typename traits::scv_type scv_type;
150
152 static const size_t numSCVF = traits::numSCVF;
153
156
158 static const size_t nsh = local_shape_fct_set_type::nsh;
159
161 static const size_t nip = 1;
162
163 public:
165
174 class SCVF
175 {
176 public:
178 static const size_t numCo = traits::NumCornersOfSCVF;
179
180 public:
181 SCVF() {}
182
184 inline size_t from() const {return From;}
185
187 inline size_t to() const {return To;}
188
190 inline const MathVector<worldDim>& normal() const {return Normal;}
191
193 inline size_t num_ip() const {return nip;}
194
196 inline const MathVector<dim>& local_ip() const {return localIP;}
197
199 inline const MathVector<worldDim>& global_ip() const {return globalIP;}
200
202 inline const MathMatrix<worldDim,dim>& JTInv() const {return JtInv;}
203
205 inline number detJ() const {return detj;}
206
208 inline size_t num_sh() const {return nsh;}
209
211 inline number shape(size_t sh) const {return vShape[sh];}
212
214 inline const number* shape_vector() const {return vShape;}
215
217 inline const MathVector<dim>& local_grad(size_t sh) const
218 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vLocalGrad[sh];}
219
221 inline const MathVector<dim>* local_grad_vector() const {return vLocalGrad;}
222
224 inline const MathVector<worldDim>& global_grad(size_t sh) const
225 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vGlobalGrad[sh];}
226
229
231 inline size_t num_corners() const {return numCo;}
232
234 inline const MathVector<dim>& local_corner(size_t co) const
235 {UG_ASSERT(co < num_corners(), "Invalid index."); return vLocPos[co];}
236
238 inline const MathVector<worldDim>& global_corner(size_t co) const
239 {UG_ASSERT(co < num_corners(), "Invalid index."); return vGloPos[co];}
240
241 private:
242 // let outer class access private members
243 friend class FV1Geometry_gen<TElem, TWorldDim, TCondensed>;
244
245 // This scvf separates the scv with the ids given in "from" and "to"
246 // The computed normal points in direction from->to
247 size_t From, To;
248
249 // The normal on the SCVF pointing (from -> to)
250 MathVector<worldDim> Normal; // normal (incl. area)
251
252 // ordering is:
253 // 1D: edgeMidPoint
254 // 2D: edgeMidPoint, CenterOfElement
255 // 3D: edgeMidPoint, Side #1, CenterOfElement, Side #2
256 MathVector<dim> vLocPos[numCo]; // local corners of scvf
257 MathVector<worldDim> vGloPos[numCo]; // global corners of scvf
258 MidID vMidID[numCo]; // dimension and id of object, that's midpoint bounds the scvf
259
260 // scvf part
261 MathVector<dim> localIP; // local integration point
262 MathVector<worldDim> globalIP; // global integration point
263
264 // shapes and derivatives
265 number vShape[nsh]; // shapes at ip
266 MathVector<dim> vLocalGrad[nsh]; // local grad at ip
267 MathVector<worldDim> vGlobalGrad[nsh]; // global grad at ip
268 MathMatrix<worldDim,dim> JtInv; // Jacobian transposed at ip
269 number detj; // Jacobian det at ip
270 };
271
273 class SCV
274 {
275 public:
277 static const size_t numCo = traits::NumCornersOfSCV;
278
279 public:
280 SCV() {};
281
283 inline number volume() const {return Vol;}
284
286 inline size_t num_corners() const {return numCo;}
287
289 inline const MathVector<dim>& local_corner(size_t co) const
290 {UG_ASSERT(co < num_corners(), "Invalid index."); return vLocPos[co];}
291
293 inline const MathVector<worldDim>& global_corner(size_t co) const
294 {UG_ASSERT(co < num_corners(), "Invalid index."); return vGloPos[co];}
295
297 inline const MathVector<dim>* local_corners() const
298 {return &vLocPos[0];}
299
302 {return &vGloPos[0];}
303
305 inline size_t node_id() const {return nodeId;}
306
308 inline size_t num_ip() const {return nip;}
309
311 inline const MathVector<dim>& local_ip() const {return vLocPos[0];}
312
314 inline const MathVector<worldDim>& global_ip() const {return vGloPos[0];}
315
317 inline const MathMatrix<worldDim,dim>& JTInv() const {return JtInv;}
318
320 inline number detJ() const {return detj;}
321
323 inline size_t num_sh() const {return nsh;}
324
326 inline number shape(size_t sh) const {return vShape[sh];}
327
329 inline const number* shape_vector() const {return vShape;}
330
332 inline const MathVector<dim>& local_grad(size_t sh) const
333 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vLocalGrad[sh];}
334
336 inline const MathVector<dim>* local_grad_vector() const {return vLocalGrad;}
337
339 inline const MathVector<worldDim>& global_grad(size_t sh) const
340 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vGlobalGrad[sh];}
341
344
345 private:
346 // let outer class access private members
347 friend class FV1Geometry_gen<TElem, TWorldDim, TCondensed>;
348
349 // node id of associated node
350 size_t nodeId;
351
352 // volume of scv
354
355 // local and global positions of this element
356 MathVector<dim> vLocPos[numCo]; // local position of node
357 MathVector<worldDim> vGloPos[numCo]; // global position of node
358 MidID midId[numCo]; // dimension and id of object, that's midpoint bounds the scv
359
360 // shapes and derivatives
361 number vShape[nsh]; // shapes at ip
362 MathVector<dim> vLocalGrad[nsh]; // local grad at ip
363 MathVector<worldDim> vGlobalGrad[nsh]; // global grad at ip
364 MathMatrix<worldDim,dim> JtInv; // Jacobian transposed at ip
365 number detj; // Jacobian det at ip
366 };
367
369 class BF
370 {
371 public:
373 static const size_t numCo = traits::NumCornersOfBF;
374
375 public:
376 BF() {}
377
379 inline size_t node_id() const {return nodeId;}
380
382 inline size_t num_ip() const {return nip;}
383
385 inline const MathVector<dim>& local_ip() const {return localIP;}
386
388 inline const MathVector<worldDim>& global_ip() const {return globalIP;}
389
391 inline const MathVector<worldDim>& normal() const {return Normal;} // includes area
392
394 inline number volume() const {return Vol;}
395
397 inline const MathMatrix<worldDim, dim>& JTInv() const {return JtInv;}
398
400 inline number detJ() const {return detj;}
401
403 inline size_t num_sh() const {return nsh;}
404
406 inline number shape(size_t sh) const
407 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vShape[sh];}
408
410 inline const number* shape_vector() const {return vShape;}
411
413 inline const MathVector<dim>& local_grad(size_t sh) const
414 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vLocalGrad[sh];}
415
417 inline const MathVector<dim>* local_grad_vector() const {return vLocalGrad;}
418
420 inline const MathVector<worldDim>& global_grad(size_t sh) const
421 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vGlobalGrad[sh];}
422
425
427 inline size_t num_corners() const {return numCo;}
428
430 inline const MathVector<dim>& local_corner(size_t co) const
431 {UG_ASSERT(co < num_corners(), "Invalid index."); return vLocPos[co];}
432
434 inline const MathVector<worldDim>& global_corner(size_t co) const
435 {UG_ASSERT(co < num_corners(), "Invalid index."); return vGloPos[co];}
436
437 private:
439 friend class FV1Geometry_gen<TElem, TWorldDim, TCondensed>;
440
441 // id of scv this bf belongs to
442 size_t nodeId;
443
444 // ordering is:
445 // 1D: edgeMidPoint
446 // 2D: edgeMidPoint, CenterOfElement
447 // 3D: edgeMidPoint, Side one, CenterOfElement, Side two
448 MathVector<dim> vLocPos[numCo]; // local corners of bf
449 MathVector<worldDim> vGloPos[numCo]; // global corners of bf
450 MidID vMidID[numCo]; // dimension and id of object, that's midpoint bounds the bf
451
452 // scvf part
453 MathVector<dim> localIP; // local integration point
454 MathVector<worldDim> globalIP; // global integration point
455 MathVector<worldDim> Normal; // normal (incl. area)
456 number Vol; // volume of bf
457
458 // shapes and derivatives
459 number vShape[nsh]; // shapes at ip
460 MathVector<dim> vLocalGrad[nsh]; // local grad at ip
461 MathVector<worldDim> vGlobalGrad[nsh]; // global grad at ip
462 MathMatrix<worldDim,dim> JtInv; // Jacobian transposed at ip
463 number detj; // Jacobian det at ip
464 };
465
466 public:
469
471 void update_local_data();
472
474 void update(GridObject* elem, const MathVector<worldDim>* vCornerCoords,
475 const ISubsetHandler* ish = NULL);
476
479 const MathVector<worldDim>* vCornerCoords,
480 const ISubsetHandler* ish = NULL);
481
483 TElem* elem() const {return m_pElem;}
484
486 const MathVector<worldDim>* corners() const {return m_vvGloMid[0];}
487
489 size_t num_scvf() const {return numSCVF;};
490
492 const SCVF& scvf(size_t i) const
493 {UG_ASSERT(i < num_scvf(), "Invalid Index."); return m_vSCVF[i];}
494
496 // do not use this method to obtain the number of shape functions,
497 // since this is NOT the same for pyramids; use num_sh() instead.
498 size_t num_scv() const {return numSCV;}
499
501 const SCV& scv(size_t i) const
502 {UG_ASSERT(i < num_scv(), "Invalid Index."); return m_vSCV[i];}
503
505 size_t num_sh() const {return nsh;};
506
508 ReferenceObjectID roid() const {return ref_elem_type::REFERENCE_OBJECT_ID;}
509
510
511 public:
513 size_t num_scvf_ips() const {return numSCVF;}
514
517
520
522 size_t num_scv_ips() const {return numSCV;}
523
526 if(ref_elem_type::REFERENCE_OBJECT_ID == ROID_PYRAMID || ref_elem_type::REFERENCE_OBJECT_ID == ROID_OCTAHEDRON)
527 return &(m_vLocSCV_IP[0]);
528 else
529 return &(m_vvLocMid[0][0]);
530 }
531
534 if(ref_elem_type::REFERENCE_OBJECT_ID == ROID_PYRAMID || ref_elem_type::REFERENCE_OBJECT_ID == ROID_OCTAHEDRON)
535 return &(m_vGlobSCV_IP[0]);
536 else
537 return &(m_vvGloMid[0][0]);
538 }
539
541 const MathVector<dim>& local_node_position(size_t nodeID) const
542 {
543 UG_ASSERT(nodeID < (size_t) ref_elem_type::numCorners, "Invalid node id.");
544 return m_vvLocMid[0][nodeID];
545 }
546
548 const MathVector<worldDim>& global_node_position(size_t nodeID) const
549 {
550 UG_ASSERT(nodeID < (size_t) ref_elem_type::numCorners, "Invalid node id.");
551 return m_vvGloMid[0][nodeID];
552 }
553
555 const MathVector<dim>* coe_local() const {return &(m_vvLocMid[dim][0]);}
556
558 const MathVector<worldDim>* coe_global() const {return &(m_vvGloMid[dim][0]);}
559
560 protected:
561 // global and local ips on SCVF
564
565 // global and local ips on SCV (only needed for Pyramid and Octahedron)
568
569 public:
571 inline void add_boundary_subset(int subsetIndex) {m_mapVectorBF[subsetIndex];}
572
574 inline void remove_boundary_subset(int subsetIndex) {m_mapVectorBF.erase(subsetIndex);}
575
577 inline void clear_boundary_subsets() {m_mapVectorBF.clear();}
578
580 inline size_t num_boundary_subsets() {return m_mapVectorBF.size();}
581
583 inline size_t num_bf() const
584 {
585 typename std::map<int, std::vector<BF> >::const_iterator it;
586 size_t num = 0;
587 for ( it=m_mapVectorBF.begin() ; it != m_mapVectorBF.end(); it++ )
588 num += (*it).second.size();
589 return num;
590 }
591
593 inline size_t num_bf(int si) const
594 {
595 typename std::map<int, std::vector<BF> >::const_iterator it;
596 it = m_mapVectorBF.find(si);
597 if(it == m_mapVectorBF.end()) return 0;
598 else return (*it).second.size();
599 }
600
602 inline const BF& bf(int si, size_t i) const
603 {
604 UG_ASSERT(i < num_bf(si), "Invalid index.");
605 typename std::map<int, std::vector<BF> >::const_iterator it;
606 it = m_mapVectorBF.find(si);
607 if(it == m_mapVectorBF.end()) UG_THROW("FVGeom: No bnd face for subset"<<si);
608 return (*it).second[i];
609 }
610
612 inline const std::vector<BF>& bf(int si) const
613 {
614 typename std::map<int, std::vector<BF> >::const_iterator it;
615 it = m_mapVectorBF.find(si);
616 if(it == m_mapVectorBF.end()) return m_vEmptyVectorBF;
617 return (*it).second;
618 }
619
620 void reset_curr_elem() {m_pElem = NULL;}
621
622 protected:
623 std::map<int, std::vector<BF> > m_mapVectorBF;
624 std::vector<BF> m_vEmptyVectorBF;
625
626 private:
628 TElem* m_pElem;
629
631 // (most objects in 1 dim, i.e. number of edges, but +1 for 1D)
632 static const int maxMid = numSCVF + 1;
633
637
640
643
646
649
652};
653
655// Dim-dependent Finite Volume Geometry
657
659
663template <int TDim, int TWorldDim = TDim>
665{
666 public:
669
670 public:
672 static const int dim = TDim;
673
675 static const int worldDim = TWorldDim;
676
678 static const bool usesHangingNodes = false;
679
681 static const bool staticLocalData = false;
682
683 public:
685 static const int order = 1;
686
688 static const size_t maxNumSCV = traits::maxNumSCV;
689
691 typedef typename traits::scv_type scv_type;
692
694 static const size_t maxNumSCVF = traits::maxNumSCVF;
695
697 static const size_t maxNSH = traits::maxNSH;
698
700 static const size_t nip = 1;
701
702 public:
704
713 class SCVF
714 {
715 public:
717 static const size_t numCo = traits::NumCornersOfSCVF;
718
719 public:
720 SCVF() {}
721
723 inline size_t from() const {return From;}
724
726 inline size_t to() const {return To;}
727
729 inline size_t num_ip() const {return nip;}
730
732 inline const MathVector<dim>& local_ip() const {return localIP;}
733
735 inline const MathVector<worldDim>& global_ip() const {return globalIP;}
736
738 inline const MathVector<worldDim>& normal() const {return Normal;}
739
741 inline const MathMatrix<worldDim,dim>& JTInv() const {return JtInv;}
742
744 inline number detJ() const {return detj;}
745
747 inline size_t num_sh() const {return numSH;}
748
750 inline number shape(size_t sh) const {return vShape[sh];}
751
753 inline const number* shape_vector() const {return vShape;}
754
756 inline const MathVector<dim>& local_grad(size_t sh) const
757 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vLocalGrad[sh];}
758
760 inline const MathVector<dim>* local_grad_vector() const {return vLocalGrad;}
761
763 inline const MathVector<worldDim>& global_grad(size_t sh) const
764 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vGlobalGrad[sh];}
765
768
770 inline size_t num_corners() const {return numCo;}
771
773 inline const MathVector<dim>& local_corner(size_t co) const
774 {UG_ASSERT(co < num_corners(), "Invalid index."); return vLocPos[co];}
775
777 inline const MathVector<worldDim>& global_corner(size_t co) const
778 {UG_ASSERT(co < num_corners(), "Invalid index."); return vGloPos[co];}
779
780 private:
781 // let outer class access private members
782 friend class DimFV1Geometry<dim, worldDim>;
783
784 // This scvf separates the scv with the ids given in "from" and "to"
785 // The computed normal points in direction from->to
786 size_t From, To;
787
788 // The normal on the SCVF pointing (from -> to)
789 MathVector<worldDim> Normal; // normal (incl. area)
790
791 // ordering is:
792 // 1D: edgeMidPoint
793 // 2D: edgeMidPoint, CenterOfElement
794 // 3D: edgeMidPoint, Side one, CenterOfElement, Side two
795 MathVector<dim> vLocPos[numCo]; // local corners of scvf
796 MathVector<worldDim> vGloPos[numCo]; // global corners of scvf
797 MidID vMidID[numCo]; // dimension and id of object, that's midpoint bounds the scvf
798
799 // scvf part
800 MathVector<dim> localIP; // local integration point
801 MathVector<worldDim> globalIP; // global integration point
802
803 // shapes and derivatives
804 size_t numSH;
805 number vShape[maxNSH]; // shapes at ip
806 MathVector<dim> vLocalGrad[maxNSH]; // local grad at ip
808 MathMatrix<worldDim,dim> JtInv; // Jacobian transposed at ip
809 number detj; // Jacobian det at ip
810 };
811
813 class SCV
814 {
815 public:
817 static const size_t numCo = traits::NumCornersOfSCV;
818
819 public:
820 SCV() {};
821
823 inline number volume() const {return Vol;}
824
826 inline size_t num_corners() const {return numCo;}
827
829 inline const MathVector<dim>& local_corner(size_t co) const
830 {UG_ASSERT(co < num_corners(), "Invalid index."); return vLocPos[co];}
831
833 inline const MathVector<worldDim>& global_corner(size_t co) const
834 {UG_ASSERT(co < num_corners(), "Invalid index."); return vGloPos[co];}
835
837 inline size_t node_id() const {return nodeId;}
838
840 inline size_t num_ip() const {return nip;}
841
843 inline const MathVector<dim>& local_ip() const {return vLocPos[0];}
844
846 inline const MathVector<worldDim>& global_ip() const {return vGloPos[0];}
847
849 inline const MathMatrix<worldDim,dim>& JTInv() const {return JtInv;}
850
852 inline number detJ() const {return detj;}
853
855 inline size_t num_sh() const {return numSH;}
856
858 inline number shape(size_t sh) const {return vShape[sh];}
859
861 inline const number* shape_vector() const {return vShape;}
862
864 inline const MathVector<dim>& local_grad(size_t sh) const
865 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vLocalGrad[sh];}
866
868 inline const MathVector<dim>* local_grad_vector() const {return vLocalGrad;}
869
871 inline const MathVector<worldDim>& global_grad(size_t sh) const
872 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vGlobalGrad[sh];}
873
876
877 private:
878 // let outer class access private members
879 friend class DimFV1Geometry<dim, worldDim>;
880
881 // node id of associated node
882 size_t nodeId;
883
884 // volume of scv
886
887 // local and global positions of this element
888 MathVector<dim> vLocPos[numCo]; // local position of node
889 MathVector<worldDim> vGloPos[numCo]; // global position of node
890 MidID vMidID[numCo]; // dimension and id of object, that's midpoint bounds the scv
891
892 // shapes and derivatives
893 size_t numSH;
894 number vShape[maxNSH]; // shapes at ip
895 MathVector<dim> vLocalGrad[maxNSH]; // local grad at ip
897 MathMatrix<worldDim,dim> JtInv; // Jacobian transposed at ip
898 number detj; // Jacobian det at ip
899 };
900
902 class BF
903 {
904 public:
906 static const size_t numCo = traits::NumCornersOfSCVF;
907
908 public:
909 BF() {}
910
912 inline size_t node_id() const {return nodeId;}
913
915 inline size_t num_ip() const {return nip;}
916
918 inline const MathVector<dim>& local_ip() const {return localIP;}
919
921 inline const MathVector<worldDim>& global_ip() const {return globalIP;}
922
924 inline const MathVector<worldDim>& normal() const {return Normal;} // includes area
925
927 inline number volume() const {return Vol;}
928
930 inline const MathMatrix<worldDim, dim>& JTInv() const {return JtInv;}
931
933 inline number detJ() const {return detj;}
934
936 inline size_t num_sh() const {return numSH;}
937
939 inline number shape(size_t sh) const
940 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vShape[sh];}
941
943 inline const number* shape_vector() const {return vShape;}
944
946 inline const MathVector<dim>& local_grad(size_t sh) const
947 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vLocalGrad[sh];}
948
950 inline const MathVector<dim>* local_grad_vector() const {return vLocalGrad;}
951
953 inline const MathVector<worldDim>& global_grad(size_t sh) const
954 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vGlobalGrad[sh];}
955
958
960 inline size_t num_corners() const {return numCo;}
961
963 inline const MathVector<dim>& local_corner(size_t co) const
964 {UG_ASSERT(co < num_corners(), "Invalid index."); return vLocPos[co];}
965
967 inline const MathVector<worldDim>& global_corner(size_t co) const
968 {UG_ASSERT(co < num_corners(), "Invalid index."); return vGloPos[co];}
969
970 private:
972 friend class DimFV1Geometry<dim, worldDim>;
973
974 // id of scv this bf belongs to
975 size_t nodeId;
976
977 // ordering is:
978 // 1D: edgeMidPoint
979 // 2D: edgeMidPoint, CenterOfElement
980 // 3D: edgeMidPoint, Side one, CenterOfElement, Side two
981 MathVector<dim> vLocPos[numCo]; // local corners of bf
982 MathVector<worldDim> vGloPos[numCo]; // global corners of bf
983 MidID vMidID[numCo]; // dimension and id of object, that's midpoint bounds the bf
984
985 // scvf part
986 MathVector<dim> localIP; // local integration point
987 MathVector<worldDim> globalIP; // global integration point
988 MathVector<worldDim> Normal; // normal (incl. area)
989 number Vol; // volume of bf
990
991 // shapes and derivatives
992 size_t numSH;
993 number vShape[maxNSH]; // shapes at ip
994 MathVector<dim> vLocalGrad[maxNSH]; // local grad at ip
996 MathMatrix<worldDim,dim> JtInv; // Jacobian transposed at ip
997 number detj; // Jacobian det at ip
998 };
999
1000 public:
1003
1005 void update(GridObject* elem, const MathVector<worldDim>* vCornerCoords,
1006 const ISubsetHandler* ish = NULL);
1007
1010 const MathVector<worldDim>* vCornerCoords,
1011 const ISubsetHandler* ish = NULL);
1012
1014 GridObject* elem() const {return m_pElem;}
1015
1017 const MathVector<worldDim>* corners() const {return m_vvGloMid[0];}
1018
1020 size_t num_scvf() const {return m_numSCVF;};
1021
1023 const SCVF& scvf(size_t i) const
1024 {UG_ASSERT(i < num_scvf(), "Invalid Index."); return m_vSCVF[i];}
1025
1027 // do not use this method to obtain the number of shape functions,
1028 // since this is NOT the same for pyramids; use num_sh() instead.
1029 size_t num_scv() const {return m_numSCV;}
1030
1032 const SCV& scv(size_t i) const
1033 {UG_ASSERT(i < num_scv(), "Invalid Index."); return m_vSCV[i];}
1034
1036 size_t num_sh() const {return m_nsh;};
1037
1038 public:
1040 size_t num_scvf_ips() const {return m_numSCVF;}
1041
1044
1047
1049 size_t num_scv_ips() const {return m_numSCV;}
1050
1054 return &(m_vLocSCV_IP[0]);
1055 else
1056 return &(m_vvLocMid[0][0]);
1057 }
1058
1062 return &(m_vGlobSCV_IP[0]);
1063 else
1064 return &(m_vvGloMid[0][0]);
1065 }
1066
1068 const MathVector<dim>& local_node_position(size_t nodeID) const
1069 {
1070 UG_ASSERT(nodeID < (size_t) maxMid, "Invalid node id.");
1071 return m_vvLocMid[0][nodeID];
1072 }
1073
1075 const MathVector<worldDim>& global_node_position(size_t nodeID) const
1076 {
1077 UG_ASSERT(nodeID < (size_t) maxMid, "Invalid node id.");
1078 return m_vvGloMid[0][nodeID];
1079 }
1080
1082 const MathVector<dim>* coe_local() const {return &(m_vvLocMid[dim][0]);}
1083
1085 const MathVector<worldDim>* coe_global() const {return &(m_vvGloMid[dim][0]);}
1086
1089
1092
1093 protected:
1094 // global and local ips on SCVF
1097
1098 // global and local ips on SCV (only needed for Pyramid and Octahedron)
1101
1102 public:
1104 inline void add_boundary_subset(int subsetIndex) {m_mapVectorBF[subsetIndex];}
1105
1107 inline void remove_boundary_subset(int subsetIndex) {m_mapVectorBF.erase(subsetIndex);}
1108
1110 inline void clear_boundary_subsets() {m_mapVectorBF.clear();}
1111
1113 inline size_t num_boundary_subsets() {return m_mapVectorBF.size();}
1114
1116 inline size_t num_bf() const
1117 {
1118 typename std::map<int, std::vector<BF> >::const_iterator it;
1119 size_t num = 0;
1120 for ( it=m_mapVectorBF.begin() ; it != m_mapVectorBF.end(); it++ )
1121 num += (*it).second.size();
1122 return num;
1123 }
1124
1126 inline size_t num_bf(int si) const
1127 {
1128 typename std::map<int, std::vector<BF> >::const_iterator it;
1129 it = m_mapVectorBF.find(si);
1130 if(it == m_mapVectorBF.end()) return 0;
1131 else return (*it).second.size();
1132 }
1133
1135 inline const BF& bf(int si, size_t i) const
1136 {
1137 typename std::map<int, std::vector<BF> >::const_iterator it;
1138 it = m_mapVectorBF.find(si);
1139 if(it == m_mapVectorBF.end()) UG_THROW("DimFV1Geom: No BndSubset "<<si);
1140 return (*it).second[i];
1141 }
1142
1144 inline const std::vector<BF>& bf(int si) const
1145 {
1146 typename std::map<int, std::vector<BF> >::const_iterator it;
1147 it = m_mapVectorBF.find(si);
1148 if(it == m_mapVectorBF.end()) return m_vEmptyVectorBF;
1149 return (*it).second;
1150 }
1151
1152 protected:
1153 std::map<int, std::vector<BF> > m_mapVectorBF;
1154 std::vector<BF> m_vEmptyVectorBF;
1155
1156 private:
1159
1162
1164 size_t m_numSCV;
1165
1168
1170 size_t m_nsh;
1171
1173 // (most objects in 1 dim, i.e. number of edges, but +1 for 1D)
1174 static const int maxMid = maxNumSCVF + 1;
1175
1179
1182
1185};
1186
1187
1189// FV1 Manifold Geometry
1191
1192template <typename TElem, int TWorldDim>
1194{
1195 public:
1196 // type of element
1197 typedef TElem elem_type;
1198
1199 // type of reference element
1201
1202 public:
1203 // order
1204 static const int order = 1;
1205
1206 // number of BoundaryFaces
1207 static const size_t m_numBF = ref_elem_type::numCorners;
1208
1209 // dimension of reference element
1210 static const int dim = ref_elem_type::dim;
1211
1212 // dimension of world
1213 static const int worldDim = TWorldDim;
1214
1215 // type of BoundaryFaces
1217
1218 // Hanging node flag: this Geometry does not support hanging nodes
1219 static const bool usesHangingNodes = false;
1220
1222 static const bool staticLocalData = true;
1223
1224 protected:
1225 struct MidID
1226 {
1227 MidID() : dim(0), id(0) {};
1228 MidID(size_t dim_, size_t id_) : dim(dim_), id(id_) {};
1229 size_t dim;
1230 size_t id;
1231 };
1232
1233 public:
1234 class BF
1235 {
1236 private:
1237 // let outer class access private members
1238 friend class FV1ManifoldGeometry<TElem, TWorldDim>;
1239
1240 // number of integration points
1241 static const size_t m_numIP = 1;
1242
1243 // max number of corners of bf
1245
1246 public:
1247 BF() {};
1248
1250 inline size_t node_id() const {return nodeId;}
1251
1253 inline size_t num_ip() const {return m_numIP;}
1254
1256 inline const MathVector<dim>& local_ip() const
1257 {return vLocPos[0];} // <-- always the vertex
1258
1260 inline const MathVector<worldDim>& global_ip() const
1261 {return vGloPos[0];} // <-- here too
1262
1264 inline number volume() const {return vol;}
1265
1267 inline size_t num_corners() const {return numCorners;}
1268
1270 inline const MathVector<dim>& local_corner(size_t i) const
1271 {UG_ASSERT(i < num_corners(), "Invalid index."); return vLocPos[i];}
1272
1274 inline const MathVector<worldDim>& global_corner(size_t i) const
1275 {UG_ASSERT(i < num_corners(), "Invalid index."); return vGloPos[i];}
1276
1278 inline size_t num_sh() const {return vShape.size();}
1279
1281 inline number shape(size_t i, size_t ip) const
1282 {UG_ASSERT(ip < num_ip(), "Invalid index"); return vShape[i];}
1283
1284 private:
1285 size_t nodeId; // id of associated node
1286
1287 // CORNERS: ordering is:
1288 // 1D: edgeMidPoint, CenterOfElement
1289 // 2D: edgeMidPoint, Side one, CenterOfElement, Side two
1290 MathVector<dim> vLocPos[numCorners]; // local position of node
1291 MathVector<worldDim> vGloPos[numCorners]; // global position of node
1292 MidID midId[numCorners]; // dimension and id of object, whose midpoint bounds the scv
1293
1294 // IPs & shapes
1295 std::vector<number> vShape; // shapes at ip
1296
1298 };
1299
1300 protected:
1302 {
1303 for (size_t i = 0; i < bf.num_corners(); ++i)
1304 {
1305 const size_t dim = bf.midId[i].dim;
1306 const size_t id = bf.midId[i].id;
1307 bf.vLocPos[i] = m_locMid[dim][id];
1308 }
1309 }
1310
1312 {
1313 for (size_t i = 0; i < bf.num_corners(); ++i)
1314 {
1315 const size_t dim = bf.midId[i].dim;
1316 const size_t id = bf.midId[i].id;
1317 bf.vGloPos[i] = m_gloMid[dim][id];
1318 }
1319 }
1320
1321 std::vector<MathVector<dim> > m_vLocBFIP;
1322 std::vector<MathVector<worldDim> > m_vGlobBFIP;
1323
1324
1325 public:
1328
1330 void update(GridObject* elem, const MathVector<worldDim>* vCornerCoords,
1331 const ISubsetHandler* ish = NULL);
1332
1334 GridObject* elem() const {return m_pElem;}
1335
1337 const MathVector<worldDim>* corners() const {return m_gloMid[0];}
1338
1340 inline size_t num_bf() const {return m_numBF;}
1341
1343 inline const BF& bf(size_t i) const
1344 {UG_ASSERT(i < num_bf(), "Invalid Index."); return m_vBF[i];}
1345
1348
1350 size_t num_bf_global_ips() const {return m_vGlobBFIP.size();}
1351
1353 const MathVector<dim>* bf_local_ips() const {return &m_vLocBFIP[0];}
1354
1356 size_t num_bf_local_ips() const {return m_vLocBFIP.size();}
1357
1358 private:
1359 // pointer to current element
1361
1362 // local and global geom object midpoints for each dimension
1365
1366 // BndFaces
1368
1369 // Reference Mapping
1371
1372 // Reference Element
1374};
1375
1376}
1377
1378#endif /* __H__UG__LIB_DISC__SPATIAL_DISC__DISC_HELPER__FV1_GEOMETRY__ */
boundary face
Definition fv1_geom.h:903
MathVector< dim > localIP
Definition fv1_geom.h:986
const MathMatrix< worldDim, dim > & JTInv() const
Transposed Inverse of Jacobian in integration point.
Definition fv1_geom.h:930
const MathVector< dim > & local_grad(size_t sh) const
value of local gradient of shape function i in integration point
Definition fv1_geom.h:946
const number * shape_vector() const
vector of local gradients in ip point
Definition fv1_geom.h:943
const MathVector< worldDim > & global_corner(size_t co) const
return global corner number i
Definition fv1_geom.h:967
const MathVector< dim > * local_grad_vector() const
vector of local gradients in ip point
Definition fv1_geom.h:950
MathVector< dim > vLocalGrad[maxNSH]
Definition fv1_geom.h:994
MidID vMidID[numCo]
Definition fv1_geom.h:983
MathVector< worldDim > Normal
Definition fv1_geom.h:988
number shape(size_t sh) const
value of shape function i in integration point
Definition fv1_geom.h:939
MathMatrix< worldDim, dim > JtInv
Definition fv1_geom.h:996
size_t node_id() const
index of SubControlVolume of the bf
Definition fv1_geom.h:912
size_t numSH
Definition fv1_geom.h:992
number vShape[maxNSH]
Definition fv1_geom.h:993
const MathVector< worldDim > & normal() const
outer normal on bf. Norm is equal to area
Definition fv1_geom.h:924
number Vol
Definition fv1_geom.h:989
MathVector< dim > vLocPos[numCo]
Definition fv1_geom.h:981
const MathVector< dim > & local_ip() const
local integration point of bf
Definition fv1_geom.h:918
MathVector< worldDim > vGloPos[numCo]
Definition fv1_geom.h:982
number detj
Definition fv1_geom.h:997
BF()
Definition fv1_geom.h:909
static const size_t numCo
Number of corners of bf.
Definition fv1_geom.h:906
size_t num_sh() const
number of shape functions
Definition fv1_geom.h:936
const MathVector< worldDim > & global_ip() const
global integration point of bf
Definition fv1_geom.h:921
const MathVector< dim > & local_corner(size_t co) const
return local corner number i
Definition fv1_geom.h:963
size_t num_ip() const
number of integration points on bf
Definition fv1_geom.h:915
MathVector< worldDim > vGlobalGrad[maxNSH]
Definition fv1_geom.h:995
MathVector< worldDim > globalIP
Definition fv1_geom.h:987
const MathVector< worldDim > * global_grad_vector() const
vector of global gradients in ip point
Definition fv1_geom.h:957
size_t nodeId
Definition fv1_geom.h:975
size_t num_corners() const
number of corners, that bound the scvf
Definition fv1_geom.h:960
const MathVector< worldDim > & global_grad(size_t sh) const
value of global gradient of shape function i in integration point
Definition fv1_geom.h:953
number volume() const
volume of bf
Definition fv1_geom.h:927
number detJ() const
Determinant of Jacobian in integration point.
Definition fv1_geom.h:933
Sub-Control Volume Face structure.
Definition fv1_geom.h:714
size_t From
Definition fv1_geom.h:786
size_t num_ip() const
number of integration points on scvf
Definition fv1_geom.h:729
MathVector< worldDim > globalIP
Definition fv1_geom.h:801
const MathVector< dim > & local_corner(size_t co) const
return local corner number i
Definition fv1_geom.h:773
MathVector< dim > localIP
Definition fv1_geom.h:800
MathVector< worldDim > Normal
Definition fv1_geom.h:789
const MathVector< dim > * local_grad_vector() const
vector of local gradients in ip point
Definition fv1_geom.h:760
number vShape[maxNSH]
Definition fv1_geom.h:805
const MathVector< worldDim > * global_grad_vector() const
vector of gloabl gradients in ip point
Definition fv1_geom.h:767
size_t to() const
index of SubControlVolume on one side of the scvf
Definition fv1_geom.h:726
size_t num_sh() const
number of shape functions
Definition fv1_geom.h:747
const MathVector< worldDim > & normal() const
normal on scvf (points direction "from"->"to"). Norm is equal to area
Definition fv1_geom.h:738
MathVector< dim > vLocPos[numCo]
Definition fv1_geom.h:795
MathVector< dim > vLocalGrad[maxNSH]
Definition fv1_geom.h:806
const MathMatrix< worldDim, dim > & JTInv() const
Transposed Inverse of Jacobian in integration point.
Definition fv1_geom.h:741
const number * shape_vector() const
vector of shape functions in ip point
Definition fv1_geom.h:753
size_t numSH
Definition fv1_geom.h:804
const MathVector< worldDim > & global_corner(size_t co) const
return glbal corner number i
Definition fv1_geom.h:777
MathVector< worldDim > vGloPos[numCo]
Definition fv1_geom.h:796
MathVector< worldDim > vGlobalGrad[maxNSH]
Definition fv1_geom.h:807
const MathVector< dim > & local_grad(size_t sh) const
value of local gradient of shape function i in integration point
Definition fv1_geom.h:756
number detj
Definition fv1_geom.h:809
size_t To
Definition fv1_geom.h:786
const MathVector< dim > & local_ip() const
local integration point of scvf
Definition fv1_geom.h:732
size_t from() const
index of SubControlVolume on one side of the scvf
Definition fv1_geom.h:723
static const size_t numCo
Number of corners of scvf.
Definition fv1_geom.h:717
number detJ() const
Determinant of Jacobian in integration point.
Definition fv1_geom.h:744
MidID vMidID[numCo]
Definition fv1_geom.h:797
MathMatrix< worldDim, dim > JtInv
Definition fv1_geom.h:808
const MathVector< worldDim > & global_ip() const
global integration point of scvf
Definition fv1_geom.h:735
number shape(size_t sh) const
value of shape function i in integration point
Definition fv1_geom.h:750
SCVF()
Definition fv1_geom.h:720
size_t num_corners() const
number of corners, that bound the scvf
Definition fv1_geom.h:770
const MathVector< worldDim > & global_grad(size_t sh) const
value of global gradient of shape function i in integration point
Definition fv1_geom.h:763
sub control volume structure
Definition fv1_geom.h:814
size_t num_ip() const
number of integration points
Definition fv1_geom.h:840
size_t num_sh() const
number of shape functions
Definition fv1_geom.h:855
const MathVector< dim > & local_grad(size_t sh) const
value of local gradient of shape function i in integration point
Definition fv1_geom.h:864
MathVector< worldDim > vGloPos[numCo]
Definition fv1_geom.h:889
static const size_t numCo
Number of corners of scv.
Definition fv1_geom.h:817
const MathVector< dim > & local_corner(size_t co) const
return local corner number i
Definition fv1_geom.h:829
const number * shape_vector() const
vector of shape functions in ip point
Definition fv1_geom.h:861
number volume() const
volume of scv
Definition fv1_geom.h:823
const MathVector< worldDim > & global_grad(size_t sh) const
value of global gradient of shape function i in integration point
Definition fv1_geom.h:871
MathVector< worldDim > vGlobalGrad[maxNSH]
Definition fv1_geom.h:896
number detj
Definition fv1_geom.h:898
const MathMatrix< worldDim, dim > & JTInv() const
Transposed Inverse of Jacobian in integration point.
Definition fv1_geom.h:849
size_t numSH
Definition fv1_geom.h:893
const MathVector< worldDim > & global_ip() const
global integration point
Definition fv1_geom.h:846
number detJ() const
Determinant of Jacobian in integration point.
Definition fv1_geom.h:852
number shape(size_t sh) const
value of shape function i in integration point
Definition fv1_geom.h:858
MidID vMidID[numCo]
Definition fv1_geom.h:890
size_t nodeId
Definition fv1_geom.h:882
const MathVector< dim > * local_grad_vector() const
vector of local gradients in ip point
Definition fv1_geom.h:868
MathMatrix< worldDim, dim > JtInv
Definition fv1_geom.h:897
MathVector< dim > vLocPos[numCo]
Definition fv1_geom.h:888
number vShape[maxNSH]
Definition fv1_geom.h:894
SCV()
Definition fv1_geom.h:820
const MathVector< worldDim > & global_corner(size_t co) const
return global corner number i
Definition fv1_geom.h:833
const MathVector< worldDim > * global_grad_vector() const
vector of global gradients in ip point
Definition fv1_geom.h:875
size_t node_id() const
node id that this scv is associated to
Definition fv1_geom.h:837
number Vol
Definition fv1_geom.h:885
size_t num_corners() const
number of corners, that bound the scv
Definition fv1_geom.h:826
const MathVector< dim > & local_ip() const
local integration point of scv
Definition fv1_geom.h:843
MathVector< dim > vLocalGrad[maxNSH]
Definition fv1_geom.h:895
Geometry and shape functions for 1st order Vertex-Centered Finite Volume.
Definition fv1_geom.h:665
GridObject * m_pElem
pointer to current element
Definition fv1_geom.h:1158
MathVector< dim > m_vLocSCV_IP[maxNumSCV]
Definition fv1_geom.h:1100
size_t num_scv() const
number of SubControlVolumes
Definition fv1_geom.h:1029
static const int order
order
Definition fv1_geom.h:685
static const size_t maxNSH
max number of shape functions
Definition fv1_geom.h:697
static const size_t maxNumSCVF
max number of SubControlVolumeFaces
Definition fv1_geom.h:694
size_t m_numSCV
current number of scv
Definition fv1_geom.h:1164
const MathVector< worldDim > * scvf_global_ips() const
returns all ips of scvf as they appear in scv loop
Definition fv1_geom.h:1046
static const int dim
dimension of reference element
Definition fv1_geom.h:672
const MathVector< worldDim > * coe_global() const
returns the global coordinates of the center of mass of the element
Definition fv1_geom.h:1085
const MathVector< worldDim > * scv_global_ips() const
returns all ips of scv as they appear in scv loop
Definition fv1_geom.h:1060
MathVector< dim > m_vLocSCVF_IP[maxNumSCVF]
Definition fv1_geom.h:1096
const SCVF & scvf(size_t i) const
const access to SubControlVolumeFace number i
Definition fv1_geom.h:1023
static const bool usesHangingNodes
Hanging node flag: this Geometry does not support hanging nodes.
Definition fv1_geom.h:678
MathVector< worldDim > m_vGlobSCVF_IP[maxNumSCVF]
Definition fv1_geom.h:1095
MathVector< dim > m_vvLocMid[dim+1][maxMid]
local and global geom object midpoints for each dimension
Definition fv1_geom.h:1177
static const int maxMid
max number of geometric objects in a dimension
Definition fv1_geom.h:1174
const MathVector< dim > * coe_local() const
returns the local coordinates of the center of mass of the element
Definition fv1_geom.h:1082
size_t num_sh() const
number of shape functions
Definition fv1_geom.h:1036
const BF & bf(int si, size_t i) const
returns the boundary face i for subsetIndex
Definition fv1_geom.h:1135
const MathVector< worldDim > * corners() const
get vector of corners for current element
Definition fv1_geom.h:1017
const MathVector< dim > * scv_local_ips() const
returns all ips of scv as they appear in scv loop
Definition fv1_geom.h:1052
MathVector< worldDim > m_vGlobSCV_IP[maxNumSCV]
Definition fv1_geom.h:1099
static const int worldDim
dimension of world
Definition fv1_geom.h:675
fv1_dim_traits< TDim, TWorldDim > traits
used traits
Definition fv1_geom.h:668
const MathVector< dim > * scvf_local_ips() const
returns all ips of scvf as they appear in scv loop
Definition fv1_geom.h:1043
size_t num_bf() const
number of all boundary faces
Definition fv1_geom.h:1116
void update(GridObject *elem, const MathVector< worldDim > *vCornerCoords, const ISubsetHandler *ish=NULL)
update data for given element
Definition fv1_geom.cpp:669
ReferenceObjectID roid() const
returns reference object id
Definition fv1_geom.h:1088
DimFV1Geometry()
construct object and initialize local values and sizes
Definition fv1_geom.h:1002
size_t num_bf(int si) const
number of boundary faces on subset 'subsetIndex'
Definition fv1_geom.h:1126
static const size_t nip
number of integration points
Definition fv1_geom.h:700
SCV m_vSCV[maxNumSCV]
SubControlVolumes.
Definition fv1_geom.h:1184
const std::vector< BF > & bf(int si) const
returns reference to vector of boundary faces for subsetIndex
Definition fv1_geom.h:1144
const MathVector< dim > & local_node_position(size_t nodeID) const
return local coords for node ID
Definition fv1_geom.h:1068
SCVF m_vSCVF[maxNumSCVF]
SubControlVolumeFaces.
Definition fv1_geom.h:1181
std::vector< BF > m_vEmptyVectorBF
Definition fv1_geom.h:1154
GridObject * elem() const
get the element
Definition fv1_geom.h:1014
void add_boundary_subset(int subsetIndex)
add subset that is interpreted as boundary subset.
Definition fv1_geom.h:1104
traits::scv_type scv_type
type of SubControlVolume
Definition fv1_geom.h:691
void update_boundary_faces(GridObject *elem, const MathVector< worldDim > *vCornerCoords, const ISubsetHandler *ish=NULL)
update boundary data for given element
Definition fv1_geom.cpp:779
size_t m_nsh
current number of shape functions
Definition fv1_geom.h:1170
size_t num_boundary_subsets()
number of registered boundary subsets
Definition fv1_geom.h:1113
void update_local(ReferenceObjectID roid)
update local data
Definition fv1_geom.cpp:578
static const bool staticLocalData
flag indicating if local data may change
Definition fv1_geom.h:681
MathVector< worldDim > m_vvGloMid[dim+1][maxMid]
Definition fv1_geom.h:1178
size_t num_scvf_ips() const
returns number of all scvf ips
Definition fv1_geom.h:1040
void remove_boundary_subset(int subsetIndex)
removes subset that is interpreted as boundary subset.
Definition fv1_geom.h:1107
ReferenceObjectID m_roid
current reference object id
Definition fv1_geom.h:1161
size_t num_scv_ips() const
returns number of all scv ips
Definition fv1_geom.h:1049
void clear_boundary_subsets()
reset all boundary subsets
Definition fv1_geom.h:1110
const SCV & scv(size_t i) const
const access to SubControlVolume number i
Definition fv1_geom.h:1032
const MathVector< worldDim > & global_node_position(size_t nodeID) const
return global coords for node ID
Definition fv1_geom.h:1075
size_t m_numSCVF
current number of scvf
Definition fv1_geom.h:1167
static const size_t maxNumSCV
number of SubControlVolumes
Definition fv1_geom.h:688
size_t num_scvf() const
number of SubControlVolumeFaces
Definition fv1_geom.h:1020
std::map< int, std::vector< BF > > m_mapVectorBF
Definition fv1_geom.h:1153
Geometry and shape functions for 1st order Vertex-Centered Finite Volume.
Definition fv1_geom.h:95
boundary face
Definition fv1_geom.h:370
MathVector< worldDim > globalIP
Definition fv1_geom.h:454
number vShape[nsh]
Definition fv1_geom.h:459
BF()
Definition fv1_geom.h:376
MathMatrix< worldDim, dim > JtInv
Definition fv1_geom.h:462
MathVector< dim > vLocalGrad[nsh]
Definition fv1_geom.h:460
number shape(size_t sh) const
value of shape function i in integration point
Definition fv1_geom.h:406
const number * shape_vector() const
vector of local gradients in ip point
Definition fv1_geom.h:410
const MathVector< worldDim > * global_grad_vector() const
vector of global gradients in ip point
Definition fv1_geom.h:424
const MathVector< worldDim > & global_corner(size_t co) const
return global corner number i
Definition fv1_geom.h:434
const MathVector< dim > & local_ip() const
local integration point of bf
Definition fv1_geom.h:385
static const size_t numCo
Number of corners of bf.
Definition fv1_geom.h:373
MathVector< worldDim > Normal
Definition fv1_geom.h:455
number detj
Definition fv1_geom.h:463
size_t num_ip() const
number of integration points on bf
Definition fv1_geom.h:382
const MathVector< worldDim > & global_grad(size_t sh) const
value of global gradient of shape function i in integration point
Definition fv1_geom.h:420
const MathVector< worldDim > & normal() const
outer normal on bf. Norm is equal to area
Definition fv1_geom.h:391
number volume() const
volume of bf
Definition fv1_geom.h:394
MathVector< worldDim > vGloPos[numCo]
Definition fv1_geom.h:449
MidID vMidID[numCo]
Definition fv1_geom.h:450
size_t node_id() const
index of SubControlVolume of the bf
Definition fv1_geom.h:379
size_t nodeId
Definition fv1_geom.h:442
const MathMatrix< worldDim, dim > & JTInv() const
Transposed Inverse of Jacobian in integration point.
Definition fv1_geom.h:397
const MathVector< dim > & local_corner(size_t co) const
return local corner number i
Definition fv1_geom.h:430
MathVector< dim > localIP
Definition fv1_geom.h:453
const MathVector< worldDim > & global_ip() const
global integration point of bf
Definition fv1_geom.h:388
size_t num_corners() const
number of corners, that bound the scvf
Definition fv1_geom.h:427
MathVector< dim > vLocPos[numCo]
Definition fv1_geom.h:448
number Vol
Definition fv1_geom.h:456
const MathVector< dim > * local_grad_vector() const
vector of local gradients in ip point
Definition fv1_geom.h:417
number detJ() const
Determinant of Jacobian in integration point.
Definition fv1_geom.h:400
MathVector< worldDim > vGlobalGrad[nsh]
Definition fv1_geom.h:461
const MathVector< dim > & local_grad(size_t sh) const
value of local gradient of shape function i in integration point
Definition fv1_geom.h:413
size_t num_sh() const
number of shape functions
Definition fv1_geom.h:403
Sub-Control Volume Face structure.
Definition fv1_geom.h:175
MathVector< dim > localIP
Definition fv1_geom.h:261
const number * shape_vector() const
vector of shape functions in ip point
Definition fv1_geom.h:214
MathVector< dim > vLocalGrad[nsh]
Definition fv1_geom.h:266
number detJ() const
Determinant of Jacobian in integration point.
Definition fv1_geom.h:205
SCVF()
Definition fv1_geom.h:181
const MathMatrix< worldDim, dim > & JTInv() const
Transposed Inverse of Jacobian in integration point.
Definition fv1_geom.h:202
size_t From
Definition fv1_geom.h:247
const MathVector< worldDim > & global_grad(size_t sh) const
value of global gradient of shape function i in integration point
Definition fv1_geom.h:224
number shape(size_t sh) const
value of shape function i in integration point
Definition fv1_geom.h:211
number vShape[nsh]
Definition fv1_geom.h:265
size_t num_sh() const
number of shape functions
Definition fv1_geom.h:208
static const size_t numCo
Number of corners of scvf.
Definition fv1_geom.h:178
MathVector< worldDim > vGlobalGrad[nsh]
Definition fv1_geom.h:267
const MathVector< dim > & local_corner(size_t co) const
return local corner number i
Definition fv1_geom.h:234
MathVector< worldDim > globalIP
Definition fv1_geom.h:262
MathVector< dim > vLocPos[numCo]
Definition fv1_geom.h:256
MathVector< worldDim > vGloPos[numCo]
Definition fv1_geom.h:257
MathVector< worldDim > Normal
Definition fv1_geom.h:250
MathMatrix< worldDim, dim > JtInv
Definition fv1_geom.h:268
number detj
Definition fv1_geom.h:269
const MathVector< worldDim > & normal() const
normal on scvf (points direction "from"->"to"). Norm is equal to area
Definition fv1_geom.h:190
const MathVector< worldDim > & global_ip() const
global integration point of scvf
Definition fv1_geom.h:199
size_t num_corners() const
number of corners, that bound the scvf
Definition fv1_geom.h:231
MidID vMidID[numCo]
Definition fv1_geom.h:258
const MathVector< dim > & local_grad(size_t sh) const
value of local gradient of shape function i in integration point
Definition fv1_geom.h:217
const MathVector< worldDim > * global_grad_vector() const
vector of global gradients in ip point
Definition fv1_geom.h:228
const MathVector< dim > * local_grad_vector() const
vector of local gradients in ip point
Definition fv1_geom.h:221
size_t num_ip() const
number of integration points on scvf
Definition fv1_geom.h:193
size_t To
Definition fv1_geom.h:247
size_t from() const
index of SubControlVolume on one side of the scvf
Definition fv1_geom.h:184
const MathVector< worldDim > & global_corner(size_t co) const
return global corner number i
Definition fv1_geom.h:238
size_t to() const
index of SubControlVolume on one side of the scvf
Definition fv1_geom.h:187
const MathVector< dim > & local_ip() const
local integration point of scvf
Definition fv1_geom.h:196
sub control volume structure
Definition fv1_geom.h:274
MidID midId[numCo]
Definition fv1_geom.h:358
const MathVector< dim > * local_grad_vector() const
vector of local gradients in ip point
Definition fv1_geom.h:336
MathVector< worldDim > vGloPos[numCo]
Definition fv1_geom.h:357
const MathVector< dim > & local_corner(size_t co) const
return local corner number i
Definition fv1_geom.h:289
MathMatrix< worldDim, dim > JtInv
Definition fv1_geom.h:364
const MathVector< worldDim > * global_corners() const
return global corners
Definition fv1_geom.h:301
number shape(size_t sh) const
value of shape function i in integration point
Definition fv1_geom.h:326
SCV()
Definition fv1_geom.h:280
const MathVector< dim > & local_ip() const
local integration point of scv
Definition fv1_geom.h:311
MathVector< dim > vLocPos[numCo]
Definition fv1_geom.h:356
static const size_t numCo
Number of corners of scvf.
Definition fv1_geom.h:277
size_t num_corners() const
number of corners, that bound the scvf
Definition fv1_geom.h:286
const MathVector< worldDim > & global_ip() const
global integration point
Definition fv1_geom.h:314
const MathVector< dim > & local_grad(size_t sh) const
value of local gradient of shape function i in integration point
Definition fv1_geom.h:332
MathVector< dim > vLocalGrad[nsh]
Definition fv1_geom.h:362
const number * shape_vector() const
vector of shape functions in ip point
Definition fv1_geom.h:329
size_t nodeId
Definition fv1_geom.h:350
const MathVector< worldDim > & global_grad(size_t sh) const
value of global gradient of shape function i in integration point
Definition fv1_geom.h:339
const MathVector< worldDim > * global_grad_vector() const
vector of global gradients in ip point
Definition fv1_geom.h:343
number detj
Definition fv1_geom.h:365
size_t num_ip() const
number of integration points
Definition fv1_geom.h:308
const MathVector< worldDim > & global_corner(size_t co) const
return global corner number i
Definition fv1_geom.h:293
number detJ() const
Determinant of Jacobian in integration point.
Definition fv1_geom.h:320
MathVector< worldDim > vGlobalGrad[nsh]
Definition fv1_geom.h:363
number vShape[nsh]
Definition fv1_geom.h:361
number Vol
Definition fv1_geom.h:353
const MathMatrix< worldDim, dim > & JTInv() const
Transposed Inverse of Jacobian in integration point.
Definition fv1_geom.h:317
size_t num_sh() const
number of shape functions
Definition fv1_geom.h:323
number volume() const
volume of scv
Definition fv1_geom.h:283
const MathVector< dim > * local_corners() const
return local corners
Definition fv1_geom.h:297
size_t node_id() const
node id that this scv is associated to
Definition fv1_geom.h:305
Geometry and shape functions for 1st order Vertex-Centered Finite Volume.
Definition fv1_geom.h:114
TElem * elem() const
get the element
Definition fv1_geom.h:483
size_t num_boundary_subsets()
number of registered boundary subsets
Definition fv1_geom.h:580
const MathVector< worldDim > * corners() const
get vector of the global coordinates of corners for current element
Definition fv1_geom.h:486
static const size_t numSCV
number of SubControlVolumes
Definition fv1_geom.h:146
TElem * m_pElem
pointer to current element
Definition fv1_geom.h:628
void reset_curr_elem()
Definition fv1_geom.h:620
const BF & bf(int si, size_t i) const
returns the boundary face i for subsetIndex
Definition fv1_geom.h:602
const local_shape_fct_set_type & m_rTrialSpace
Shape function set.
Definition fv1_geom.h:651
FV1Geometry_gen()
construct object and initialize local values and sizes
Definition fv1_geom.cpp:251
const MathVector< dim > * scv_local_ips() const
returns all ips of scv as they appear in scv loop
Definition fv1_geom.h:525
const std::vector< BF > & bf(int si) const
returns reference to vector of boundary faces for subsetIndex
Definition fv1_geom.h:612
ReferenceMapping< ref_elem_type, worldDim > m_mapping
Reference Mapping.
Definition fv1_geom.h:645
void update_boundary_faces(GridObject *elem, const MathVector< worldDim > *vCornerCoords, const ISubsetHandler *ish=NULL)
update boundary data for given element
Definition fv1_geom.cpp:468
size_t num_bf() const
number of all boundary faces
Definition fv1_geom.h:583
LagrangeP1< ref_elem_type > local_shape_fct_set_type
type of Shape function used
Definition fv1_geom.h:155
void remove_boundary_subset(int subsetIndex)
removes subset that is interpreted as boundary subset.
Definition fv1_geom.h:574
MathVector< dim > m_vLocSCV_IP[numSCV]
Definition fv1_geom.h:567
MathVector< worldDim > m_vGlobSCVF_IP[numSCVF]
Definition fv1_geom.h:562
const MathVector< worldDim > & global_node_position(size_t nodeID) const
return global coords for node ID
Definition fv1_geom.h:548
static const bool staticLocalData
flag indicating if local data may change
Definition fv1_geom.h:136
static const int order
order
Definition fv1_geom.h:143
static const int worldDim
dimension of world
Definition fv1_geom.h:130
SCV m_vSCV[numSCV]
SubControlVolumes.
Definition fv1_geom.h:642
void clear_boundary_subsets()
reset all boundary subsets
Definition fv1_geom.h:577
size_t num_bf(int si) const
number of boundary faces on subset 'subsetIndex'
Definition fv1_geom.h:593
static const int dim
dimension of reference element
Definition fv1_geom.h:127
size_t num_scv() const
number of SubControlVolumes
Definition fv1_geom.h:498
static const bool condensed_scvf_ips
whether the scheme shifts the scvf ip's to midpoints of the edges
Definition fv1_geom.h:139
MathVector< worldDim > m_vGlobSCV_IP[numSCV]
Definition fv1_geom.h:566
fv1_traits< ref_elem_type, TWorldDim > traits
used traits
Definition fv1_geom.h:123
size_t num_scvf_ips() const
returns number of all scvf ips
Definition fv1_geom.h:513
ReferenceObjectID roid() const
returns reference object id
Definition fv1_geom.h:508
static const int maxMid
max number of geom objects in all dimensions
Definition fv1_geom.h:632
const SCVF & scvf(size_t i) const
const access to SubControlVolumeFace number i
Definition fv1_geom.h:492
reference_element_traits< TElem >::reference_element_type ref_elem_type
type of reference element
Definition fv1_geom.h:120
std::vector< BF > m_vEmptyVectorBF
Definition fv1_geom.h:624
const MathVector< dim > & local_node_position(size_t nodeID) const
return local coords for node ID
Definition fv1_geom.h:541
const MathVector< worldDim > * scvf_global_ips() const
returns all ips of scvf as they appear in scv loop
Definition fv1_geom.h:519
size_t num_scvf() const
number of SubControlVolumeFaces
Definition fv1_geom.h:489
const MathVector< worldDim > * scv_global_ips() const
returns all ips of scv as they appear in scv loop
Definition fv1_geom.h:533
static const bool usesHangingNodes
Hanging node flag: this Geometry does not support hanging nodes.
Definition fv1_geom.h:133
size_t num_scv_ips() const
returns number of all scv ips
Definition fv1_geom.h:522
static const size_t numSCVF
number of SubControlVolumeFaces
Definition fv1_geom.h:152
SCVF m_vSCVF[numSCVF]
SubControlVolumeFaces.
Definition fv1_geom.h:639
TElem elem_type
type of element
Definition fv1_geom.h:117
void add_boundary_subset(int subsetIndex)
add subset that is interpreted as boundary subset.
Definition fv1_geom.h:571
MathVector< dim > m_vLocSCVF_IP[numSCVF]
Definition fv1_geom.h:563
const MathVector< dim > * scvf_local_ips() const
returns all ips of scvf as they appear in scv loop
Definition fv1_geom.h:516
void update_local_data()
update local data
Definition fv1_geom.cpp:260
MathVector< worldDim > m_vvGloMid[dim+1][maxMid]
Definition fv1_geom.h:636
const MathVector< worldDim > * coe_global() const
returns the global coordinates of the center of mass of the element
Definition fv1_geom.h:558
void update(GridObject *elem, const MathVector< worldDim > *vCornerCoords, const ISubsetHandler *ish=NULL)
update data for given element
Definition fv1_geom.cpp:331
traits::scv_type scv_type
type of SubControlVolume
Definition fv1_geom.h:149
size_t num_sh() const
number of shape functions
Definition fv1_geom.h:505
static const size_t nsh
number of shape functions
Definition fv1_geom.h:158
static const size_t nip
number of integration points
Definition fv1_geom.h:161
std::map< int, std::vector< BF > > m_mapVectorBF
Definition fv1_geom.h:623
const ref_elem_type & m_rRefElem
Reference Element.
Definition fv1_geom.h:648
MathVector< dim > m_vvLocMid[dim+1][maxMid]
local and global geom object midpoints for each dimension
Definition fv1_geom.h:635
const MathVector< dim > * coe_local() const
returns the local coordinates of the center of mass of the element
Definition fv1_geom.h:555
const SCV & scv(size_t i) const
const access to SubControlVolume number i
Definition fv1_geom.h:501
Geometry and shape functions for 1st order Vertex-Centered Finite Volume.
Definition fv1_geom.h:76
Definition fv1_geom.h:1235
MidID midId[numCorners]
Definition fv1_geom.h:1292
const MathVector< worldDim > & global_corner(size_t i) const
return global position of corner number i
Definition fv1_geom.h:1274
MathVector< worldDim > vGloPos[numCorners]
Definition fv1_geom.h:1291
BF()
Definition fv1_geom.h:1247
size_t num_ip() const
number of integration points
Definition fv1_geom.h:1253
size_t num_sh() const
number of shape functions
Definition fv1_geom.h:1278
size_t nodeId
Definition fv1_geom.h:1285
const MathVector< dim > & local_corner(size_t i) const
return local position of corner number i
Definition fv1_geom.h:1270
size_t num_corners() const
number of corners, that bound the bf
Definition fv1_geom.h:1267
const MathVector< worldDim > & global_ip() const
global integration point
Definition fv1_geom.h:1260
number shape(size_t i, size_t ip) const
value of shape function i in integration point
Definition fv1_geom.h:1281
number volume() const
volume of bf
Definition fv1_geom.h:1264
const MathVector< dim > & local_ip() const
local integration point of bf
Definition fv1_geom.h:1256
static const size_t numCorners
Definition fv1_geom.h:1244
size_t node_id() const
node id that this bf is associated to
Definition fv1_geom.h:1250
number vol
Definition fv1_geom.h:1297
static const size_t m_numIP
Definition fv1_geom.h:1241
MathVector< dim > vLocPos[numCorners]
Definition fv1_geom.h:1290
std::vector< number > vShape
Definition fv1_geom.h:1295
Definition fv1_geom.h:1194
const MathVector< worldDim > * bf_global_ips() const
returns all ips of scvf as they appear in scv loop
Definition fv1_geom.h:1347
static const bool usesHangingNodes
Definition fv1_geom.h:1219
static const int dim
Definition fv1_geom.h:1210
static const size_t m_numBF
Definition fv1_geom.h:1207
void copy_local_corners(BF &bf)
Definition fv1_geom.h:1301
size_t num_bf() const
number of BoundaryFaces
Definition fv1_geom.h:1340
ReferenceMapping< ref_elem_type, worldDim > m_rMapping
Definition fv1_geom.h:1370
static const bool staticLocalData
flag indicating if local data may change
Definition fv1_geom.h:1222
static const int worldDim
Definition fv1_geom.h:1213
size_t num_bf_local_ips() const
returns number of all BF ips
Definition fv1_geom.h:1356
MathVector< dim > m_locMid[dim+1][m_numBF]
Definition fv1_geom.h:1363
GridObject * m_pElem
Definition fv1_geom.h:1360
const BF & bf(size_t i) const
const access to Boundary Face number i
Definition fv1_geom.h:1343
void update(GridObject *elem, const MathVector< worldDim > *vCornerCoords, const ISubsetHandler *ish=NULL)
update data for given element
Definition fv1_geom.cpp:999
static const int order
Definition fv1_geom.h:1204
BF m_vBF[m_numBF]
Definition fv1_geom.h:1367
TElem elem_type
Definition fv1_geom.h:1197
const ref_elem_type & m_rRefElem
Definition fv1_geom.h:1373
fv1_traits< ref_elem_type, dim >::scv_type bf_type
Definition fv1_geom.h:1216
reference_element_traits< TElem >::reference_element_type ref_elem_type
Definition fv1_geom.h:1200
std::vector< MathVector< dim > > m_vLocBFIP
Definition fv1_geom.h:1321
void copy_global_corners(BF &bf)
Definition fv1_geom.h:1311
FV1ManifoldGeometry()
constructor
Definition fv1_geom.cpp:903
size_t num_bf_global_ips() const
returns number of all BF ips
Definition fv1_geom.h:1350
GridObject * elem() const
get the element
Definition fv1_geom.h:1334
std::vector< MathVector< worldDim > > m_vGlobBFIP
Definition fv1_geom.h:1322
const MathVector< worldDim > * corners() const
get vector of corners for current element
Definition fv1_geom.h:1337
MathVector< worldDim > m_gloMid[dim+1][m_numBF]
Definition fv1_geom.h:1364
const MathVector< dim > * bf_local_ips() const
returns all ips of BF as they appear in scv loop
Definition fv1_geom.h:1353
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 subset_handler_interface.h:223
Lagrange shape functions of first order.
Definition lagrangep1.h:55
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_mapping.h:65
#define UG_ASSERT(expr, msg)
Definition assert.h:70
#define UG_THROW(msg)
Definition error.h:57
double number
Definition types.h:124
the ug namespace
ReferenceObjectID
these ids are used to identify the shape of a geometric object.
Definition grid_base_objects.h:74
@ ROID_PYRAMID
Definition grid_base_objects.h:83
@ ROID_UNKNOWN
Definition grid_base_objects.h:75
@ ROID_OCTAHEDRON
Definition grid_base_objects.h:84
Definition fv1_geom.h:1226
size_t dim
Definition fv1_geom.h:1229
size_t id
Definition fv1_geom.h:1230
MidID(size_t dim_, size_t id_)
Definition fv1_geom.h:1228
MidID()
Definition fv1_geom.h:1227
helper class to store dimension and id of a midpoint of a sub-element
Definition fv_geom_base.h:41
Traits for Finite Volumes for a generic element of the fixed dimensionalities.
Definition fv_util.h:1325
Traits for Finite Volumes (dummy implementation, s. the instantiations below)
Definition fv_util.h:178
void scv_type
Definition fv_util.h:204
traits for reference elements
Definition reference_element_traits.h:48