Loading [MathJax]/extensions/tex2jax.js
ug4
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
volume_util_impl.hpp
Go to the documentation of this file.
1/*
2 * Copyright (c) 2010-2015: G-CSC, Goethe University Frankfurt
3 * Author: Sebastian Reiter
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__LIB_GRID__VOLUME_UTIL_IMPL__
34#define __H__LIB_GRID__VOLUME_UTIL_IMPL__
35
36#include <algorithm>
37#include <vector>
38#include "lib_grid/lg_base.h"
44
45namespace ug
46{
48template <class TAAPos>
49bool
50ContainsPoint(Volume* vol, const vector3& p, TAAPos aaPos)
51{
52// iterate over face descriptors of the sides and check whether the point
53// lies inside or outside
55 vector3 n, dir;
56
57// to minimize rounding errors we'll compare against a small-constant which is
58// relative to the first edge of the examined volume.
60 vol->edge_desc(0, ed);
61 number len = EdgeLength(&ed, aaPos);
62
63 // the constant should be relative to the same geometric measure as what it is
64 // compared against later on, i.e. length*area, since otherwise problems arise
65 // with geometries scaled to very small extensions;
66 // which is why I changed sqrt(lenSq) to lenSq^1.5 (mbreit, 2015-05-11)
67 const number locSmall = len * len * len * SMALL;
68
69 for(size_t i = 0; i < vol->num_faces(); ++i){
70 vol->face_desc(i, fd);
71 CalculateNormalNoNormalize(n, &fd, aaPos);
72 VecSubtract(dir, aaPos[fd.vertex(0)], p);
73
74 if(VecDot(dir, n) < -locSmall)
75 return false;
76 }
77 return true;
78}
79
81// PointIsInsideTetrahedron
82inline bool
85{
86 return PointIsInsideTetrahedron(v, aaPos[tet->vertex(0)], aaPos[tet->vertex(1)],
87 aaPos[tet->vertex(2)], aaPos[tet->vertex(3)]);
88}
89
91// CalculateCenter
92template<class TVertexPositionAttachmentAccessor>
93typename TVertexPositionAttachmentAccessor::ValueType
94CalculateCenter(const VolumeVertices* vol, TVertexPositionAttachmentAccessor& aaPosVRT)
95{
96 typename TVertexPositionAttachmentAccessor::ValueType v;
97// init v with 0.
98 VecSet(v, 0);
99
100 size_t numVrts = vol->num_vertices();
102
103// sum up
104 for(size_t i = 0; i < numVrts; ++i)
105 {
106 VecAdd(v, v, aaPosVRT[vrts[i]]);
107 }
108
109// average
110 if(numVrts > 0)
111 VecScale(v, v, 1./(number)numVrts);
112
113 return v;
114}
115
117template<class TAAPosVRT, class TAAWeightVRT>
118UG_API
119typename TAAPosVRT::ValueType
120CalculateCenter(const VolumeVertices* vol, TAAPosVRT& aaPos, TAAWeightVRT& aaWeight)
121{
122 typename TAAPosVRT::ValueType v;
123 typedef typename TAAWeightVRT::ValueType weight_t;
124
125// init v with 0.
126 VecSet(v, 0);
127
128 size_t numVrts = vol->num_vertices();
130
131// sum up
132 weight_t totalWeight = 0;
133 for(size_t i = 0; i < numVrts; ++i)
134 {
135 weight_t w = aaWeight[vrts[i]];
136 VecScaleAppend(v, w, aaPos[vrts[i]]);
137 totalWeight += w;
138 }
139
140// average
141 if(totalWeight != 0)
142 VecScale(v, v, 1./(number)totalWeight);
143
144 return v;
145}
146
147
149template <class TElem>
151 public:
152 CmpVrtsByHash(TElem* e) : m_e(e) {};
153 bool operator () (int i0, int i1) {
154 return m_e->vertex(i0)->get_hash_value() < m_e->vertex(i1)->get_hash_value();
155 }
156 private:
157 TElem* m_e;
158};
159
160
161template <class TVolIter>
163 Grid& grid,
164 TVolIter volsBegin,
165 TVolIter volsEnd)
166{
167 using namespace std;
168
169// first we'll collect all quadrilaterals that are connected to selected
170// volumes. Avoid 'grid.mark' here, since it may be used by 'grid.associated_elements'
172 vector<Face*> quads;
173
174 for(TVolIter iv = volsBegin; iv != volsEnd; ++iv){
175 grid.associated_elements(faces, *iv);
176 for_each_in_vec(Face* f, faces){
177 if(f->num_vertices() == 4)
178 quads.push_back(f);
179 }end_for;
180 }
181
182// remove double entries
183 grid.begin_marking();
184 size_t offset = 0;
185
186 for(size_t i = 0; i + offset < quads.size();){
187 if(offset > 0)
188 quads[i] = quads[i + offset];
189
190 if(!grid.is_marked(quads[i])){
191 grid.mark(quads[i]);
192 ++i;
193 }
194 else{
195 ++offset;
196 }
197 }
198
199 if(offset > 0)
200 quads.resize(quads.size() - offset);
201
202 grid.end_marking();
203
204 for_each_in_vec(Face* f, quads){
205//todo in a parallel environment, global id's should be compared here
206 CmpVrtsByHash<Face> cmp(f);
207 // get the smallest vertex of the face
208 int smallest = 0;
209 for(int i = 1; i < 4; ++i){
210 if(cmp(i, smallest))
211 smallest = i;
212 }
213
214 int i0 = smallest;
215 int i1 = (smallest + 1) % 4;
216 int i2 = (smallest + 2) % 4;
217 int i3 = (smallest + 3) % 4;
218 grid.create<Triangle>(TriangleDescriptor(f->vertex(i0), f->vertex(i1), f->vertex(i2)), f);
219 grid.create<Triangle>(TriangleDescriptor(f->vertex(i2), f->vertex(i3), f->vertex(i0)), f);
220 }end_for;
221
222
223// now convert the given volume-elements
228 HEX_RULES_MAX_NUM_CONVERT_TO_TETS_INDS_OUT__considered_to_be_highest_among_prism_pyra_and_hex);
229 static const int arrayLen = hex_rules::MAX_NUM_CONVERT_TO_TETS_INDS_OUT;
230 int inds[arrayLen];
231
232 vector<Volume*> volsToErase;
233 for(TVolIter iv = volsBegin; iv != volsEnd; ++iv){
234 Volume* vol = *iv;
235 const ReferenceObjectID roid = vol->reference_object_id();
236 CmpVrtsByHash<Volume> cmp(vol);
237 size_t numEntries = 0;
238
239 switch(roid){
240 case ROID_PYRAMID:
241 numEntries = pyra_rules::ConvertToTetrahedra(inds, cmp);
242 break;
243
244 case ROID_PRISM:
245 numEntries = prism_rules::ConvertToTetrahedra(inds, cmp);
246 break;
247
248 case ROID_HEXAHEDRON:
249 UG_THROW("ConvertToTetrahedra for hexahedra not yet implemented!");
250 break;
251 default:
252 break;
253 }
254
255 if(numEntries > 0){
256 volsToErase.push_back(vol);
257 size_t i = 0;
259 while(i < numEntries){
260 int goid = inds[i++];
262 "Only tetrahedra may result from ConvertToTetrahedra");
263 int i0 = inds[i++];
264 int i1 = inds[i++];
265 int i2 = inds[i++];
266 int i3 = inds[i++];
267 grid.create<Tetrahedron>(
268 TetrahedronDescriptor(vrts[i0], vrts[i1], vrts[i2], vrts[i3]),
269 vol);
270 }
271 }
272 }
273
274// finally erase all unused volumes and quadrilaterals
275 for_each_in_vec(Volume* v, volsToErase){
276 grid.erase(v);
277 }end_for;
278
280 for_each_in_vec(Face* f, quads){
281 grid.associated_elements(assVols, f);
282 if(assVols.empty())
283 grid.erase(f);
284 }end_for;
285}
286
287}// end of namespace
288
289#endif
parameterString p
Can be used to compare vertices of their grids through their hash-value.
Definition volume_util_impl.hpp:150
TElem * m_e
Definition volume_util_impl.hpp:157
CmpVrtsByHash(TElem *e)
Definition volume_util_impl.hpp:152
bool operator()(int i0, int i1)
Definition volume_util_impl.hpp:153
Can be used to store information about an edge and to construct an edge.
Definition grid_base_objects.h:464
Can be queried for the edges and vertices of a face.
Definition grid_base_objects.h:684
virtual Vertex * vertex(size_t index) const
Definition grid_base_objects.h:701
Faces are 2-dimensional objects.
Definition grid_base_objects.h:510
virtual size_t num_vertices() const
Definition grid_base_objects.h:488
virtual Vertex * vertex(size_t index) const
Definition grid_base_objects.h:486
Manages the elements of a grid and their interconnection.
Definition grid.h:132
void end_marking()
ends a marking sequence. Call this method when you're done with marking.
Definition grid.cpp:1285
bool is_marked(GridObject *obj) const
returns true if the object is marked, false if not.
Definition grid_impl.hpp:843
void mark(GridObject *obj)
marks the object. Calls are only valid between calls to Grid::begin_marking and Grid::end_marking.
Definition grid_impl.hpp:773
void erase(GridObject *geomObj)
Definition grid.cpp:459
void begin_marking()
begin marking.
Definition grid.cpp:1262
void associated_elements(traits< Vertex >::secure_container &elemsOut, TElem *e)
Puts all elements of type TAss which are contained in 'e' or which contain 'e' into elemsOut.
Definition grid_impl.hpp:466
geometry_traits< TGeomObj >::iterator create(GridObject *pParent=NULL)
create a custom element.
Definition grid_impl.hpp:69
Vertex *const * ConstVertexArray
Definition grid_base_objects.h:319
only used to initialize a tetrahedron. for all other tasks you should use VolumeDescripor.
Definition grid_objects_3d.h:68
the most simple volume-element.
Definition grid_objects_3d.h:91
virtual Vertex * vertex(size_t index) const
Definition grid_objects_3d.h:106
only used to initialize a triangle. for all other tasks you should use FaceDescriptor.
Definition grid_objects_2d.h:78
the most simple form of a face
Definition grid_objects_2d.h:174
Volumes are 3-dimensional objects.
Definition grid_base_objects.h:754
virtual ReferenceObjectID reference_object_id() const
Definition grid_base_objects.h:927
virtual FaceDescriptor face_desc(int index) const
Definition grid_base_objects.h:784
virtual EdgeDescriptor edge_desc(int index) const
Definition grid_base_objects.h:780
virtual uint num_faces() const
Definition grid_base_objects.h:786
holds the vertices of a Volume or a VolumeDescriptor
Definition grid_base_objects.h:723
virtual ConstVertexArray vertices() const
Definition grid_base_objects.h:728
virtual size_t num_vertices() const
Definition grid_base_objects.h:729
int CalculateNormalNoNormalize(vector3 &vNormOut, Grid &grid, Edge *e, Grid::VertexAttachmentAccessor< APosition > &aaPos, Grid::FaceAttachmentAccessor< ANormal > *paaNormFACE)
Calculates the normal of the given edge.
Definition edge_util.cpp:358
UG_API number EdgeLength(const EdgeVertices *e, TAAPosVRT &aaPos)
Calculates the length of the given edge.
Definition edge_util_impl.hpp:80
UG_API bool ContainsPoint(const EdgeVertices *e, const vector_t &p, TAAPos aaPos)
Returns true if the given point lies on the given edge.
Definition edge_util_impl.hpp:513
void ConvertToTetrahedra(Grid &grid, TVolIter volsBegin, TVolIter volsEnd)
Converts all volumes in the given range to tetrahedra.
Definition volume_util_impl.hpp:162
#define UG_THROW(msg)
Definition error.h:57
#define UG_STATIC_ASSERT(expr, msg)
Checks an expression at compile-time and raises a compile-error if the expression equals 0.
Definition static_assert.h:63
#define UG_API
Definition ug_config.h:65
#define UG_COND_THROW(cond, msg)
UG_COND_THROW(cond, msg) : performs a UG_THROW(msg) if cond == true.
Definition error.h:61
double number
Definition types.h:124
void CalculateCenter(vector_t &centerOut, const vector_t *pointSet, size_t numPoints)
calculates the center of a point-set
Definition math_util_impl.hpp:98
bool PointIsInsideTetrahedron(const vector_t &v, const vector_t &v0, const vector_t &v1, const vector_t &v2, const vector_t &v3)
Returns true if the point lies inside or on the boundary of a tetrahedron.
Definition math_util_impl.hpp:1110
void VecScaleAppend(vector_t &vOut, typename vector_t::value_type s1, const vector_t &v1)
Scales a Vector and adds it to a second vector.
Definition math_vector_functions_common_impl.hpp:126
void VecSet(vector_t &vInOut, typename vector_t::value_type s)
Set each vector component to scalar (componentwise)
Definition math_vector_functions_common_impl.hpp:539
void VecSubtract(vector_t &vOut, const vector_t &v1, const vector_t &v2)
subtracts v2 from v1 and stores the result in a vOut
Definition math_vector_functions_common_impl.hpp:226
void VecAdd(vector_t &vOut, const vector_t &v1, const vector_t &v2)
adds two MathVector<N>s and stores the result in a third one
Definition math_vector_functions_common_impl.hpp:185
void VecScale(vector_t &vOut, const vector_t &v, typename vector_t::value_type s)
scales a MathVector<N>
Definition math_vector_functions_common_impl.hpp:252
vector_t::value_type VecDot(const vector_t &v1, const vector_t &v2)
returns the dot-product of two vector_ts
Definition math_vector_functions_common_impl.hpp:385
Definition smart_pointer.h:814
const int MAX_NUM_CONVERT_TO_TETS_INDS_OUT
Definition hexahedron_rules.h:51
int ConvertToTetrahedra(int *newIndsOut, TCmp cmp)
fills an array of integers describing tetrahedra that shall replace the prism
Definition prism_rules_impl.h:43
const int MAX_NUM_CONVERT_TO_TETS_INDS_OUT
Definition prism_rules.h:51
int ConvertToTetrahedra(int *newIndsOut, TCmp cmp)
fills an array of integers describing tetrahedra that shall replace the pyramid
Definition pyramid_rules_impl.h:43
const int MAX_NUM_CONVERT_TO_TETS_INDS_OUT
Definition pyramid_rules.h:52
the ug namespace
@ GOID_TETRAHEDRON
Definition grid_object_ids.h:43
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_PRISM
Definition grid_base_objects.h:82
@ ROID_HEXAHEDRON
Definition grid_base_objects.h:81
const number SMALL
Definition math_constants.h:41
PointerConstArray< Face * > secure_container
Definition grid.h:146
#define for_each_in_vec(_vfeDecl, _vfeVec)
Allows iteration over all members of an std::vector compatible type.
Definition vec_for_each.h:52
#define end_for
Allows iteration over all members of an std::vector compatible type.
Definition vec_for_each.h:56