Loading [MathJax]/extensions/tex2jax.js
ug4
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
fv1ib_geom.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2014-2015: G-CSC, Goethe University Frankfurt
3 * Author: Susanne Höllbacher
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 FV1IB_GEOM_H_
34#define FV1IB_GEOM_H_
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
63
67template < typename TElem, int TWorldDim>
69{
70public:
72 typedef TElem elem_type;
73
76
79
80public:
82 static const int dim = ref_elem_type::dim;
83
85 static const int worldDim = TWorldDim;
86
88 static const bool usesHangingNodes = false;
89
91 static const bool staticLocalData = true;
92
93public:
95 static const int order = 1;
96
98 static const size_t numSCV = (ref_elem_type::REFERENCE_OBJECT_ID != ROID_PYRAMID)
99 ? ref_elem_type::numCorners : 8;
100
102 typedef typename traits::scv_type scv_type;
103
105 static const size_t numSCVF = (ref_elem_type::REFERENCE_OBJECT_ID != ROID_PYRAMID)
106 ? ref_elem_type::numEdges : 12;
107
110
112 static const size_t nsh = local_shape_fct_set_type::nsh;
113
115 static const size_t nip = 1;
116
117public:
119
128 class SCVF
129 {
130 public:
132 static const size_t numCo = traits::NumCornersOfSCVF;
133
134 public:
135 SCVF() {}
136
138 inline size_t from() const {return From;}
139
141 inline size_t to() const {return To;}
142
144 inline const MathVector<worldDim>& normal() const {return Normal;}
145
147 inline size_t num_ip() const {return nip;}
148
150 inline const MathVector<dim>& local_ip() const {return localIP;}
151
153 inline const MathVector<worldDim>& global_ip() const {return globalIP;}
154
156 inline const MathMatrix<worldDim,dim>& JTInv() const {return JtInv;}
157
159 inline number detJ() const {return detj;}
160
162 inline size_t num_sh() const {return nsh;}
163
165 inline number shape(size_t sh) const {return vShape[sh];}
166
168 inline const number* shape_vector() const {return vShape;}
169
171 inline const MathVector<dim>& local_grad(size_t sh) const
172 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vLocalGrad[sh];}
173
175 inline const MathVector<dim>* local_grad_vector() const {return vLocalGrad;}
176
178 inline const MathVector<worldDim>& global_grad(size_t sh) const
179 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vGlobalGrad[sh];}
180
183
185 inline size_t num_corners() const {return numCo;}
186
188 inline const MathVector<dim>& local_corner(size_t co) const
189 {UG_ASSERT(co < num_corners(), "Invalid index."); return vLocPos[co];}
190
192 inline const MathVector<worldDim>& global_corner(size_t co) const
193 {UG_ASSERT(co < num_corners(), "Invalid index."); return vGloPos[co];}
194
195 protected:
196 // let outer class access private members
197 friend class FV1IBGeometry<TElem, TWorldDim>;
198
199 // This scvf separates the scv with the ids given in "from" and "to"
200 // The computed normal points in direction from->to
201 size_t From, To;
202
203 // The normal on the SCVF pointing (from -> to)
204 MathVector<worldDim> Normal; // normal (incl. area)
205
206 // ordering is:
207 // 1D: edgeMidPoint
208 // 2D: edgeMidPoint, CenterOfElement
209 // 3D: edgeMidPoint, Side one, CenterOfElement, Side two
210 MathVector<dim> vLocPos[numCo]; // local corners of scvf
211 MathVector<worldDim> vGloPos[numCo]; // global corners of scvf
212 MidID vMidID[numCo]; // dimension and id of object, that's midpoint bounds the scvf
213
214 // scvf part
215 MathVector<dim> localIP; // local integration point
216 MathVector<worldDim> globalIP; // global integration point
217
218 // shapes and derivatives
219 number vShape[nsh]; // shapes at ip
220 MathVector<dim> vLocalGrad[nsh]; // local grad at ip
221 MathVector<worldDim> vGlobalGrad[nsh]; // global grad at ip
222 MathMatrix<worldDim,dim> JtInv; // Jacobian transposed at ip
223 number detj; // Jacobian det at ip
224 };
225
227 class SCV
228 {
229 public:
231 static const size_t numCo = traits::NumCornersOfSCV;
232
233 public:
234 SCV() {};
235
237 inline number volume() const {return Vol;}
238
240 inline size_t num_corners() const {return numCo;}
241
243 inline const MathVector<dim>& local_corner(size_t co) const
244 {UG_ASSERT(co < num_corners(), "Invalid index."); return vLocPos[co];}
245
247 inline const MathVector<worldDim>& global_corner(size_t co) const
248 {UG_ASSERT(co < num_corners(), "Invalid index."); return vGloPos[co];}
249
251 inline const MathVector<dim>* local_corners() const
252 {return &vLocPos[0];}
253
256 {return &vGloPos[0];}
257
259 inline size_t node_id() const {return nodeId;}
260
262 inline size_t num_ip() const {return nip;}
263
265 inline const MathVector<dim>& local_ip() const {return vLocPos[0];}
266
268 inline const MathVector<worldDim>& global_ip() const {return vGloPos[0];}
269
271 inline const MathMatrix<worldDim,dim>& JTInv() const {return JtInv;}
272
274 inline number detJ() const {return detj;}
275
277 inline size_t num_sh() const {return nsh;}
278
280 inline number shape(size_t sh) const {return vShape[sh];}
281
283 inline const number* shape_vector() const {return vShape;}
284
286 inline const MathVector<dim>& local_grad(size_t sh) const
287 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vLocalGrad[sh];}
288
290 inline const MathVector<dim>* local_grad_vector() const {return vLocalGrad;}
291
293 inline const MathVector<worldDim>& global_grad(size_t sh) const
294 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vGlobalGrad[sh];}
295
298
299 private:
300 // let outer class access private members
301 friend class FV1IBGeometry<TElem, TWorldDim>;
302
303 // node id of associated node
304 size_t nodeId;
305
306 // volume of scv
308
309 // local and global positions of this element
310 MathVector<dim> vLocPos[numCo]; // local position of node
311 MathVector<worldDim> vGloPos[numCo]; // global position of node
312 MidID midId[numCo]; // dimension and id of object, that's midpoint bounds the scv
313
314 // shapes and derivatives
315 number vShape[nsh]; // shapes at ip
316 MathVector<dim> vLocalGrad[nsh]; // local grad at ip
317 MathVector<worldDim> vGlobalGrad[nsh]; // global grad at ip
318 MathMatrix<worldDim,dim> JtInv; // Jacobian transposed at ip
319 number detj; // Jacobian det at ip
320 };
321
323 class BF
324 {
325 public:
327 static const size_t numCo = traits::NumCornersOfSCVF;
328
329 public:
330 BF() {}
331
333 inline size_t node_id() const {return nodeId;}
334
336 inline size_t num_ip() const {return nip;}
337
339 inline const MathVector<dim>& local_ip() const {return localIP;}
340
342 inline const MathVector<worldDim>& global_ip() const {return globalIP;}
343
345 inline const MathVector<worldDim>& normal() const {return Normal;} // includes area
346
348 inline number volume() const {return Vol;}
349
351 inline const MathMatrix<worldDim, dim>& JTInv() const {return JtInv;}
352
354 inline number detJ() const {return detj;}
355
357 inline size_t num_sh() const {return nsh;}
358
360 inline number shape(size_t sh) const
361 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vShape[sh];}
362
364 inline const number* shape_vector() const {return vShape;}
365
367 inline const MathVector<dim>& local_grad(size_t sh) const
368 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vLocalGrad[sh];}
369
371 inline const MathVector<dim>* local_grad_vector() const {return vLocalGrad;}
372
374 inline const MathVector<worldDim>& global_grad(size_t sh) const
375 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vGlobalGrad[sh];}
376
379
381 inline size_t num_corners() const {return numCo;}
382
384 inline const MathVector<dim>& local_corner(size_t co) const
385 {UG_ASSERT(co < num_corners(), "Invalid index."); return vLocPos[co];}
386
388 inline const MathVector<worldDim>& global_corner(size_t co) const
389 {UG_ASSERT(co < num_corners(), "Invalid index."); return vGloPos[co];}
390
391 private:
393 friend class FV1IBGeometry<TElem, TWorldDim>;
394
395 // id of scv this bf belongs to
396 size_t nodeId;
397
398 // ordering is:
399 // 1D: edgeMidPoint
400 // 2D: edgeMidPoint, CenterOfElement
401 // 3D: edgeMidPoint, Side one, CenterOfElement, Side two
402 MathVector<dim> vLocPos[numCo]; // local corners of bf
403 MathVector<worldDim> vGloPos[numCo]; // global corners of bf
404 MidID vMidID[numCo]; // dimension and id of object, that's midpoint bounds the bf
405
406 // scvf part
407 MathVector<dim> localIP; // local integration point
408 MathVector<worldDim> globalIP; // global integration point
409 MathVector<worldDim> Normal; // normal (incl. area)
410 number Vol; // volume of bf
411
412 // shapes and derivatives
413 number vShape[nsh]; // shapes at ip
414 MathVector<dim> vLocalGrad[nsh]; // local grad at ip
415 MathVector<worldDim> vGlobalGrad[nsh]; // global grad at ip
416 MathMatrix<worldDim,dim> JtInv; // Jacobian transposed at ip
417 number detj; // Jacobian det at ip
418 };
419
420
421 public:
424
426 void adapt(GridObject* elem, const MathVector<TWorldDim>* vCornerCoords,
427 const ISubsetHandler* ish = NULL);
428
430 void adapt_normals(GridObject* elem, const MathVector<TWorldDim>* vCornerCoords,
431 const ISubsetHandler* ish = NULL);
432
434 void adapt_integration_points(GridObject* elem, const MathVector<TWorldDim>* vCornerCoords,
435 const ISubsetHandler* ish = NULL);
436
438 void update_local_data();
439
441 void update(GridObject* elem, const MathVector<worldDim>* vCornerCoords,
442 const ISubsetHandler* ish = NULL);
443
446 const MathVector<worldDim>* vCornerCoords,
447 const ISubsetHandler* ish = NULL);
448
450 const MathVector<worldDim>* corners() const {return m_vvGloMid[0];}
451
453 inline size_t num_scvf() const {return numSCVF;};
454
456 inline const SCVF& scvf(size_t i) const
457 {UG_ASSERT(i < num_scvf(), "Invalid Index."); return m_vSCVF[i];}
458
460 // do not use this method to obtain the number of shape functions,
461 // since this is NOT the same for pyramids; use num_sh() instead.
462 inline size_t num_scv() const {return numSCV;}
463
465 inline const SCV& scv(size_t i) const
466 {UG_ASSERT(i < num_scv(), "Invalid Index."); return m_vSCV[i];}
467
469 inline size_t num_sh() const {return nsh;};
470
471 public:
473 size_t num_scvf_ips() const {return numSCVF;}
474
477
480
482 size_t num_scv_ips() const {return numSCV;}
483
485 const MathVector<dim>* scv_local_ips() const {return &(m_vvLocMid[0][0]);}
486
488 const MathVector<worldDim>* scv_global_ips() const {return &(m_vvGloMid[0][0]);}
489
490
491 protected:
492 // global and local ips on SCVF
495
496 public:
498 inline void add_boundary_subset(int subsetIndex) {m_mapVectorBF[subsetIndex];}
499
501 inline void remove_boundary_subset(int subsetIndex) {m_mapVectorBF.erase(subsetIndex);}
502
504 inline void clear_boundary_subsets() {m_mapVectorBF.clear();}
505
507 inline size_t num_boundary_subsets() {return m_mapVectorBF.size();}
508
510 inline size_t num_bf() const
511 {
512 typename std::map<int, std::vector<BF> >::const_iterator it;
513 size_t num = 0;
514 for ( it=m_mapVectorBF.begin() ; it != m_mapVectorBF.end(); it++ )
515 num += (*it).second.size();
516 return num;
517 }
518
520 inline size_t num_bf(int si) const
521 {
522 typename std::map<int, std::vector<BF> >::const_iterator it;
523 it = m_mapVectorBF.find(si);
524 if(it == m_mapVectorBF.end()) return 0;
525 else return (*it).second.size();
526 }
527
529 inline const BF& bf(int si, size_t i) const
530 {
531 UG_ASSERT(i < num_bf(si), "Invalid index.");
532 typename std::map<int, std::vector<BF> >::const_iterator it;
533 it = m_mapVectorBF.find(si);
534 if(it == m_mapVectorBF.end()) UG_THROW("FVGeom: No bnd face for subset"<<si);
535 return (*it).second[i];
536 }
537
539 inline const std::vector<BF>& bf(int si) const
540 {
541 typename std::map<int, std::vector<BF> >::const_iterator it;
542 it = m_mapVectorBF.find(si);
543 if(it == m_mapVectorBF.end()) return m_vEmptyVectorBF;
544 return (*it).second;
545 }
546
547 protected:
548 std::map<int, std::vector<BF> > m_mapVectorBF;
549 std::vector<BF> m_vEmptyVectorBF;
550
551 protected:
553 TElem* m_pElem;
554
556 // (most objects in 1 dim, i.e. number of edges, but +1 for 1D)
557 static const int maxMid = numSCVF + 1;
558
562
565
568
571
574
577
578 public:
580 MathVector<worldDim> vSCVF_PosOnEdge[numSCVF]; // global position of scvf on edges
581 MathVector<worldDim> m_MidPoint; // global position of the midpoint, connecting each entry of PosOnEdges to a scvf
582
583
584};
585
587// Dim-dependent Finite Volume Geometry
589
591
596template <int TDim, int TWorldDim = TDim>
598{
599 public:
602
603 public:
605 static const int dim = TDim;
606
608 static const int worldDim = TWorldDim;
609
611 static const bool usesHangingNodes = false;
612
614 static const bool staticLocalData = false;
615
616 public:
618 static const int order = 1;
619
621 static const size_t maxNumSCV = traits::maxNumSCV;
622
624 typedef typename traits::scv_type scv_type;
625
627 static const size_t maxNumSCVF = traits::maxNumSCVF;
628
630 static const size_t maxNSH = traits::maxNSH;
631
633 static const size_t nip = 1;
634
635 public:
637
646 class SCVF
647 {
648 public:
650 static const size_t numCo = traits::NumCornersOfSCVF;
651
652 public:
653 SCVF() {}
654
656 inline size_t from() const {return From;}
657
659 inline size_t to() const {return To;}
660
662 inline size_t num_ip() const {return nip;}
663
665 inline const MathVector<dim>& local_ip() const {return localIP;}
666
668 inline const MathVector<worldDim>& global_ip() const {return globalIP;}
669
671 inline const MathVector<worldDim>& normal() const {return Normal;}
672
674 inline const MathMatrix<worldDim,dim>& JTInv() const {return JtInv;}
675
677 inline number detJ() const {return detj;}
678
680 inline size_t num_sh() const {return numSH;}
681
683 inline number shape(size_t sh) const {return vShape[sh];}
684
686 inline const number* shape_vector() const {return vShape;}
687
689 inline const MathVector<dim>& local_grad(size_t sh) const
690 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vLocalGrad[sh];}
691
693 inline const MathVector<dim>* local_grad_vector() const {return vLocalGrad;}
694
696 inline const MathVector<worldDim>& global_grad(size_t sh) const
697 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vGlobalGrad[sh];}
698
701
703 inline size_t num_corners() const {return numCo;}
704
706 inline const MathVector<dim>& local_corner(size_t co) const
707 {UG_ASSERT(co < num_corners(), "Invalid index."); return vLocPos[co];}
708
710 inline const MathVector<worldDim>& global_corner(size_t co) const
711 {UG_ASSERT(co < num_corners(), "Invalid index."); return vGloPos[co];}
712
713 private:
714 // let outer class access private members
715 friend class DimFV1IBGeometry<dim, worldDim>;
716
717 // This scvf separates the scv with the ids given in "from" and "to"
718 // The computed normal points in direction from->to
719 size_t From, To;
720
721 // The normal on the SCVF pointing (from -> to)
722 MathVector<worldDim> Normal; // normal (incl. area)
723
724 // ordering is:
725 // 1D: edgeMidPoint
726 // 2D: edgeMidPoint, CenterOfElement
727 // 3D: edgeMidPoint, Side one, CenterOfElement, Side two
728 MathVector<dim> vLocPos[numCo]; // local corners of scvf
729 MathVector<worldDim> vGloPos[numCo]; // global corners of scvf
730 MidID vMidID[numCo]; // dimension and id of object, that's midpoint bounds the scvf
731
732 // scvf part
733 MathVector<dim> localIP; // local integration point
734 MathVector<worldDim> globalIP; // global integration point
735
736 // shapes and derivatives
737 size_t numSH;
738 number vShape[maxNSH]; // shapes at ip
739 MathVector<dim> vLocalGrad[maxNSH]; // local grad at ip
741 MathMatrix<worldDim,dim> JtInv; // Jacobian transposed at ip
742 number detj; // Jacobian det at ip
743 };
744
746 class SCV
747 {
748 public:
750 static const size_t numCo = traits::NumCornersOfSCV;
751
752 public:
753 SCV() {};
754
756 inline number volume() const {return Vol;}
757
759 inline size_t num_corners() const {return numCo;}
760
762 inline const MathVector<dim>& local_corner(size_t co) const
763 {UG_ASSERT(co < num_corners(), "Invalid index."); return vLocPos[co];}
764
766 inline const MathVector<worldDim>& global_corner(size_t co) const
767 {UG_ASSERT(co < num_corners(), "Invalid index."); return vGloPos[co];}
768
770 inline size_t node_id() const {return nodeId;}
771
773 inline size_t num_ip() const {return nip;}
774
776 inline const MathVector<dim>& local_ip() const {return vLocPos[0];}
777
779 inline const MathVector<worldDim>& global_ip() const {return vGloPos[0];}
780
782 inline const MathMatrix<worldDim,dim>& JTInv() const {return JtInv;}
783
785 inline number detJ() const {return detj;}
786
788 inline size_t num_sh() const {return numSH;}
789
791 inline number shape(size_t sh) const {return vShape[sh];}
792
794 inline const number* shape_vector() const {return vShape;}
795
797 inline const MathVector<dim>& local_grad(size_t sh) const
798 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vLocalGrad[sh];}
799
801 inline const MathVector<dim>* local_grad_vector() const {return vLocalGrad;}
802
804 inline const MathVector<worldDim>& global_grad(size_t sh) const
805 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vGlobalGrad[sh];}
806
809
810 private:
811 // let outer class access private members
812 friend class DimFV1IBGeometry<dim, worldDim>;
813
814 // node id of associated node
815 size_t nodeId;
816
817 // volume of scv
819
820 // local and global positions of this element
821 MathVector<dim> vLocPos[numCo]; // local position of node
822 MathVector<worldDim> vGloPos[numCo]; // global position of node
823 MidID vMidID[numCo]; // dimension and id of object, that's midpoint bounds the scv
824
825 // shapes and derivatives
826 size_t numSH;
827 number vShape[maxNSH]; // shapes at ip
828 MathVector<dim> vLocalGrad[maxNSH]; // local grad at ip
830 MathMatrix<worldDim,dim> JtInv; // Jacobian transposed at ip
831 number detj; // Jacobian det at ip
832 };
833
835 class BF
836 {
837 public:
839 static const size_t numCo = traits::NumCornersOfSCVF;
840
841 public:
842 BF() {}
843
845 inline size_t node_id() const {return nodeId;}
846
848 inline size_t num_ip() const {return nip;}
849
851 inline const MathVector<dim>& local_ip() const {return localIP;}
852
854 inline const MathVector<worldDim>& global_ip() const {return globalIP;}
855
857 inline const MathVector<worldDim>& normal() const {return Normal;} // includes area
858
860 inline number volume() const {return Vol;}
861
863 inline const MathMatrix<worldDim, dim>& JTInv() const {return JtInv;}
864
866 inline number detJ() const {return detj;}
867
869 inline size_t num_sh() const {return numSH;}
870
872 inline number shape(size_t sh) const
873 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vShape[sh];}
874
876 inline const number* shape_vector() const {return vShape;}
877
879 inline const MathVector<dim>& local_grad(size_t sh) const
880 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vLocalGrad[sh];}
881
883 inline const MathVector<dim>* local_grad_vector() const {return vLocalGrad;}
884
886 inline const MathVector<worldDim>& global_grad(size_t sh) const
887 {UG_ASSERT(sh < num_sh(), "Invalid index"); return vGlobalGrad[sh];}
888
891
893 inline size_t num_corners() const {return numCo;}
894
896 inline const MathVector<dim>& local_corner(size_t co) const
897 {UG_ASSERT(co < num_corners(), "Invalid index."); return vLocPos[co];}
898
900 inline const MathVector<worldDim>& global_corner(size_t co) const
901 {UG_ASSERT(co < num_corners(), "Invalid index."); return vGloPos[co];}
902
903 private:
905 friend class DimFV1IBGeometry<dim, worldDim>;
906
907 // id of scv this bf belongs to
908 size_t nodeId;
909
910 // ordering is:
911 // 1D: edgeMidPoint
912 // 2D: edgeMidPoint, CenterOfElement
913 // 3D: edgeMidPoint, Side one, CenterOfElement, Side two
914 MathVector<dim> vLocPos[numCo]; // local corners of bf
915 MathVector<worldDim> vGloPos[numCo]; // global corners of bf
916 MidID vMidID[numCo]; // dimension and id of object, that's midpoint bounds the bf
917
918 // scvf part
919 MathVector<dim> localIP; // local integration point
920 MathVector<worldDim> globalIP; // global integration point
921 MathVector<worldDim> Normal; // normal (incl. area)
922 number Vol; // volume of bf
923
924 // shapes and derivatives
925 size_t numSH;
926 number vShape[maxNSH]; // shapes at ip
927 MathVector<dim> vLocalGrad[maxNSH]; // local grad at ip
929 MathMatrix<worldDim,dim> JtInv; // Jacobian transposed at ip
930 number detj; // Jacobian det at ip
931 };
932
933 public:
936
938 void adapt(GridObject* elem, const MathVector<TWorldDim>* vCornerCoords,
939 const ISubsetHandler* ish = NULL);
940
942 void adapt_normals(GridObject* elem, const MathVector<TWorldDim>* vCornerCoords,
943 const ISubsetHandler* ish = NULL);
944
946 void adapt_integration_points(GridObject* elem, const MathVector<TWorldDim>* vCornerCoords,
947 const ISubsetHandler* ish = NULL);
948
950 void update_local_data();
951
953 void update(GridObject* elem, const MathVector<worldDim>* vCornerCoords,
954 const ISubsetHandler* ish = NULL);
955
958 const MathVector<worldDim>* vCornerCoords,
959 const ISubsetHandler* ish = NULL);
960
962 const MathVector<worldDim>* corners() const {return m_vvGloMid[0];}
963
965 inline size_t num_scvf() const {return m_numSCVF;};
966
968 inline const SCVF& scvf(size_t i) const
969 {UG_ASSERT(i < num_scvf(), "Invalid Index."); return m_vSCVF[i];}
970
972 // do not use this method to obtain the number of shape functions,
973 // since this is NOT the same for pyramids; use num_sh() instead.
974 inline size_t num_scv() const {return m_numSCV;}
975
977 inline const SCV& scv(size_t i) const
978 {UG_ASSERT(i < num_scv(), "Invalid Index."); return m_vSCV[i];}
979
981 inline size_t num_sh() const {return m_nsh;};
982
983 public:
985 size_t num_scvf_ips() const {return m_numSCVF;}
986
989
992
994 size_t num_scv_ips() const {return m_numSCV;}
995
997 const MathVector<dim>* scv_local_ips() const {return &(m_vvLocMid[0][0]);}
998
1000 const MathVector<worldDim>* scv_global_ips() const {return &(m_vvGloMid[0][0]);}
1001
1002
1003 protected:
1004 // global and local ips on SCVF
1007
1008 public:
1010 inline void add_boundary_subset(int subsetIndex) {m_mapVectorBF[subsetIndex];}
1011
1013 inline void remove_boundary_subset(int subsetIndex) {m_mapVectorBF.erase(subsetIndex);}
1014
1016 inline void clear_boundary_subsets() {m_mapVectorBF.clear();}
1017
1019 inline size_t num_boundary_subsets() {return m_mapVectorBF.size();}
1020
1022 inline size_t num_bf() const
1023 {
1024 typename std::map<int, std::vector<BF> >::const_iterator it;
1025 size_t num = 0;
1026 for ( it=m_mapVectorBF.begin() ; it != m_mapVectorBF.end(); it++ )
1027 num += (*it).second.size();
1028 return num;
1029 }
1030
1032 inline size_t num_bf(int si) const
1033 {
1034 typename std::map<int, std::vector<BF> >::const_iterator it;
1035 it = m_mapVectorBF.find(si);
1036 if(it == m_mapVectorBF.end()) return 0;
1037 else return (*it).second.size();
1038 }
1039
1041 inline const BF& bf(int si, size_t i) const
1042 {
1043 typename std::map<int, std::vector<BF> >::const_iterator it;
1044 it = m_mapVectorBF.find(si);
1045 if(it == m_mapVectorBF.end()) UG_THROW("DimFV1IBGeom: No BndSubset "<<si);
1046 return (*it).second[i];
1047 }
1048
1050 inline const std::vector<BF>& bf(int si) const
1051 {
1052 typename std::map<int, std::vector<BF> >::const_iterator it;
1053 it = m_mapVectorBF.find(si);
1054 if(it == m_mapVectorBF.end()) return m_vEmptyVectorBF;
1055 return (*it).second;
1056 }
1057
1058 protected:
1059 std::map<int, std::vector<BF> > m_mapVectorBF;
1060 std::vector<BF> m_vEmptyVectorBF;
1061
1062 private:
1065
1068
1070 size_t m_numSCV;
1071
1074
1076 size_t m_nsh;
1077
1079 // (most objects in 1 dim, i.e. number of edges, but +1 for 1D)
1080 static const int maxMid = maxNumSCVF + 1;
1081
1085
1088
1091
1092 public:
1094 MathVector<worldDim> vSCVF_PosOnEdge[maxNumSCVF]; // global position of scvf on edges
1095 MathVector<worldDim> m_MidPoint; // global position of the midpoint, connecting each entry of PosOnEdges to a scvf
1096
1097};
1098
1099
1100
1101#include "fv1ib_geom_impl.h"
1102
1103} // end ug namespace
1104
1105#endif /* FV1IB_GEOM_H_ */
boundary face
Definition fv1ib_geom.h:836
MathVector< dim > vLocPos[numCo]
Definition fv1ib_geom.h:914
size_t nodeId
Definition fv1ib_geom.h:908
size_t node_id() const
index of SubControlVolume of the bf
Definition fv1ib_geom.h:845
const MathVector< worldDim > * global_grad_vector() const
vector of global gradients in ip point
Definition fv1ib_geom.h:890
const MathMatrix< worldDim, dim > & JTInv() const
Transposed Inverse of Jacobian in integration point.
Definition fv1ib_geom.h:863
number Vol
Definition fv1ib_geom.h:922
number detj
Definition fv1ib_geom.h:930
const MathVector< worldDim > & global_grad(size_t sh) const
value of global gradient of shape function i in integration point
Definition fv1ib_geom.h:886
const MathVector< worldDim > & global_ip() const
global integration point of bf
Definition fv1ib_geom.h:854
size_t num_ip() const
number of integration points on bf
Definition fv1ib_geom.h:848
MathVector< worldDim > vGlobalGrad[maxNSH]
Definition fv1ib_geom.h:928
MathVector< dim > vLocalGrad[maxNSH]
Definition fv1ib_geom.h:927
size_t numSH
Definition fv1ib_geom.h:925
BF()
Definition fv1ib_geom.h:842
number shape(size_t sh) const
value of shape function i in integration point
Definition fv1ib_geom.h:872
size_t num_sh() const
number of shape functions
Definition fv1ib_geom.h:869
MathVector< worldDim > globalIP
Definition fv1ib_geom.h:920
const MathVector< dim > & local_corner(size_t co) const
return local corner number i
Definition fv1ib_geom.h:896
MathVector< worldDim > Normal
Definition fv1ib_geom.h:921
size_t num_corners() const
number of corners, that bound the scvf
Definition fv1ib_geom.h:893
static const size_t numCo
Number of corners of bf.
Definition fv1ib_geom.h:839
number vShape[maxNSH]
Definition fv1ib_geom.h:926
const MathVector< dim > & local_ip() const
local integration point of bf
Definition fv1ib_geom.h:851
number volume() const
volume of bf
Definition fv1ib_geom.h:860
MidID vMidID[numCo]
Definition fv1ib_geom.h:916
MathVector< dim > localIP
Definition fv1ib_geom.h:919
const MathVector< dim > & local_grad(size_t sh) const
value of local gradient of shape function i in integration point
Definition fv1ib_geom.h:879
number detJ() const
Determinant of Jacobian in integration point.
Definition fv1ib_geom.h:866
MathVector< worldDim > vGloPos[numCo]
Definition fv1ib_geom.h:915
const MathVector< worldDim > & normal() const
outer normal on bf. Norm is equal to area
Definition fv1ib_geom.h:857
const MathVector< worldDim > & global_corner(size_t co) const
return global corner number i
Definition fv1ib_geom.h:900
const number * shape_vector() const
vector of local gradients in ip point
Definition fv1ib_geom.h:876
MathMatrix< worldDim, dim > JtInv
Definition fv1ib_geom.h:929
const MathVector< dim > * local_grad_vector() const
vector of local gradients in ip point
Definition fv1ib_geom.h:883
Sub-Control Volume Face structure.
Definition fv1ib_geom.h:647
MathVector< worldDim > Normal
Definition fv1ib_geom.h:722
MathVector< dim > vLocalGrad[maxNSH]
Definition fv1ib_geom.h:739
size_t num_ip() const
number of integration points on scvf
Definition fv1ib_geom.h:662
size_t numSH
Definition fv1ib_geom.h:737
number vShape[maxNSH]
Definition fv1ib_geom.h:738
number shape(size_t sh) const
value of shape function i in integration point
Definition fv1ib_geom.h:683
MathVector< dim > vLocPos[numCo]
Definition fv1ib_geom.h:728
const MathVector< worldDim > & global_corner(size_t co) const
return glbal corner number i
Definition fv1ib_geom.h:710
const MathVector< worldDim > * global_grad_vector() const
vector of gloabl gradients in ip point
Definition fv1ib_geom.h:700
MathVector< dim > localIP
Definition fv1ib_geom.h:733
MathVector< worldDim > vGlobalGrad[maxNSH]
Definition fv1ib_geom.h:740
size_t from() const
index of SubControlVolume on one side of the scvf
Definition fv1ib_geom.h:656
size_t From
Definition fv1ib_geom.h:719
const MathVector< dim > & local_grad(size_t sh) const
value of local gradient of shape function i in integration point
Definition fv1ib_geom.h:689
const MathVector< worldDim > & global_ip() const
global integration point of scvf
Definition fv1ib_geom.h:668
number detj
Definition fv1ib_geom.h:742
size_t To
Definition fv1ib_geom.h:719
MathVector< worldDim > globalIP
Definition fv1ib_geom.h:734
MathVector< worldDim > vGloPos[numCo]
Definition fv1ib_geom.h:729
MidID vMidID[numCo]
Definition fv1ib_geom.h:730
size_t num_sh() const
number of shape functions
Definition fv1ib_geom.h:680
size_t num_corners() const
number of corners, that bound the scvf
Definition fv1ib_geom.h:703
const MathVector< dim > & local_corner(size_t co) const
return local corner number i
Definition fv1ib_geom.h:706
static const size_t numCo
Number of corners of scvf.
Definition fv1ib_geom.h:650
MathMatrix< worldDim, dim > JtInv
Definition fv1ib_geom.h:741
SCVF()
Definition fv1ib_geom.h:653
const number * shape_vector() const
vector of shape functions in ip point
Definition fv1ib_geom.h:686
const MathVector< worldDim > & global_grad(size_t sh) const
value of global gradient of shape function i in integration point
Definition fv1ib_geom.h:696
const MathVector< dim > * local_grad_vector() const
vector of local gradients in ip point
Definition fv1ib_geom.h:693
const MathMatrix< worldDim, dim > & JTInv() const
Transposed Inverse of Jacobian in integration point.
Definition fv1ib_geom.h:674
const MathVector< worldDim > & normal() const
normal on scvf (points direction "from"->"to"). Norm is equal to area
Definition fv1ib_geom.h:671
number detJ() const
Determinant of Jacobian in integration point.
Definition fv1ib_geom.h:677
const MathVector< dim > & local_ip() const
local integration point of scvf
Definition fv1ib_geom.h:665
size_t to() const
index of SubControlVolume on one side of the scvf
Definition fv1ib_geom.h:659
sub control volume structure
Definition fv1ib_geom.h:747
const number * shape_vector() const
vector of shape functions in ip point
Definition fv1ib_geom.h:794
size_t numSH
Definition fv1ib_geom.h:826
number detj
Definition fv1ib_geom.h:831
size_t nodeId
Definition fv1ib_geom.h:815
const MathVector< worldDim > & global_grad(size_t sh) const
value of global gradient of shape function i in integration point
Definition fv1ib_geom.h:804
size_t node_id() const
node id that this scv is associated to
Definition fv1ib_geom.h:770
const MathVector< dim > * local_grad_vector() const
vector of local gradients in ip point
Definition fv1ib_geom.h:801
MathVector< worldDim > vGloPos[numCo]
Definition fv1ib_geom.h:822
SCV()
Definition fv1ib_geom.h:753
number vShape[maxNSH]
Definition fv1ib_geom.h:827
MidID vMidID[numCo]
Definition fv1ib_geom.h:823
const MathVector< worldDim > * global_grad_vector() const
vector of global gradients in ip point
Definition fv1ib_geom.h:808
const MathVector< dim > & local_corner(size_t co) const
return local corner number i
Definition fv1ib_geom.h:762
static const size_t numCo
Number of corners of scv.
Definition fv1ib_geom.h:750
const MathVector< worldDim > & global_ip() const
global integration point
Definition fv1ib_geom.h:779
const MathVector< dim > & local_grad(size_t sh) const
value of local gradient of shape function i in integration point
Definition fv1ib_geom.h:797
MathMatrix< worldDim, dim > JtInv
Definition fv1ib_geom.h:830
number shape(size_t sh) const
value of shape function i in integration point
Definition fv1ib_geom.h:791
number detJ() const
Determinant of Jacobian in integration point.
Definition fv1ib_geom.h:785
number volume() const
volume of scv
Definition fv1ib_geom.h:756
MathVector< dim > vLocalGrad[maxNSH]
Definition fv1ib_geom.h:828
MathVector< worldDim > vGlobalGrad[maxNSH]
Definition fv1ib_geom.h:829
const MathVector< worldDim > & global_corner(size_t co) const
return global corner number i
Definition fv1ib_geom.h:766
size_t num_sh() const
number of shape functions
Definition fv1ib_geom.h:788
size_t num_ip() const
number of integration points
Definition fv1ib_geom.h:773
number Vol
Definition fv1ib_geom.h:818
size_t num_corners() const
number of corners, that bound the scv
Definition fv1ib_geom.h:759
MathVector< dim > vLocPos[numCo]
Definition fv1ib_geom.h:821
const MathVector< dim > & local_ip() const
local integration point of scv
Definition fv1ib_geom.h:776
const MathMatrix< worldDim, dim > & JTInv() const
Transposed Inverse of Jacobian in integration point.
Definition fv1ib_geom.h:782
Geometry and shape functions for 1st order Vertex-Centered Finite Volume.
Definition fv1ib_geom.h:598
MathVector< worldDim > m_vGlobSCVF_IP[maxNumSCVF]
Definition fv1ib_geom.h:1005
size_t num_scvf_ips() const
returns number of all scvf ips
Definition fv1ib_geom.h:985
const MathVector< dim > * scvf_local_ips() const
returns all ips of scvf as they appear in scv loop
Definition fv1ib_geom.h:988
size_t num_boundary_subsets()
number of registered boundary subsets
Definition fv1ib_geom.h:1019
DimFV1IBGeometry()
construct object and initialize local values and sizes
Definition fv1ib_geom.h:935
ReferenceObjectID m_roid
current reference object id
Definition fv1ib_geom.h:1067
void add_boundary_subset(int subsetIndex)
add subset that is interpreted as boundary subset.
Definition fv1ib_geom.h:1010
size_t num_bf() const
number of all boundary faces
Definition fv1ib_geom.h:1022
size_t m_nsh
current number of shape functions
Definition fv1ib_geom.h:1076
SCV m_vSCV[maxNumSCV]
SubControlVolumes.
Definition fv1ib_geom.h:1090
MathVector< worldDim > m_vvGloMid[dim+1][maxMid]
Definition fv1ib_geom.h:1084
const SCV & scv(size_t i) const
const access to SubControlVolume number i
Definition fv1ib_geom.h:977
size_t m_numSCV
current number of scv
Definition fv1ib_geom.h:1070
std::map< int, std::vector< BF > > m_mapVectorBF
Definition fv1ib_geom.h:1059
size_t m_numSCVF
current number of scvf
Definition fv1ib_geom.h:1073
void adapt_integration_points(GridObject *elem, const MathVector< TWorldDim > *vCornerCoords, const ISubsetHandler *ish=NULL)
adapt integration points for elements cut by the inner boundary
Definition fv1ib_geom_impl.h:921
void update_local_data()
update local data
Definition fv1ib_geom_impl.h:951
const MathVector< dim > * scv_local_ips() const
returns all ips of scv as they appear in scv loop
Definition fv1ib_geom.h:997
static const int order
order
Definition fv1ib_geom.h:618
static const bool staticLocalData
flag indicating if local data may change
Definition fv1ib_geom.h:614
size_t num_scvf() const
number of SubControlVolumeFaces
Definition fv1ib_geom.h:965
size_t num_bf(int si) const
number of boundary faces on subset 'subsetIndex'
Definition fv1ib_geom.h:1032
const BF & bf(int si, size_t i) const
returns the boundary face i for subsetIndex
Definition fv1ib_geom.h:1041
void remove_boundary_subset(int subsetIndex)
removes subset that is interpreted as boundary subset.
Definition fv1ib_geom.h:1013
void update(GridObject *elem, const MathVector< worldDim > *vCornerCoords, const ISubsetHandler *ish=NULL)
update data for given element
Definition fv1ib_geom_impl.h:1057
const MathVector< worldDim > * scvf_global_ips() const
returns all ips of scvf as they appear in scv loop
Definition fv1ib_geom.h:991
static const int worldDim
dimension of world
Definition fv1ib_geom.h:608
MathVector< worldDim > m_MidPoint
Definition fv1ib_geom.h:1095
static const size_t nip
number of integration points
Definition fv1ib_geom.h:633
SCVF m_vSCVF[maxNumSCVF]
SubControlVolumeFaces.
Definition fv1ib_geom.h:1087
size_t num_scv() const
number of SubControlVolumes
Definition fv1ib_geom.h:974
MathVector< dim > m_vLocSCVF_IP[maxNumSCVF]
Definition fv1ib_geom.h:1006
const MathVector< worldDim > * corners() const
get vector of corners for current element
Definition fv1ib_geom.h:962
std::vector< BF > m_vEmptyVectorBF
Definition fv1ib_geom.h:1060
GridObject * m_pElem
pointer to current element
Definition fv1ib_geom.h:1064
static const size_t maxNSH
max number of shape functions
Definition fv1ib_geom.h:630
void adapt(GridObject *elem, const MathVector< TWorldDim > *vCornerCoords, const ISubsetHandler *ish=NULL)
adapt integration points for elements cut by the inner boundary
Definition fv1ib_geom_impl.h:858
MathVector< worldDim > vSCVF_PosOnEdge[maxNumSCVF]
new members for adaptions
Definition fv1ib_geom.h:1094
void adapt_normals(GridObject *elem, const MathVector< TWorldDim > *vCornerCoords, const ISubsetHandler *ish=NULL)
adapt normals for elements cut by the inner boundary
Definition fv1ib_geom_impl.h:930
MathVector< dim > m_vvLocMid[dim+1][maxMid]
local and global geom object midpoints for each dimension
Definition fv1ib_geom.h:1083
const MathVector< worldDim > * scv_global_ips() const
returns all ips of scv as they appear in scv loop
Definition fv1ib_geom.h:1000
traits::scv_type scv_type
type of SubControlVolume
Definition fv1ib_geom.h:624
size_t num_sh() const
number of shape functions
Definition fv1ib_geom.h:981
static const size_t maxNumSCV
number of SubControlVolumes
Definition fv1ib_geom.h:621
static const int maxMid
max number of geometric objects in a dimension
Definition fv1ib_geom.h:1080
static const size_t maxNumSCVF
max number of SubControlVolumeFaces
Definition fv1ib_geom.h:627
const std::vector< BF > & bf(int si) const
returns reference to vector of boundary faces for subsetIndex
Definition fv1ib_geom.h:1050
fv1_dim_traits< TDim, TWorldDim > traits
used traits
Definition fv1ib_geom.h:601
static const int dim
dimension of reference element
Definition fv1ib_geom.h:605
static const bool usesHangingNodes
Hanging node flag: this Geometry does not support hanging nodes.
Definition fv1ib_geom.h:611
const SCVF & scvf(size_t i) const
const access to SubControlVolumeFace number i
Definition fv1ib_geom.h:968
void update_boundary_faces(GridObject *elem, const MathVector< worldDim > *vCornerCoords, const ISubsetHandler *ish=NULL)
update boundary data for given element
Definition fv1ib_geom_impl.h:1169
size_t num_scv_ips() const
returns number of all scv ips
Definition fv1ib_geom.h:994
void clear_boundary_subsets()
reset all boundary subsets
Definition fv1ib_geom.h:1016
boundary face
Definition fv1ib_geom.h:324
const MathVector< dim > & local_grad(size_t sh) const
value of local gradient of shape function i in integration point
Definition fv1ib_geom.h:367
size_t num_corners() const
number of corners, that bound the scvf
Definition fv1ib_geom.h:381
const MathVector< dim > * local_grad_vector() const
vector of local gradients in ip point
Definition fv1ib_geom.h:371
const MathVector< dim > & local_ip() const
local integration point of bf
Definition fv1ib_geom.h:339
MathVector< worldDim > vGloPos[numCo]
Definition fv1ib_geom.h:403
size_t num_sh() const
number of shape functions
Definition fv1ib_geom.h:357
const MathVector< worldDim > & normal() const
outer normal on bf. Norm is equal to area
Definition fv1ib_geom.h:345
MathMatrix< worldDim, dim > JtInv
Definition fv1ib_geom.h:416
MathVector< worldDim > Normal
Definition fv1ib_geom.h:409
const MathVector< dim > & local_corner(size_t co) const
return local corner number i
Definition fv1ib_geom.h:384
MathVector< dim > vLocalGrad[nsh]
Definition fv1ib_geom.h:414
BF()
Definition fv1ib_geom.h:330
number detJ() const
Determinant of Jacobian in integration point.
Definition fv1ib_geom.h:354
const MathVector< worldDim > & global_ip() const
global integration point of bf
Definition fv1ib_geom.h:342
const MathMatrix< worldDim, dim > & JTInv() const
Transposed Inverse of Jacobian in integration point.
Definition fv1ib_geom.h:351
number shape(size_t sh) const
value of shape function i in integration point
Definition fv1ib_geom.h:360
const MathVector< worldDim > & global_grad(size_t sh) const
value of global gradient of shape function i in integration point
Definition fv1ib_geom.h:374
MidID vMidID[numCo]
Definition fv1ib_geom.h:404
static const size_t numCo
Number of corners of bf.
Definition fv1ib_geom.h:327
number detj
Definition fv1ib_geom.h:417
const MathVector< worldDim > & global_corner(size_t co) const
return global corner number i
Definition fv1ib_geom.h:388
MathVector< dim > localIP
Definition fv1ib_geom.h:407
size_t nodeId
Definition fv1ib_geom.h:396
size_t num_ip() const
number of integration points on bf
Definition fv1ib_geom.h:336
size_t node_id() const
index of SubControlVolume of the bf
Definition fv1ib_geom.h:333
number volume() const
volume of bf
Definition fv1ib_geom.h:348
MathVector< worldDim > globalIP
Definition fv1ib_geom.h:408
number Vol
Definition fv1ib_geom.h:410
number vShape[nsh]
Definition fv1ib_geom.h:413
const number * shape_vector() const
vector of local gradients in ip point
Definition fv1ib_geom.h:364
MathVector< dim > vLocPos[numCo]
Definition fv1ib_geom.h:402
const MathVector< worldDim > * global_grad_vector() const
vector of global gradients in ip point
Definition fv1ib_geom.h:378
MathVector< worldDim > vGlobalGrad[nsh]
Definition fv1ib_geom.h:415
Sub-Control Volume Face structure.
Definition fv1ib_geom.h:129
const MathVector< dim > & local_ip() const
local integration point of scvf
Definition fv1ib_geom.h:150
size_t num_corners() const
number of corners, that bound the scvf
Definition fv1ib_geom.h:185
number detJ() const
Determinant of Jacobian in integration point.
Definition fv1ib_geom.h:159
MathVector< dim > vLocPos[numCo]
Definition fv1ib_geom.h:210
MathVector< worldDim > vGlobalGrad[nsh]
Definition fv1ib_geom.h:221
const MathVector< worldDim > & global_ip() const
global integration point of scvf
Definition fv1ib_geom.h:153
number shape(size_t sh) const
value of shape function i in integration point
Definition fv1ib_geom.h:165
const MathVector< dim > & local_corner(size_t co) const
return local corner number i
Definition fv1ib_geom.h:188
MathVector< worldDim > globalIP
Definition fv1ib_geom.h:216
const MathVector< dim > * local_grad_vector() const
vector of local gradients in ip point
Definition fv1ib_geom.h:175
MathVector< dim > localIP
Definition fv1ib_geom.h:215
MidID vMidID[numCo]
Definition fv1ib_geom.h:212
const MathMatrix< worldDim, dim > & JTInv() const
Transposed Inverse of Jacobian in integration point.
Definition fv1ib_geom.h:156
const number * shape_vector() const
vector of shape functions in ip point
Definition fv1ib_geom.h:168
size_t To
Definition fv1ib_geom.h:201
number detj
Definition fv1ib_geom.h:223
const MathVector< worldDim > & global_corner(size_t co) const
return global corner number i
Definition fv1ib_geom.h:192
size_t to() const
index of SubControlVolume on one side of the scvf
Definition fv1ib_geom.h:141
MathVector< worldDim > vGloPos[numCo]
Definition fv1ib_geom.h:211
size_t from() const
index of SubControlVolume on one side of the scvf
Definition fv1ib_geom.h:138
const MathVector< worldDim > & normal() const
normal on scvf (points direction "from"->"to"). Norm is equal to area
Definition fv1ib_geom.h:144
number vShape[nsh]
Definition fv1ib_geom.h:219
size_t From
Definition fv1ib_geom.h:201
size_t num_sh() const
number of shape functions
Definition fv1ib_geom.h:162
MathVector< worldDim > Normal
Definition fv1ib_geom.h:204
const MathVector< worldDim > & global_grad(size_t sh) const
value of global gradient of shape function i in integration point
Definition fv1ib_geom.h:178
const MathVector< dim > & local_grad(size_t sh) const
value of local gradient of shape function i in integration point
Definition fv1ib_geom.h:171
MathVector< dim > vLocalGrad[nsh]
Definition fv1ib_geom.h:220
static const size_t numCo
Number of corners of scvf.
Definition fv1ib_geom.h:132
MathMatrix< worldDim, dim > JtInv
Definition fv1ib_geom.h:222
SCVF()
Definition fv1ib_geom.h:135
size_t num_ip() const
number of integration points on scvf
Definition fv1ib_geom.h:147
const MathVector< worldDim > * global_grad_vector() const
vector of global gradients in ip point
Definition fv1ib_geom.h:182
sub control volume structure
Definition fv1ib_geom.h:228
size_t num_ip() const
number of integration points
Definition fv1ib_geom.h:262
MidID midId[numCo]
Definition fv1ib_geom.h:312
MathVector< dim > vLocPos[numCo]
Definition fv1ib_geom.h:310
MathVector< worldDim > vGloPos[numCo]
Definition fv1ib_geom.h:311
number Vol
Definition fv1ib_geom.h:307
MathMatrix< worldDim, dim > JtInv
Definition fv1ib_geom.h:318
const MathVector< worldDim > & global_corner(size_t co) const
return global corner number i
Definition fv1ib_geom.h:247
const MathMatrix< worldDim, dim > & JTInv() const
Transposed Inverse of Jacobian in integration point.
Definition fv1ib_geom.h:271
number volume() const
volume of scv
Definition fv1ib_geom.h:237
const MathVector< dim > & local_ip() const
local integration point of scv
Definition fv1ib_geom.h:265
size_t node_id() const
node id that this scv is associated to
Definition fv1ib_geom.h:259
size_t nodeId
Definition fv1ib_geom.h:304
const MathVector< dim > & local_corner(size_t co) const
return local corner number i
Definition fv1ib_geom.h:243
MathVector< dim > vLocalGrad[nsh]
Definition fv1ib_geom.h:316
MathVector< worldDim > vGlobalGrad[nsh]
Definition fv1ib_geom.h:317
static const size_t numCo
Number of corners of scvf.
Definition fv1ib_geom.h:231
const MathVector< worldDim > & global_ip() const
global integration point
Definition fv1ib_geom.h:268
const MathVector< dim > * local_grad_vector() const
vector of local gradients in ip point
Definition fv1ib_geom.h:290
number detj
Definition fv1ib_geom.h:319
number shape(size_t sh) const
value of shape function i in integration point
Definition fv1ib_geom.h:280
size_t num_sh() const
number of shape functions
Definition fv1ib_geom.h:277
size_t num_corners() const
number of corners, that bound the scvf
Definition fv1ib_geom.h:240
SCV()
Definition fv1ib_geom.h:234
const MathVector< worldDim > * global_grad_vector() const
vector of global gradients in ip point
Definition fv1ib_geom.h:297
const MathVector< dim > * local_corners() const
return local corners
Definition fv1ib_geom.h:251
const MathVector< dim > & local_grad(size_t sh) const
value of local gradient of shape function i in integration point
Definition fv1ib_geom.h:286
const number * shape_vector() const
vector of shape functions in ip point
Definition fv1ib_geom.h:283
number detJ() const
Determinant of Jacobian in integration point.
Definition fv1ib_geom.h:274
const MathVector< worldDim > & global_grad(size_t sh) const
value of global gradient of shape function i in integration point
Definition fv1ib_geom.h:293
const MathVector< worldDim > * global_corners() const
return global corners
Definition fv1ib_geom.h:255
number vShape[nsh]
Definition fv1ib_geom.h:315
Geometry and shape functions for 1st order Vertex-Centered Finite Volume.
Definition fv1ib_geom.h:69
static const int maxMid
max number of geom objects in all dimensions
Definition fv1ib_geom.h:557
size_t num_sh() const
number of shape functions
Definition fv1ib_geom.h:469
MathVector< worldDim > m_vGlobSCVF_IP[numSCVF]
Definition fv1ib_geom.h:493
void remove_boundary_subset(int subsetIndex)
removes subset that is interpreted as boundary subset.
Definition fv1ib_geom.h:501
const MathVector< worldDim > * scv_global_ips() const
returns all ips of scv as they appear in scv loop
Definition fv1ib_geom.h:488
size_t num_boundary_subsets()
number of registered boundary subsets
Definition fv1ib_geom.h:507
LagrangeP1< ref_elem_type > local_shape_fct_set_type
type of Shape function used
Definition fv1ib_geom.h:109
MathVector< worldDim > m_vvGloMid[dim+1][maxMid]
Definition fv1ib_geom.h:561
void clear_boundary_subsets()
reset all boundary subsets
Definition fv1ib_geom.h:504
fv1_traits< ref_elem_type, TWorldDim > traits
used traits
Definition fv1ib_geom.h:78
const MathVector< dim > * scv_local_ips() const
returns all ips of scv as they appear in scv loop
Definition fv1ib_geom.h:485
const std::vector< BF > & bf(int si) const
returns reference to vector of boundary faces for subsetIndex
Definition fv1ib_geom.h:539
const SCV & scv(size_t i) const
const access to SubControlVolume number i
Definition fv1ib_geom.h:465
static const size_t nip
number of integration points
Definition fv1ib_geom.h:115
traits::scv_type scv_type
type of SubControlVolume
Definition fv1ib_geom.h:102
SCVF m_vSCVF[numSCVF]
SubControlVolumeFaces.
Definition fv1ib_geom.h:564
ReferenceMapping< ref_elem_type, worldDim > m_mapping
Reference Mapping.
Definition fv1ib_geom.h:570
const BF & bf(int si, size_t i) const
returns the boundary face i for subsetIndex
Definition fv1ib_geom.h:529
static const size_t numSCV
number of SubControlVolumes
Definition fv1ib_geom.h:98
void add_boundary_subset(int subsetIndex)
add subset that is interpreted as boundary subset.
Definition fv1ib_geom.h:498
size_t num_bf(int si) const
number of boundary faces on subset 'subsetIndex'
Definition fv1ib_geom.h:520
static const int order
order
Definition fv1ib_geom.h:95
FV1IBGeometry()
construct object and initialize local values and sizes
Definition fv1ib_geom_impl.h:552
static const size_t nsh
number of shape functions
Definition fv1ib_geom.h:112
const local_shape_fct_set_type & m_rTrialSpace
Shape function set.
Definition fv1ib_geom.h:576
static const bool usesHangingNodes
Hanging node flag: this Geometry does not support hanging nodes.
Definition fv1ib_geom.h:88
reference_element_traits< TElem >::reference_element_type ref_elem_type
type of reference element
Definition fv1ib_geom.h:75
void update_local_data()
update local data
Definition fv1ib_geom_impl.h:561
MathVector< dim > m_vvLocMid[dim+1][maxMid]
local and global geom object midpoints for each dimension
Definition fv1ib_geom.h:560
size_t num_bf() const
number of all boundary faces
Definition fv1ib_geom.h:510
size_t num_scvf_ips() const
returns number of all scvf ips
Definition fv1ib_geom.h:473
MathVector< worldDim > vSCVF_PosOnEdge[numSCVF]
new members for adaptions
Definition fv1ib_geom.h:580
size_t num_scv() const
number of SubControlVolumes
Definition fv1ib_geom.h:462
const MathVector< dim > * scvf_local_ips() const
returns all ips of scvf as they appear in scv loop
Definition fv1ib_geom.h:476
void update(GridObject *elem, const MathVector< worldDim > *vCornerCoords, const ISubsetHandler *ish=NULL)
update data for given element
Definition fv1ib_geom_impl.h:643
void adapt(GridObject *elem, const MathVector< TWorldDim > *vCornerCoords, const ISubsetHandler *ish=NULL)
adapt integration points for elements cut by the inner boundary
Definition fv1ib_geom_impl.h:47
const ref_elem_type & m_rRefElem
Reference Element.
Definition fv1ib_geom.h:573
void update_boundary_faces(GridObject *elem, const MathVector< worldDim > *vCornerCoords, const ISubsetHandler *ish=NULL)
update boundary data for given element
Definition fv1ib_geom_impl.h:741
static const int dim
dimension of reference element
Definition fv1ib_geom.h:82
static const bool staticLocalData
flag indicating if local data may change
Definition fv1ib_geom.h:91
MathVector< worldDim > m_MidPoint
Definition fv1ib_geom.h:581
MathVector< dim > m_vLocSCVF_IP[numSCVF]
Definition fv1ib_geom.h:494
std::map< int, std::vector< BF > > m_mapVectorBF
Definition fv1ib_geom.h:548
size_t num_scv_ips() const
returns number of all scv ips
Definition fv1ib_geom.h:482
const MathVector< worldDim > * scvf_global_ips() const
returns all ips of scvf as they appear in scv loop
Definition fv1ib_geom.h:479
static const size_t numSCVF
number of SubControlVolumeFaces
Definition fv1ib_geom.h:105
size_t num_scvf() const
number of SubControlVolumeFaces
Definition fv1ib_geom.h:453
static const int worldDim
dimension of world
Definition fv1ib_geom.h:85
const SCVF & scvf(size_t i) const
const access to SubControlVolumeFace number i
Definition fv1ib_geom.h:456
SCV m_vSCV[numSCV]
SubControlVolumes.
Definition fv1ib_geom.h:567
void adapt_normals(GridObject *elem, const MathVector< TWorldDim > *vCornerCoords, const ISubsetHandler *ish=NULL)
adapt normals for elements cut by the inner boundary
Definition fv1ib_geom_impl.h:120
TElem * m_pElem
pointer to current element
Definition fv1ib_geom.h:553
void adapt_integration_points(GridObject *elem, const MathVector< TWorldDim > *vCornerCoords, const ISubsetHandler *ish=NULL)
adapt integration points for elements cut by the inner boundary
Definition fv1ib_geom_impl.h:111
std::vector< BF > m_vEmptyVectorBF
Definition fv1ib_geom.h:549
const MathVector< worldDim > * corners() const
get vector of the global coordinates of corners for current element
Definition fv1ib_geom.h:450
TElem elem_type
type of element
Definition fv1ib_geom.h:72
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
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
traits for reference elements
Definition reference_element_traits.h:48