Loading [MathJax]/extensions/tex2jax.js
ug4
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
fv1ib_geom_impl.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/*
34 * copy of 'fv1_geom.h' with two additional methods:
35 * 'adapt_normals()' and 'adapt_integration_points()'
36 */
37
38
40 #include "fv1ib_geom.h"
44
45template <typename TElem, int TWorldDim>
47 adapt(GridObject* elem, const MathVector<TWorldDim>* vCornerCoords,
48 const ISubsetHandler* ish)
49{
50
51 for(size_t ip = 0; ip < num_scvf(); ++ip)
52 UG_LOG("#### vSCVF_PosOnEdge = " << vSCVF_PosOnEdge[ip] << "\n");
53
54// compute global informations for scvf
55 for(size_t i = 0; i < num_scvf(); ++i)
56 {
57 UG_LOG("ip's: " << m_vSCVF[i].globalIP << "\n");
58 UG_LOG("normals: " << m_vSCVF[i].Normal << "\n");
59
60 }
61
62// compute global informations for scvf
63 for(size_t i = 0; i < num_scvf(); ++i)
64 {
65 //VecScaleAdd(m_vSCVF[i].globalIP, vSCVF_PosOnEdge[i], 0.5, m_MidPoint, 0.5);
66 for(size_t d = 0; d < m_vSCVF[i].globalIP.size(); ++d)
67 m_vSCVF[i].globalIP[d] = 0.5*(vSCVF_PosOnEdge[i][d] + m_MidPoint[d]);
68
69 //traits::NormalOnSCVF(m_vSCVF[i].Normal, vSCVF_PosOnEdge[i], m_MidPoint);
70 // 'NormalOnSCVF' braucht Felder im 2. und 3. Parameter
71 //traits::NormalOnSCVF(m_vSCVF[i].Normal, vSCVF_PosOnEdge, vSCVF_PosOnEdge);
72 MathVector<worldDim> buffer_comp;
73 VecSubtract(buffer_comp, vSCVF_PosOnEdge[i], m_MidPoint);
74
75 MathVector<worldDim> buffer_copy;
76 buffer_copy[0] = -buffer_comp[1];
77 buffer_copy[1] = buffer_comp[0];
78
79 // check the orientation of the newly computet Normal:
80 // -> has to be in direction from -> to, i.e. VecDot(Normal_alt, Normal_neu) > 0
81 if ( VecDot(m_vSCVF[i].Normal, buffer_copy) > 0 )
82 {
83 UG_LOG("Positiv: normals: " << m_vSCVF[i].Normal << "\n");
84 UG_LOG("Positiv: normals: " << buffer_copy << "\n");
85
86 m_vSCVF[i].Normal[0] = buffer_copy[0];
87 m_vSCVF[i].Normal[1] = buffer_copy[1];
88 }
89 else
90 {
91 UG_LOG("Negativ: normals: " << m_vSCVF[i].Normal << "\n");
92 UG_LOG("Negativ: normals: " << buffer_copy << "\n");
93
94 m_vSCVF[i].Normal[0] = -buffer_copy[0];
95 m_vSCVF[i].Normal[1] = -buffer_copy[1];
96 }
97
98
99 UG_LOG("***ip's: " << m_vSCVF[i].globalIP << "\n");
100 UG_LOG("***normals: " << m_vSCVF[i].Normal << "\n");
101 UG_LOG("***buffer_copy: " << buffer_copy << "\n");
102
103 }
104
105 //UG_THROW("### ok done...\n");
106}
107
108
109template <typename TElem, int TWorldDim>
112 const ISubsetHandler* ish)
113{
114 UG_THROW("adapt_integration_points: ok done...\n");
115
116}
117
118template <typename TElem, int TWorldDim>
120adapt_normals(GridObject* elem, const MathVector<TWorldDim>* vCornerCoords,
121 const ISubsetHandler* ish)
122{
125
126
127 for(size_t i = 0; i < num_scvf(); ++i)
128 {
129 UG_LOG("ip's: " << m_vSCVF[i].globalIP << "\n");
130 UG_LOG("normals: " << m_vSCVF[i].Normal << "\n");
131
132 }
133
134 UG_THROW("stop hier...\n");
135
136
137}
138
139
140
146template <int dim, typename TRefElem, int maxMid>
147static void ComputeMidPoints(const TRefElem& rRefElem,
148 const MathVector<dim> vCorner[],
149 MathVector<dim> vvMid[][maxMid])
150{
151// compute local midpoints for all geometric objects with 0 < d <= dim
152 for(int d = 1; d <= dim; ++d)
153 {
154 // loop geometric objects of dimension d
155 for(size_t i = 0; i < rRefElem.num(d); ++i)
156 {
157 // set first node
158 const size_t coID0 = rRefElem.id(d, i, 0, 0);
159 vvMid[d][i] = vCorner[coID0];
160
161 // add corner coordinates of the corners of the geometric object
162 for(size_t j = 1; j < rRefElem.num(d, i, 0); ++j)
163 {
164 const size_t coID = rRefElem.id(d, i, 0, j);
165 vvMid[d][i] += vCorner[coID];
166 }
167
168 // scale for correct averaging
169 vvMid[d][i] *= 1./(rRefElem.num(d, i, 0));
170 }
171 }
172
173 // for PYRAMIDS: add midpoints of imaginary faces, edges and volumes
174 // resulting from the division into two tetrahedra alongside x==y
175 if (rRefElem.roid() == ROID_PYRAMID)
176 {
177 // diagonal 2->0, diagonal 0->2
178 VecScaleAdd(vvMid[1][rRefElem.num(1)], 0.5, vCorner[2], 0.5, vCorner[0]);
179 VecScaleAdd(vvMid[1][rRefElem.num(1)+1], 0.5, vCorner[0], 0.5, vCorner[2]);
180
181 // subface 0,1,2; subface 0,2,3; face 0,4,2; face 0,2,4
182 vvMid[2][rRefElem.num(2)] = vCorner[0];
183 vvMid[2][rRefElem.num(2)] += vCorner[1];
184 vvMid[2][rRefElem.num(2)] += vCorner[2];
185 vvMid[2][rRefElem.num(2)] *= 1.0/3.0;
186
187 vvMid[2][rRefElem.num(2)+1] = vCorner[0];
188 vvMid[2][rRefElem.num(2)+1] += vCorner[2];
189 vvMid[2][rRefElem.num(2)+1] += vCorner[3];
190 vvMid[2][rRefElem.num(2)+1] *= 1.0/3.0;
191
192 vvMid[2][rRefElem.num(2)+2] = vCorner[0];
193 vvMid[2][rRefElem.num(2)+2] += vCorner[4];
194 vvMid[2][rRefElem.num(2)+2] += vCorner[2];
195 vvMid[2][rRefElem.num(2)+2] *= 1.0/3.0;
196
197 vvMid[2][rRefElem.num(2)+3] = vCorner[0];
198 vvMid[2][rRefElem.num(2)+3] += vCorner[2];
199 vvMid[2][rRefElem.num(2)+3] += vCorner[4];
200 vvMid[2][rRefElem.num(2)+3] *= 1.0/3.0;
201
202 // subvolume 0,1,2,4; subvolume 0,2,3,4
203
204 vvMid[3][rRefElem.num(3)] = vCorner[0];
205 vvMid[3][rRefElem.num(3)] += vCorner[1];
206 vvMid[3][rRefElem.num(3)] += vCorner[2];
207 vvMid[3][rRefElem.num(3)] += vCorner[4];
208 vvMid[3][rRefElem.num(3)] *= 0.25;
209
210 vvMid[3][rRefElem.num(3)+1] = vCorner[0];
211 vvMid[3][rRefElem.num(3)+1] += vCorner[2];
212 vvMid[3][rRefElem.num(3)+1] += vCorner[3];
213 vvMid[3][rRefElem.num(3)+1] += vCorner[4];
214 vvMid[3][rRefElem.num(3)+1] *= 0.25;
215 }
216}
217
221template <typename TRefElem>
222static void ComputeSCVFMidID(const TRefElem& rRefElem,
223 MidID vMidID[], int i)
224{
225 static const int dim = TRefElem::dim;
226
227 if (rRefElem.roid() != ROID_PYRAMID)
228 {
229 // set mid ids
230 {
231 // start at edge midpoint
232 vMidID[0] = MidID(1,i);
233
234 // loop up dimension
235 if(dim == 2)
236 {
237 vMidID[1] = MidID(dim, 0); // center of element
238 }
239 else if (dim == 3)
240 {
241 vMidID[1] = MidID(2, rRefElem.id(1, i, 2, 0)); // side 0
242 vMidID[2] = MidID(dim, 0); // center of element
243 vMidID[3] = MidID(2, rRefElem.id(1, i, 2, 1)); // side 1
244 }
245 }
246 }
247 // pyramid here
248 else
249 {
250 switch (i)
251 {
252 // scvf of edge 0
253 case 0: vMidID[0] = MidID(1,0); // edge 0
254 vMidID[1] = MidID(2,5); // subface 0/0
255 vMidID[2] = MidID(3,1); // subvolume 0/0
256 vMidID[3] = MidID(2,1); // face 1
257 break;
258 // scvf of edge 1
259 case 1: vMidID[0] = MidID(1,1); // edge 1
260 vMidID[1] = MidID(2,5); // subface 0/0
261 vMidID[2] = MidID(3,1); // subvolume 0/0
262 vMidID[3] = MidID(2,2); // face 2
263 break;
264 // scvf of diagonal 2->0
265 case 2: vMidID[0] = MidID(1,8); // diagonal 2->0
266 vMidID[1] = MidID(2,5); // subface 0/0
267 vMidID[2] = MidID(3,1); // subvolume 0/0
268 vMidID[3] = MidID(2,7); // face 0,4,2
269 break;
270 // scvf of edge 4 in subvolume 0/0
271 case 3: vMidID[0] = MidID(1,4); // edge 4
272 vMidID[1] = MidID(2,1); // face 1
273 vMidID[2] = MidID(3,1); // subvolume 0/0
274 vMidID[3] = MidID(2,7); // face 0,4,2
275 break;
276 // scvf of edge 5
277 case 4: vMidID[0] = MidID(1,5); // edge 5
278 vMidID[1] = MidID(2,2); // face 2
279 vMidID[2] = MidID(3,1); // subvolume 0/0
280 vMidID[3] = MidID(2,1); // face 1
281 break;
282 // scvf of edge 6 in subvolume 0/0
283 case 5: vMidID[0] = MidID(1,6); // edge 6
284 vMidID[1] = MidID(2,7); // face 0,4,2
285 vMidID[2] = MidID(3,1); // subvolume 0/0
286 vMidID[3] = MidID(2,2); // face 2
287 break;
288 // edge 0->2
289 case 6: vMidID[0] = MidID(1,9); // edge 0->2
290 vMidID[1] = MidID(2,6); // subface 1/0
291 vMidID[2] = MidID(3,2); // subvolume 1/0
292 vMidID[3] = MidID(2,8); // face 0,2,4
293 break;
294 // scvf of edge 2
295 case 7: vMidID[0] = MidID(1,2); // edge 2
296 vMidID[1] = MidID(2,6); // subface 1/0
297 vMidID[2] = MidID(3,2); // subvolume 1/0
298 vMidID[3] = MidID(2,3); // face 3
299 break;
300 // scvf of edge 3
301 case 8: vMidID[0] = MidID(1,3); // edge 3
302 vMidID[1] = MidID(2,6); // subface 1/0
303 vMidID[2] = MidID(3,2); // subvolume 1/0
304 vMidID[3] = MidID(2,4); // face 4
305 break;
306 // scvf of edge 4 in subvolume 1/0
307 case 9: vMidID[0] = MidID(1,4); // edge 4
308 vMidID[1] = MidID(2,8); // face 0,2,4
309 vMidID[2] = MidID(3,2); // subvolume 1/0
310 vMidID[3] = MidID(2,4); // face 4
311 break;
312 // scvf of edge 6 in subvolume 1/0
313 case 10:vMidID[0] = MidID(1,6); // edge 6
314 vMidID[1] = MidID(2,3); // face 3
315 vMidID[2] = MidID(3,2); // subvolume 1/0
316 vMidID[1] = MidID(2,8); // face 0,2,4
317 break;
318 // scvf of edge 7
319 case 11:vMidID[0] = MidID(1,7); // edge 7
320 vMidID[3] = MidID(2,4); // face 4
321 vMidID[2] = MidID(3,2); // subvolume 1/0
322 vMidID[1] = MidID(2,3); // face 3
323 break;
324 default:UG_THROW("Pyramid only has 12 SCVFs (no. 0-11), but requested no. " << i << ".");
325 break;
326 }
327 }
328}
329
333template <typename TRefElem>
334static void ComputeSCVMidID(const TRefElem& rRefElem,
335 MidID vMidID[], int i)
336{
337 static const int dim = TRefElem::dim;
338
339 if (rRefElem.roid() != ROID_PYRAMID)
340 {
341 if(dim == 1)
342 {
343 vMidID[0] = MidID(0, i); // set node as corner of scv
344 vMidID[1] = MidID(dim, 0); // center of element
345 }
346 else if(dim == 2)
347 {
348 vMidID[0] = MidID(0, i); // set node as corner of scv
349 vMidID[1] = MidID(1, rRefElem.id(0, i, 1, 0)); // edge 1
350 vMidID[2] = MidID(dim, 0); // center of element
351 vMidID[3] = MidID(1, rRefElem.id(0, i, 1, 1)); // edge 2
352 }
353 else if(dim == 3)
354 {
355 vMidID[0] = MidID(0, i); // set node as corner of scv
356 vMidID[1] = MidID(1, rRefElem.id(0, i, 1, 1)); // edge 1
357 vMidID[2] = MidID(2, rRefElem.id(0, i, 2, 0)); // face 0
358 vMidID[3] = MidID(1, rRefElem.id(0, i, 1, 0)); // edge 0
359 vMidID[4] = MidID(1, rRefElem.id(0, i, 1, 2)); // edge 2
360 vMidID[5] = MidID(2, rRefElem.id(0, i, 2, 2)); // face 2
361 vMidID[6] = MidID(dim, 0); // center of element
362 vMidID[7] = MidID(2, rRefElem.id(0, i, 2, 1)); // face 1
363 }
364 else {UG_THROW("Dimension higher that 3 not implemented.");}
365 }
366 // pyramid here
367 else
368 {
369 switch (i)
370 {
371 // scv of corner 0 in subvolume 0/0
372 case 0: vMidID[0] = MidID(0,0); // corner 0
373 vMidID[1] = MidID(1,0); // edge 0
374 vMidID[2] = MidID(2,5); // subface 0/0
375 vMidID[3] = MidID(1,8); // edge 2->0
376 vMidID[4] = MidID(1,4); // edge 4
377 vMidID[5] = MidID(2,1); // face 1
378 vMidID[6] = MidID(3,1); // subvolume 0/0
379 vMidID[7] = MidID(2,7); // face 0,4,2
380 break;
381 // scv of corner 1
382 case 1: vMidID[0] = MidID(0,1); // corner 1
383 vMidID[1] = MidID(1,1); // edge 1
384 vMidID[2] = MidID(2,5); // subface 0/0
385 vMidID[3] = MidID(1,0); // edge 0
386 vMidID[4] = MidID(1,5); // edge 5
387 vMidID[5] = MidID(2,2); // face 2
388 vMidID[6] = MidID(3,1); // subvolume 0/0
389 vMidID[7] = MidID(2,1); // face 1
390 break;
391 // scv of corner 2 in subvolume 0/0
392 case 2: vMidID[0] = MidID(0,2); // corner 2
393 vMidID[1] = MidID(1,8); // edge 2->0
394 vMidID[2] = MidID(2,5); // subface 0/0
395 vMidID[3] = MidID(1,1); // edge 1
396 vMidID[4] = MidID(1,6); // edge 6
397 vMidID[5] = MidID(2,7); // face 0,4,2
398 vMidID[6] = MidID(3,1); // subvolume 0/0
399 vMidID[7] = MidID(2,2); // face 2
400 break;
401 // scv of corner 4 in subvolume 0/0
402 case 3: vMidID[0] = MidID(0,4); // corner 4
403 vMidID[1] = MidID(1,5); // edge 5
404 vMidID[2] = MidID(2,1); // face 1
405 vMidID[3] = MidID(1,4); // edge 4
406 vMidID[4] = MidID(1,6); // edge 6
407 vMidID[5] = MidID(2,2); // face 2
408 vMidID[6] = MidID(3,1); // subvolume 0/0
409 vMidID[7] = MidID(2,7); // face 0,4,2
410 break;
411 // scv of corner 0 in subvolume 1/0
412 case 4: vMidID[0] = MidID(0,0); // corner 0
413 vMidID[1] = MidID(1,9); // edge 0->2
414 vMidID[2] = MidID(2,6); // subface 1/0
415 vMidID[3] = MidID(1,3); // edge 3
416 vMidID[4] = MidID(1,4); // edge 4
417 vMidID[5] = MidID(2,8); // face 0,2,4
418 vMidID[6] = MidID(3,2); // subvolume 1/0
419 vMidID[7] = MidID(2,4); // face 4
420 break;
421 // scv of corner 2 in subvolume 1/0
422 case 5: vMidID[0] = MidID(0,2); // corner 2
423 vMidID[1] = MidID(1,2); // edge 2
424 vMidID[2] = MidID(2,6); // subface 1/0
425 vMidID[3] = MidID(1,9); // edge 0->2
426 vMidID[4] = MidID(1,6); // edge 6
427 vMidID[5] = MidID(2,3); // face 3
428 vMidID[6] = MidID(3,2); // subvolume 1/0
429 vMidID[7] = MidID(2,8); // face 0,2,4
430 break;
431 // scv of corner 3
432 case 6: vMidID[0] = MidID(0,3); // corner 3
433 vMidID[1] = MidID(1,3); // edge 3
434 vMidID[2] = MidID(2,6); // subface 1/0
435 vMidID[3] = MidID(1,2); // edge 2
436 vMidID[4] = MidID(1,7); // edge 7
437 vMidID[5] = MidID(2,4); // face 4
438 vMidID[6] = MidID(3,2); // subvolume 1/0
439 vMidID[7] = MidID(2,3); // face 3
440 break;
441 // scv of corner 4 in subvolume 1/0
442 case 7: vMidID[0] = MidID(0,4); // corner 4
443 vMidID[1] = MidID(1,6); // edge 6
444 vMidID[2] = MidID(2,8); // face 0,2,4
445 vMidID[3] = MidID(1,4); // edge 4
446 vMidID[4] = MidID(1,7); // edge 7
447 vMidID[5] = MidID(2,3); // face 3
448 vMidID[6] = MidID(3,2); // subvolume 1/0
449 vMidID[7] = MidID(2,4); // face 4
450 break;
451 default:UG_THROW("Pyramid only has 8 SCVs (no. 0-7), but requested no. " << i << ".");
452 break;
453 }
454 }
455}
456
460template <typename TRefElem>
461static void ComputeBFMidID(const TRefElem& rRefElem, int side,
462 MidID vMidID[], int co)
463{
464 static const int dim = TRefElem::dim;
465
466 if (rRefElem.roid() != ROID_PYRAMID || side != 0)
467 {
468 // number of corners of side
469 const int coOfSide = rRefElem.num(dim-1, side, 0);
470
471 // set mid ids
472 if(dim == 2)
473 {
474 vMidID[co%2] = MidID(0, rRefElem.id(1, side, 0, co)); // corner of side
475 vMidID[(co+1)%2] = MidID(1, side); // side midpoint
476 }
477 else if (dim == 3)
478 {
479 vMidID[0] = MidID(0, rRefElem.id(2, side, 0, co)); // corner of side
480 vMidID[1] = MidID(1, rRefElem.id(2, side, 1, co)); // edge co
481 vMidID[2] = MidID(2, side); // side midpoint
482 vMidID[3] = MidID(1, rRefElem.id(2, side, 1, (co -1 + coOfSide)%coOfSide)); // edge co-1
483 }
484 }
485 // bottom side of pyramid here
486 else
487 {
488 switch (co)
489 {
490 // bf of corner 0 in subface 0/0
491 case 0: vMidID[0] = MidID(0,0); // corner 0
492 vMidID[1] = MidID(1,8); // edge 2->0
493 vMidID[2] = MidID(2,5); // subface 0/0
494 vMidID[3] = MidID(1,0); // edge 0
495 break;
496 // bf of corner 1
497 case 1: vMidID[0] = MidID(0,1); // corner 1
498 vMidID[1] = MidID(1,0); // edge 0
499 vMidID[2] = MidID(2,5); // subface 0/0
500 vMidID[3] = MidID(1,1); // edge 1
501 break;
502 // bf of corner 2 in subvolume 0/0
503 case 2: vMidID[0] = MidID(0,2); // corner 2
504 vMidID[1] = MidID(1,1); // edge 1
505 vMidID[2] = MidID(2,5); // subface 0/0
506 vMidID[3] = MidID(1,8); // edge 2->0
507 break;
508 // bf of corner 0 in subvolume 1/0
509 case 3: vMidID[0] = MidID(0,0); // corner 0
510 vMidID[1] = MidID(1,3); // edge 3
511 vMidID[2] = MidID(2,6); // subface 1/0
512 vMidID[3] = MidID(1,9); // edge 0->2
513 break;
514 // bf of corner 2 in subvolume 1/0
515 case 4: vMidID[0] = MidID(0,2); // corner 2
516 vMidID[1] = MidID(1,9); // edge 0->2
517 vMidID[2] = MidID(2,6); // subface 1/0
518 vMidID[3] = MidID(1,2); // edge 2
519 break;
520 // bf of corner 3
521 case 5: vMidID[0] = MidID(0,3); // corner 3
522 vMidID[1] = MidID(1,2); // edge 2
523 vMidID[2] = MidID(2,6); // subface 1/0
524 vMidID[3] = MidID(1,3); // edge 3
525 break;
526 default:UG_THROW("Pyramid only has 6 BFs on bottom side (no. 0-5), but requested no. " << co << ".");
527 break;
528 }
529 }
530}
531
532template <int dim, int maxMid>
533static void CopyCornerByMidID(MathVector<dim> vCorner[],
534 const MidID vMidID[],
535 MathVector<dim> vvMidPos[][maxMid],
536 const size_t numCo)
537{
538 for(size_t i = 0; i < numCo; ++i)
539 {
540 const size_t d = vMidID[i].dim;
541 const size_t id = vMidID[i].id;
542 vCorner[i] = vvMidPos[d][id];
543 }
544}
545
547// FV1 Geometry for Reference Element Type
549
550template <typename TElem, int TWorldDim>
553 : m_pElem(NULL), m_rRefElem(Provider<ref_elem_type>::get()),
554 m_rTrialSpace(Provider<local_shape_fct_set_type>::get())
555{
557}
558
559template <typename TElem, int TWorldDim>
562{
563// set corners of element as local centers of nodes
564 for(size_t i = 0; i < m_rRefElem.num(0); ++i)
565 m_vvLocMid[0][i] = m_rRefElem.corner(i);
566
567// compute local midpoints
568 ComputeMidPoints<dim, ref_elem_type, maxMid>(m_rRefElem, m_vvLocMid[0], m_vvLocMid);
569
570// set up local information for SubControlVolumeFaces (scvf)
571 for(size_t i = 0; i < num_scvf(); ++i)
572 {
573
574 // this scvf separates the given nodes
575 if (m_rRefElem.REFERENCE_OBJECT_ID != ROID_PYRAMID)
576 {
577 m_vSCVF[i].From = m_rRefElem.id(1, i, 0, 0);
578 m_vSCVF[i].To = m_rRefElem.id(1, i, 0, 1);
579 }
580 // special case pyramid (scvf not mappable by edges)
581 else
582 {
583 // map according to order defined in ComputeSCVFMidID
584 m_vSCVF[i].From = ((i>6 && i%3) ? (i%3)+1 : i%3);
585 m_vSCVF[i].To = i%6 > 2 ? 4 : ((i+1)%3 + (i>5 && i<8 ? 1 : 0));
586 }
587
588 // compute mid ids of the scvf
589 ComputeSCVFMidID(m_rRefElem, m_vSCVF[i].vMidID, i);
590
591 // copy local corners of scvf
592 CopyCornerByMidID<dim, maxMid>(m_vSCVF[i].vLocPos, m_vSCVF[i].vMidID, m_vvLocMid, SCVF::numCo);
593
594 // integration point
595 AveragePositions(m_vSCVF[i].localIP, m_vSCVF[i].vLocPos, SCVF::numCo);
596 }
597
598// set up local informations for SubControlVolumes (scv)
599// each scv is associated to one corner of the element
600 for(size_t i = 0; i < num_scv(); ++i)
601 {
602 // store associated node
603 if (m_rRefElem.REFERENCE_OBJECT_ID != ROID_PYRAMID)
604 {
605 m_vSCV[i].nodeId = i;
606 }
607 // special case pyramid (scv not mappable by corners)
608 else
609 {
610 // map according to order defined in ComputeSCVMidID
611 m_vSCV[i].nodeId = i<3 ? i : (i<5 ? (i+1)%5 : i-3);
612 }
613
614 // compute mid ids scv
615 ComputeSCVMidID(m_rRefElem, m_vSCV[i].midId, i);
616
617 // copy local corners of scv
618 CopyCornerByMidID<dim, maxMid>(m_vSCV[i].vLocPos, m_vSCV[i].midId, m_vvLocMid, m_vSCV[i].num_corners());
619 }
620
621
622// compute Shapes and Derivatives
623 for(size_t i = 0; i < num_scvf(); ++i)
624 {
625 m_rTrialSpace.shapes(&(m_vSCVF[i].vShape[0]), m_vSCVF[i].local_ip());
626 m_rTrialSpace.grads(&(m_vSCVF[i].vLocalGrad[0]), m_vSCVF[i].local_ip());
627 }
628
629 for(size_t i = 0; i < num_scv(); ++i)
630 {
631 m_rTrialSpace.shapes(&(m_vSCV[i].vShape[0]), m_vSCV[i].local_ip());
632 m_rTrialSpace.grads(&(m_vSCV[i].vLocalGrad[0]), m_vSCV[i].local_ip());
633 }
634
635// copy ip positions in a list for Sub Control Volumes Faces (SCVF)
636 for(size_t i = 0; i < num_scvf(); ++i)
637 m_vLocSCVF_IP[i] = scvf(i).local_ip();
638}
639
641template <typename TElem, int TWorldDim>
643update(GridObject* elem, const MathVector<worldDim>* vCornerCoords, const ISubsetHandler* ish)
644{
645
646 UG_ASSERT(dynamic_cast<TElem*>(elem) != NULL, "Wrong element type.");
647 TElem* pElem = static_cast<TElem*>(elem);
648
649// if already update for this element, do nothing
650 if(m_pElem == pElem) return; else m_pElem = pElem;
651
652// remember global position of nodes
653 for(size_t i = 0; i < m_rRefElem.num(0); ++i)
654 m_vvGloMid[0][i] = vCornerCoords[i];
655
656// compute global midpoints
657 ComputeMidPoints<worldDim, ref_elem_type, maxMid>(m_rRefElem, m_vvGloMid[0], m_vvGloMid);
658
659// compute global informations for scvf
660 for(size_t i = 0; i < num_scvf(); ++i)
661 {
662 // copy local corners of scvf
663 CopyCornerByMidID<worldDim, maxMid>(m_vSCVF[i].vGloPos, m_vSCVF[i].vMidID, m_vvGloMid, SCVF::numCo);
664
665 // integration point
666 AveragePositions(m_vSCVF[i].globalIP, m_vSCVF[i].vGloPos, SCVF::numCo);
667
668 // normal on scvf
669 traits::NormalOnSCVF(m_vSCVF[i].Normal, m_vSCVF[i].vGloPos, m_vvGloMid[0]);
670 }
671
672// compute size of scv
673 for(size_t i = 0; i < num_scv(); ++i)
674 {
675 // copy global corners
676 CopyCornerByMidID<worldDim, maxMid>(m_vSCV[i].vGloPos, m_vSCV[i].midId, m_vvGloMid, m_vSCV[i].num_corners());
677
678 // compute volume of scv
679 m_vSCV[i].Vol = ElementSize<scv_type, worldDim>(m_vSCV[i].vGloPos);
680 }
681
682// Shapes and Derivatives
683 m_mapping.update(vCornerCoords);
684
685// if mapping is linear, compute jacobian only once and copy
687 {
689 m_mapping.jacobian_transposed_inverse(JtInv, m_vSCVF[0].local_ip());
690 const number detJ = m_mapping.sqrt_gram_det(m_vSCVF[0].local_ip());
691
692 for(size_t i = 0; i < num_scvf(); ++i)
693 {
694 m_vSCVF[i].JtInv = JtInv;
695 m_vSCVF[i].detj = detJ;
696 }
697
698 for(size_t i = 0; i < num_scv(); ++i)
699 {
700 m_vSCV[i].JtInv = JtInv;
701 m_vSCV[i].detj = detJ;
702 }
703 }
704// else compute jacobian for each integration point
705 else
706 {
707 for(size_t i = 0; i < num_scvf(); ++i)
708 {
709 m_mapping.jacobian_transposed_inverse(m_vSCVF[i].JtInv, m_vSCVF[i].local_ip());
710 m_vSCVF[i].detj = m_mapping.sqrt_gram_det(m_vSCVF[i].local_ip());
711 }
712 for(size_t i = 0; i < num_scv(); ++i)
713 {
714 m_mapping.jacobian_transposed_inverse(m_vSCV[i].JtInv, m_vSCV[i].local_ip());
715 m_vSCV[i].detj = m_mapping.sqrt_gram_det(m_vSCV[i].local_ip());
716 }
717 }
718
719// compute global gradients
720 for(size_t i = 0; i < num_scvf(); ++i)
721 for(size_t sh = 0 ; sh < scvf(i).num_sh(); ++sh)
722 MatVecMult(m_vSCVF[i].vGlobalGrad[sh], m_vSCVF[i].JtInv, m_vSCVF[i].vLocalGrad[sh]);
723
724 for(size_t i = 0; i < num_scv(); ++i)
725 for(size_t sh = 0 ; sh < scv(i).num_sh(); ++sh)
726 MatVecMult(m_vSCV[i].vGlobalGrad[sh], m_vSCV[i].JtInv, m_vSCV[i].vLocalGrad[sh]);
727
728// Copy ip pos in list for SCVF
729 for(size_t i = 0; i < num_scvf(); ++i)
730 m_vGlobSCVF_IP[i] = scvf(i).global_ip();
731
732// if no boundary subsets required, return
733 if(num_boundary_subsets() == 0 || ish == NULL) return;
734 else update_boundary_faces(pElem, vCornerCoords, ish);
735
736
737}
738
739template <typename TElem, int TWorldDim>
741update_boundary_faces(GridObject* elem, const MathVector<worldDim>* vCornerCoords, const ISubsetHandler* ish)
742{
743 UG_ASSERT(dynamic_cast<TElem*>(elem) != NULL, "Wrong element type.");
744 TElem* pElem = static_cast<TElem*>(elem);
745
746// get grid
747 Grid& grid = *(ish->grid());
748
749// vector of subset indices of side
750 std::vector<int> vSubsetIndex;
751
752// get subset indices for sides (i.e. edge in 2d, faces in 3d)
753 if(dim == 1) {
754 std::vector<Vertex*> vVertex;
755 CollectVertices(vVertex, grid, pElem);
756 vSubsetIndex.resize(vVertex.size());
757 for(size_t i = 0; i < vVertex.size(); ++i)
758 vSubsetIndex[i] = ish->get_subset_index(vVertex[i]);
759 }
760 if(dim == 2) {
761 std::vector<Edge*> vEdges;
762 CollectEdgesSorted(vEdges, grid, pElem);
763 vSubsetIndex.resize(vEdges.size());
764 for(size_t i = 0; i < vEdges.size(); ++i)
765 vSubsetIndex[i] = ish->get_subset_index(vEdges[i]);
766 }
767 if(dim == 3) {
768 std::vector<Face*> vFaces;
769 CollectFacesSorted(vFaces, grid, pElem);
770 vSubsetIndex.resize(vFaces.size());
771 for(size_t i = 0; i < vFaces.size(); ++i)
772 vSubsetIndex[i] = ish->get_subset_index(vFaces[i]);
773 }
774
775// loop requested subset
776 typename std::map<int, std::vector<BF> >::iterator it;
777 for (it=m_mapVectorBF.begin() ; it != m_mapVectorBF.end(); ++it)
778 {
779 // get subset index
780 const int bndIndex = (*it).first;
781
782 // get vector of BF for element
783 std::vector<BF>& vBF = (*it).second;
784
785 // clear vector
786 vBF.clear();
787
788 // current number of bf
789 size_t curr_bf = 0;
790
791 // loop sides of element
792 for(size_t side = 0; side < vSubsetIndex.size(); ++side)
793 {
794 // skip non boundary sides
795 if(vSubsetIndex[side] != bndIndex) continue;
796
797 // number of corners of side (special case bottom side pyramid)
798 const int coOfSide = (m_rRefElem.REFERENCE_OBJECT_ID != ROID_PYRAMID || side != 0)
799 ? m_rRefElem.num(dim-1, side, 0) : m_rRefElem.num(dim-1, side, 0) + 2;
800
801 // resize vector
802 vBF.resize(curr_bf + coOfSide);
803
804 // loop corners
805 for(int co = 0; co < coOfSide; ++co)
806 {
807 // get current bf
808 BF& bf = vBF[curr_bf];
809
810 // set node id == scv this bf belongs to
811 if (m_rRefElem.REFERENCE_OBJECT_ID != ROID_PYRAMID || side != 0)
812 bf.nodeId = m_rRefElem.id(dim-1, side, 0, co);
813 else
814 {
815 // map according to order defined in ComputeBFMidID
816 bf.nodeId = m_rRefElem.id(dim-1, side, 0, (co % 3) + (co>3 ? 1 : 0));
817 }
818
819 // Compute MidID for BF
820 ComputeBFMidID(m_rRefElem, side, bf.vMidID, co);
821
822 // copy corners of bf
823 CopyCornerByMidID<dim, maxMid>(bf.vLocPos, bf.vMidID, m_vvLocMid, BF::numCo);
824 CopyCornerByMidID<worldDim, maxMid>(bf.vGloPos, bf.vMidID, m_vvGloMid, BF::numCo);
825
826 // integration point
827 AveragePositions(bf.localIP, bf.vLocPos, BF::numCo);
828 AveragePositions(bf.globalIP, bf.vGloPos, BF::numCo);
829
830 // normal on scvf
831 traits::NormalOnSCVF(bf.Normal, bf.vGloPos, m_vvGloMid[0]);
832
833 // compute volume
834 bf.Vol = VecTwoNorm(bf.Normal);
835
836 m_rTrialSpace.shapes(&(bf.vShape[0]), bf.localIP);
837 m_rTrialSpace.grads(&(bf.vLocalGrad[0]), bf.localIP);
838
839 m_mapping.jacobian_transposed_inverse(bf.JtInv, bf.localIP);
840 bf.detj = m_mapping.sqrt_gram_det(bf.localIP);
841
842 for(size_t sh = 0 ; sh < bf.num_sh(); ++sh)
843 MatVecMult(bf.vGlobalGrad[sh], bf.JtInv, bf.vLocalGrad[sh]);
844
845 // increase curr_bf
846 ++curr_bf;
847 }
848 }
849 }
850}
851
852
854// Dim-dependent Finite Volume Geometry
856template <int TDim, int TWorldDim>
858adapt(GridObject* elem, const MathVector<TWorldDim>* vCornerCoords,
859 const ISubsetHandler* ish)
860{
861
862 for(size_t ip = 0; ip < num_scvf(); ++ip)
863 UG_LOG("Dim- vSCVF_PosOnEdge = " << vSCVF_PosOnEdge[ip] << "\n");
864
865// compute global informations for scvf
866 for(size_t i = 0; i < num_scvf(); ++i)
867 {
868 UG_LOG("Dim-ip's: " << m_vSCVF[i].globalIP << "\n");
869 UG_LOG("Dim-normals: " << m_vSCVF[i].Normal << "\n");
870
871 }
872
873// compute global informations for scvf
874 for(size_t i = 0; i < num_scvf(); ++i)
875 {
876 //VecScaleAdd(m_vSCVF[i].globalIP, vSCVF_PosOnEdge[i], 0.5, m_MidPoint, 0.5);
877 for(size_t d = 0; d < m_vSCVF[i].globalIP.size(); ++d)
878 m_vSCVF[i].globalIP[d] = 0.5*(vSCVF_PosOnEdge[i][d] + m_MidPoint[d]);
879
880 //traits::NormalOnSCVF(m_vSCVF[i].Normal, vSCVF_PosOnEdge[i], m_MidPoint);
881 // 'NormalOnSCVF' braucht Felder im 2. und 3. Parameter
882 //traits::NormalOnSCVF(m_vSCVF[i].Normal, vSCVF_PosOnEdge, vSCVF_PosOnEdge);
883 MathVector<worldDim> buffer_comp;
884 VecSubtract(buffer_comp, vSCVF_PosOnEdge[i], m_MidPoint);
885
886 MathVector<worldDim> buffer_copy;
887 buffer_copy[0] = -buffer_comp[1];
888 buffer_copy[1] = buffer_comp[0];
889
890 // check the orientation of the newly computet Normal:
891 // -> has to be in direction from -> to, i.e. VecDot(Normal_alt, Normal_neu) > 0
892 if ( VecDot(m_vSCVF[i].Normal, buffer_copy) > 0 )
893 {
894 UG_LOG("Dim-Positiv: normals: " << m_vSCVF[i].Normal << "\n");
895 UG_LOG("Dim-Positiv: normals: " << buffer_copy << "\n");
896
897 m_vSCVF[i].Normal[0] = buffer_copy[0];
898 m_vSCVF[i].Normal[1] = buffer_copy[1];
899 }
900 else
901 {
902 UG_LOG("Dim-Negativ: normals: " << m_vSCVF[i].Normal << "\n");
903 UG_LOG("Dim-Negativ: normals: " << buffer_copy << "\n");
904
905 m_vSCVF[i].Normal[0] = -buffer_copy[0];
906 m_vSCVF[i].Normal[1] = -buffer_copy[1];
907 }
908
909
910 UG_LOG("Dim-***ip's: " << m_vSCVF[i].globalIP << "\n");
911 UG_LOG("Dim-***normals: " << m_vSCVF[i].Normal << "\n");
912 UG_LOG("Dim-***buffer_copy: " << buffer_copy << "\n");
913
914 }
915
916 //UG_THROW("### ok done...\n");
917}
918
919template <int TDim, int TWorldDim>
922 const ISubsetHandler* ish)
923{
924 UG_THROW("adapt_integration_points: ok done...\n");
925
926}
927
928template <int TDim, int TWorldDim>
930adapt_normals(GridObject* elem, const MathVector<TWorldDim>* vCornerCoords,
931 const ISubsetHandler* ish)
932{
935
936
937 for(size_t i = 0; i < num_scvf(); ++i)
938 {
939 UG_LOG("ip's: " << m_vSCVF[i].globalIP << "\n");
940 UG_LOG("normals: " << m_vSCVF[i].Normal << "\n");
941
942 }
943
944 UG_THROW("stop hier...\n");
945
946
947}
948
949template <int TDim, int TWorldDim>
952{
953// get reference element
954 try{
955 const DimReferenceElement<dim>& rRefElem
956 = ReferenceElementProvider::get<dim>(m_roid);
957
958// set corners of element as local centers of nodes
959 for(size_t i = 0; i < rRefElem.num(0); ++i)
960 m_vvLocMid[0][i] = rRefElem.corner(i);
961
962// compute local midpoints
963 ComputeMidPoints<dim, DimReferenceElement<dim>, maxMid>
964 (rRefElem, m_vvLocMid[0], m_vvLocMid);
965
966// set number of scvf / scv of this roid
967 m_numSCV = (m_roid != ROID_PYRAMID) ? rRefElem.num(0) : 8; // number of corners
968 m_numSCVF = (m_roid != ROID_PYRAMID) ? rRefElem.num(1) : 12; // number of edges
969
970// set up local informations for SubControlVolumeFaces (scvf)
971// each scvf is associated to one edge of the element
972 for(size_t i = 0; i < num_scvf(); ++i)
973 {
974 // this scvf separates the given nodes
975 if (m_roid != ROID_PYRAMID)
976 {
977 m_vSCVF[i].From = rRefElem.id(1, i, 0, 0);
978 m_vSCVF[i].To = rRefElem.id(1, i, 0, 1);
979 }
980 // special case pyramid (scvf not mappable by edges)
981 else
982 {
983 // map according to order defined in ComputeSCVFMidID
984 m_vSCVF[i].From = ((i>6 && i%3) ? (i%3)+1 : i%3);
985 m_vSCVF[i].To = i%6 > 2 ? 4 : ((i+1)%3 + (i>5 && i<8 ? 1 : 0));
986 }
987
988
989 // compute mid ids of the scvf
990 ComputeSCVFMidID(rRefElem, m_vSCVF[i].vMidID, i);
991
992 // copy local corners of scvf
993 CopyCornerByMidID<dim, maxMid>(m_vSCVF[i].vLocPos, m_vSCVF[i].vMidID, m_vvLocMid, SCVF::numCo);
994
995 // integration point
996 AveragePositions(m_vSCVF[i].localIP, m_vSCVF[i].vLocPos, SCVF::numCo);
997 }
998
999// set up local informations for SubControlVolumes (scv)
1000// each scv is associated to one corner of the element
1001 for(size_t i = 0; i < num_scv(); ++i)
1002 {
1003 // store associated node
1004 if (m_roid != ROID_PYRAMID)
1005 {
1006 m_vSCV[i].nodeId = i;
1007 }
1008 // special case pyramid (scv not mappable by corners)
1009 else
1010 {
1011 // map according to order defined in ComputeSCVMidID
1012 m_vSCV[i].nodeId = i<3 ? i : (i<5 ? (i+1)%5 : i-3);
1013 }
1014
1015 // compute mid ids scv
1016 ComputeSCVMidID(rRefElem, m_vSCV[i].vMidID, i);
1017
1018 // copy local corners of scv
1019 CopyCornerByMidID<dim, maxMid>(m_vSCV[i].vLocPos, m_vSCV[i].vMidID, m_vvLocMid, m_vSCV[i].num_corners());
1020 }
1021
1023 // Shapes and Derivatives
1025
1026 const LocalShapeFunctionSet<dim>& TrialSpace =
1027 LocalFiniteElementProvider::get<dim>(m_roid, LFEID(LFEID::LAGRANGE, dim, 1));
1028
1029 m_nsh = TrialSpace.num_sh();
1030
1031 for(size_t i = 0; i < num_scvf(); ++i)
1032 {
1033 m_vSCVF[i].numSH = TrialSpace.num_sh();
1034 TrialSpace.shapes(&(m_vSCVF[i].vShape[0]), m_vSCVF[i].localIP);
1035 TrialSpace.grads(&(m_vSCVF[i].vLocalGrad[0]), m_vSCVF[i].localIP);
1036 }
1037
1038 for(size_t i = 0; i < num_scv(); ++i)
1039 {
1040 m_vSCV[i].numSH = TrialSpace.num_sh();
1041 TrialSpace.shapes(&(m_vSCV[i].vShape[0]), m_vSCV[i].vLocPos[0]);
1042 TrialSpace.grads(&(m_vSCV[i].vLocalGrad[0]), m_vSCV[i].vLocPos[0]);
1043 }
1044
1045 }
1046 UG_CATCH_THROW("DimFV1IBGeometry: update failed.");
1047
1048// copy ip positions in a list for Sub Control Volumes Faces (SCVF)
1049 for(size_t i = 0; i < num_scvf(); ++i)
1050 m_vLocSCVF_IP[i] = scvf(i).local_ip();
1051}
1052
1053
1055template <int TDim, int TWorldDim>
1057update(GridObject* pElem, const MathVector<worldDim>* vCornerCoords, const ISubsetHandler* ish)
1058{
1059// If already update for this element, do nothing
1060 if(m_pElem == pElem) return; else m_pElem = pElem;
1061
1062// refresh local data, if different roid given
1063 if(m_roid != pElem->reference_object_id())
1064 {
1065 // remember new roid
1066 m_roid = (ReferenceObjectID) pElem->reference_object_id();
1067
1068 // update local data
1069 update_local_data();
1070 }
1071
1072// get reference element
1073 try{
1074 const DimReferenceElement<dim>& rRefElem
1075 = ReferenceElementProvider::get<dim>(m_roid);
1076
1077// remember global position of nodes
1078 for(size_t i = 0; i < rRefElem.num(0); ++i)
1079 m_vvGloMid[0][i] = vCornerCoords[i];
1080
1081// compute local midpoints
1082 ComputeMidPoints<worldDim, DimReferenceElement<dim>, maxMid>(rRefElem, m_vvGloMid[0], m_vvGloMid);
1083
1084// compute global informations for scvf
1085 for(size_t i = 0; i < num_scvf(); ++i)
1086 {
1087 // copy local corners of scvf
1088 CopyCornerByMidID<worldDim, maxMid>(m_vSCVF[i].vGloPos, m_vSCVF[i].vMidID, m_vvGloMid, SCVF::numCo);
1089
1090 // integration point
1091 AveragePositions(m_vSCVF[i].globalIP, m_vSCVF[i].vGloPos, SCVF::numCo);
1092
1093 // normal on scvf
1094 traits::NormalOnSCVF(m_vSCVF[i].Normal, m_vSCVF[i].vGloPos, m_vvGloMid[0]);
1095 }
1096
1097// compute size of scv
1098 for(size_t i = 0; i < num_scv(); ++i)
1099 {
1100 // copy global corners
1101 CopyCornerByMidID<worldDim, maxMid>(m_vSCV[i].vGloPos, m_vSCV[i].vMidID, m_vvGloMid, m_vSCV[i].num_corners());
1102
1103 // compute volume of scv
1104 m_vSCV[i].Vol = ElementSize<scv_type, worldDim>(m_vSCV[i].vGloPos);
1105 }
1106
1107// get reference mapping
1108 DimReferenceMapping<dim, worldDim>& rMapping = ReferenceMappingProvider::get<dim, worldDim>(m_roid);
1109 rMapping.update(vCornerCoords);
1110
1111 //\todo compute with on virt. call
1112// compute jacobian for linear mapping
1113 if(rMapping.is_linear())
1114 {
1116 rMapping.jacobian_transposed_inverse(JtInv, m_vSCVF[0].local_ip());
1117 const number detJ = rMapping.sqrt_gram_det(m_vSCVF[0].local_ip());
1118
1119 for(size_t i = 0; i < num_scvf(); ++i)
1120 {
1121 m_vSCVF[i].JtInv = JtInv;
1122 m_vSCVF[i].detj = detJ;
1123 }
1124
1125 for(size_t i = 0; i < num_scv(); ++i)
1126 {
1127 m_vSCV[i].JtInv = JtInv;
1128 m_vSCV[i].detj = detJ;
1129 }
1130 }
1131// else compute jacobian for each integration point
1132 else
1133 {
1134 for(size_t i = 0; i < num_scvf(); ++i)
1135 {
1136 rMapping.jacobian_transposed_inverse(m_vSCVF[i].JtInv, m_vSCVF[i].local_ip());
1137 m_vSCVF[i].detj = rMapping.sqrt_gram_det(m_vSCVF[i].local_ip());
1138 }
1139 for(size_t i = 0; i < num_scv(); ++i)
1140 {
1141 rMapping.jacobian_transposed_inverse(m_vSCV[i].JtInv, m_vSCV[i].local_ip());
1142 m_vSCV[i].detj = rMapping.sqrt_gram_det(m_vSCV[i].local_ip());
1143 }
1144 }
1145
1146// compute global gradients
1147 for(size_t i = 0; i < num_scvf(); ++i)
1148 for(size_t sh = 0; sh < scvf(i).num_sh(); ++sh)
1149 MatVecMult(m_vSCVF[i].vGlobalGrad[sh], m_vSCVF[i].JtInv, m_vSCVF[i].vLocalGrad[sh]);
1150
1151 for(size_t i = 0; i < num_scv(); ++i)
1152 for(size_t sh = 0; sh < scv(i).num_sh(); ++sh)
1153 MatVecMult(m_vSCV[i].vGlobalGrad[sh], m_vSCV[i].JtInv, m_vSCV[i].vLocalGrad[sh]);
1154
1155// copy ip points in list (SCVF)
1156 for(size_t i = 0; i < num_scvf(); ++i)
1157 m_vGlobSCVF_IP[i] = scvf(i).global_ip();
1158
1159 }
1160 UG_CATCH_THROW("DimFV1IBGeometry: update failed.");
1161
1162// if no boundary subsets required, return
1163 if(num_boundary_subsets() == 0 || ish == NULL) return;
1164 else update_boundary_faces(pElem, vCornerCoords, ish);
1165}
1166
1167template <int TDim, int TWorldDim>
1169update_boundary_faces(GridObject* pElem, const MathVector<worldDim>* vCornerCoords, const ISubsetHandler* ish)
1170{
1171// get grid
1172 Grid& grid = *(ish->grid());
1173
1174// vector of subset indices of side
1175 std::vector<int> vSubsetIndex;
1176
1177// get subset indices for sides (i.e. edge in 2d, faces in 3d)
1178 if(dim == 1) {
1179 std::vector<Vertex*> vVertex;
1180 CollectVertices(vVertex, grid, pElem);
1181 vSubsetIndex.resize(vVertex.size());
1182 for(size_t i = 0; i < vVertex.size(); ++i)
1183 vSubsetIndex[i] = ish->get_subset_index(vVertex[i]);
1184 }
1185 if(dim == 2) {
1186 std::vector<Edge*> vEdges;
1187 CollectEdgesSorted(vEdges, grid, pElem);
1188 vSubsetIndex.resize(vEdges.size());
1189 for(size_t i = 0; i < vEdges.size(); ++i)
1190 vSubsetIndex[i] = ish->get_subset_index(vEdges[i]);
1191 }
1192 if(dim == 3) {
1193 std::vector<Face*> vFaces;
1194 CollectFacesSorted(vFaces, grid, pElem);
1195 vSubsetIndex.resize(vFaces.size());
1196 for(size_t i = 0; i < vFaces.size(); ++i)
1197 vSubsetIndex[i] = ish->get_subset_index(vFaces[i]);
1198 }
1199
1200 try{
1201 const DimReferenceElement<dim>& rRefElem
1202 = ReferenceElementProvider::get<dim>(m_roid);
1203
1204 DimReferenceMapping<dim, worldDim>& rMapping = ReferenceMappingProvider::get<dim, worldDim>(m_roid);
1205 rMapping.update(vCornerCoords);
1206
1207 const LocalShapeFunctionSet<dim>& TrialSpace =
1208 LocalFiniteElementProvider::get<dim>(m_roid, LFEID(LFEID::LAGRANGE, dim, 1));
1209
1210// loop requested subset
1211 typename std::map<int, std::vector<BF> >::iterator it;
1212 for (it=m_mapVectorBF.begin() ; it != m_mapVectorBF.end(); ++it)
1213 {
1214 // get subset index
1215 const int bndIndex = (*it).first;
1216
1217 // get vector of BF for element
1218 std::vector<BF>& vBF = (*it).second;
1219
1220 // clear vector
1221 vBF.clear();
1222
1223 // current number of bf
1224 size_t curr_bf = 0;
1225
1226 // loop sides of element
1227 for(size_t side = 0; side < vSubsetIndex.size(); ++side)
1228 {
1229 // skip non boundary sides
1230 if(vSubsetIndex[side] != bndIndex) continue;
1231
1232 // number of corners of side (special case bottom side pyramid)
1233 const int coOfSide = (pElem->reference_object_id() != ROID_PYRAMID || side != 0)
1234 ? rRefElem.num(dim-1, side, 0) : rRefElem.num(dim-1, side, 0) + 2;
1235 // resize vector
1236 vBF.resize(curr_bf + coOfSide);
1237
1238 // loop corners
1239 for(int co = 0; co < coOfSide; ++co)
1240 {
1241 // get current bf
1242 BF& bf = vBF[curr_bf];
1243
1244 // set node id == scv this bf belongs to
1245 if (pElem->reference_object_id() != ROID_PYRAMID || side != 0)
1246 bf.nodeId = rRefElem.id(dim-1, side, 0, co);
1247 else
1248 {
1249 // map according to order defined in ComputeBFMidID
1250 bf.nodeId = rRefElem.id(dim-1, side, 0, (co % 3) + (co>3 ? 1 : 0));
1251 }
1252
1253 // Compute MidID for BF
1254 ComputeBFMidID(rRefElem, side, bf.vMidID, co);
1255
1256 // copy corners of bf
1257 CopyCornerByMidID<dim, maxMid>(bf.vLocPos, bf.vMidID, m_vvLocMid, BF::numCo);
1258 CopyCornerByMidID<worldDim, maxMid>(bf.vGloPos, bf.vMidID, m_vvGloMid, BF::numCo);
1259
1260 // integration point
1261 AveragePositions(bf.localIP, bf.vLocPos, BF::numCo);
1262 AveragePositions(bf.globalIP, bf.vGloPos, BF::numCo);
1263
1264 // normal on scvf
1265 traits::NormalOnSCVF(bf.Normal, bf.vGloPos, m_vvGloMid[0]);
1266
1267 // compute volume
1268 bf.Vol = VecTwoNorm(bf.Normal);
1269
1270 // compute shapes and grads
1271 bf.numSH = TrialSpace.num_sh();
1272 TrialSpace.shapes(&(bf.vShape[0]), bf.localIP);
1273 TrialSpace.grads(&(bf.vLocalGrad[0]), bf.localIP);
1274
1275 // get reference mapping
1276 rMapping.jacobian_transposed_inverse(bf.JtInv, bf.localIP);
1277 bf.detj = rMapping.sqrt_gram_det(bf.localIP);
1278
1279 // compute global gradients
1280 for(size_t sh = 0 ; sh < bf.num_sh(); ++sh)
1281 MatVecMult(bf.vGlobalGrad[sh], bf.JtInv, bf.vLocalGrad[sh]);
1282
1283 // increase curr_bf
1284 ++curr_bf;
1285 }
1286 }
1287 }
1288
1289 }
1290 UG_CATCH_THROW("DimFV1IBGeometry: update failed.");
1291}
1292
1294// FV1IBGeometry
1295
1296template class FV1IBGeometry<RegularEdge, 1>;
1297template class FV1IBGeometry<RegularEdge, 2>;
1298template class FV1IBGeometry<RegularEdge, 3>;
1299
1300template class FV1IBGeometry<Triangle, 2>;
1301template class FV1IBGeometry<Triangle, 3>;
1302
1303template class FV1IBGeometry<Quadrilateral, 2>;
1304template class FV1IBGeometry<Quadrilateral, 3>;
1305
1306template class FV1IBGeometry<Tetrahedron, 3>;
1307template class FV1IBGeometry<Prism, 3>;
1308template class FV1IBGeometry<Pyramid, 3>;
1309template class FV1IBGeometry<Hexahedron, 3>;
1310
1311
1312
1314// DimFV1IBGeometry
1315template class DimFV1IBGeometry<1, 1>;
1316template class DimFV1IBGeometry<1, 2>;
1317template class DimFV1IBGeometry<1, 3>;
1318
1319template class DimFV1IBGeometry<2, 2>;
1320template class DimFV1IBGeometry<2, 3>;
1321
1322template class DimFV1IBGeometry<3, 3>;
1323
1324
void VecScaleAdd(vector_t &dest, double alpha1, const vector_t &v1, double alpha2, const vector_t &v2)
calculates dest = alpha1*v1 + alpha2*v2
boundary face
Definition fv1ib_geom.h:836
MathVector< dim > vLocPos[numCo]
Definition fv1ib_geom.h:914
size_t nodeId
Definition fv1ib_geom.h:908
number Vol
Definition fv1ib_geom.h:922
number detj
Definition fv1ib_geom.h:930
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
size_t num_sh() const
number of shape functions
Definition fv1ib_geom.h:869
MathVector< worldDim > globalIP
Definition fv1ib_geom.h:920
MathVector< worldDim > Normal
Definition fv1ib_geom.h:921
number vShape[maxNSH]
Definition fv1ib_geom.h:926
MidID vMidID[numCo]
Definition fv1ib_geom.h:916
MathVector< dim > localIP
Definition fv1ib_geom.h:919
MathVector< worldDim > vGloPos[numCo]
Definition fv1ib_geom.h:915
MathMatrix< worldDim, dim > JtInv
Definition fv1ib_geom.h:929
Geometry and shape functions for 1st order Vertex-Centered Finite Volume.
Definition fv1ib_geom.h:598
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
void update(GridObject *elem, const MathVector< worldDim > *vCornerCoords, const ISubsetHandler *ish=NULL)
update data for given element
Definition fv1ib_geom_impl.h:1057
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
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
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
dimension dependent base class for reference elements
Definition reference_element.h:183
const MathVector< dim > & corner(size_t i) const
coordinates of reference corner (i = 0 ... num(0))
Definition reference_element.h:192
virtual base class for reference mappings
Definition reference_mapping_provider.h:53
virtual number sqrt_gram_det(const MathVector< dim > &locPos) const =0
returns the determinate of the jacobian
virtual number jacobian_transposed_inverse(MathMatrix< worldDim, dim > &JTInv, const MathVector< dim > &locPos) const =0
returns transposed of the inverse of the jacobian and returns sqrt of gram determinante
virtual void update(const MathVector< worldDim > *vCornerCoord)=0
refresh mapping for new set of corners
virtual bool is_linear() const =0
returns if mapping is affine
boundary face
Definition fv1ib_geom.h:324
MathVector< worldDim > vGloPos[numCo]
Definition fv1ib_geom.h:403
size_t num_sh() const
number of shape functions
Definition fv1ib_geom.h:357
MathMatrix< worldDim, dim > JtInv
Definition fv1ib_geom.h:416
MathVector< worldDim > Normal
Definition fv1ib_geom.h:409
MathVector< dim > vLocalGrad[nsh]
Definition fv1ib_geom.h:414
MidID vMidID[numCo]
Definition fv1ib_geom.h:404
number detj
Definition fv1ib_geom.h:417
MathVector< dim > localIP
Definition fv1ib_geom.h:407
size_t nodeId
Definition fv1ib_geom.h:396
MathVector< worldDim > globalIP
Definition fv1ib_geom.h:408
number Vol
Definition fv1ib_geom.h:410
number vShape[nsh]
Definition fv1ib_geom.h:413
MathVector< dim > vLocPos[numCo]
Definition fv1ib_geom.h:402
MathVector< worldDim > vGlobalGrad[nsh]
Definition fv1ib_geom.h:415
Geometry and shape functions for 1st order Vertex-Centered Finite Volume.
Definition fv1ib_geom.h:69
LagrangeP1< ref_elem_type > local_shape_fct_set_type
type of Shape function used
Definition fv1ib_geom.h:109
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
void update(GridObject *elem, const MathVector< worldDim > *vCornerCoords, const ISubsetHandler *ish=NULL)
update data for given element
Definition fv1ib_geom_impl.h:643
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
Manages the elements of a grid and their interconnection.
Definition grid.h:132
The base class for all geometric objects, such as vertices, edges, faces, volumes,...
Definition grid_base_objects.h:157
virtual ReferenceObjectID reference_object_id() const =0
Definition subset_handler_interface.h:223
Grid * grid() const
returns a pointer to the grid on which the subset-handler works.
Definition subset_handler_interface.cpp:304
int get_subset_index(GridObject *elem) const
Definition subset_handler_interface.cpp:560
Identifier for Local Finite Elements.
Definition local_finite_element_id.h:98
@ LAGRANGE
Definition local_finite_element_id.h:104
virtual size_t num_sh() const
Definition local_dof_set.cpp:46
virtual base class for local shape function sets
Definition local_shape_function_set.h:70
virtual void grads(grad_type *vGrad, const MathVector< dim > &x) const =0
returns all gradients evaluated at a point
virtual void shapes(shape_type *vShape, const MathVector< dim > &x) const =0
returns all shape functions evaluated at a point
A class for fixed size, dense matrices.
Definition math_matrix.h:63
a mathematical Vector with N entries.
Definition math_vector.h:97
Provider, holding a single instance of an object.
Definition provider.h:48
size_t num(int dim) const
returns the number of geometric objects of dim
Definition reference_element.h:95
int id(int dim_i, size_t i, int dim_j, size_t j) const
id of object j in dimension dim_j of obj i in dimension dim_i
Definition reference_element.h:127
Definition reference_mapping.h:65
static void ComputeSCVMidID(const TRefElem &rRefElem, MidID vMidID[], int i)
Definition fv1ib_geom_impl.h:334
static void ComputeBFMidID(const TRefElem &rRefElem, int side, MidID vMidID[], int co)
Definition fv1ib_geom_impl.h:461
static void ComputeSCVFMidID(const TRefElem &rRefElem, MidID vMidID[], int i)
Definition fv1ib_geom_impl.h:222
static void ComputeMidPoints(const TRefElem &rRefElem, const MathVector< dim > vCorner[], MathVector< dim > vvMid[][maxMid])
Definition fv1ib_geom_impl.h:147
static void CopyCornerByMidID(MathVector< dim > vCorner[], const MidID vMidID[], MathVector< dim > vvMidPos[][maxMid], const size_t numCo)
Definition fv1ib_geom_impl.h:533
void CollectEdgesSorted(vector< Edge * > &vEdgesOut, Grid &grid, Vertex *v, bool clearContainer)
Collects all edges of a vertex, thus, none.
Definition grid_util.cpp:205
void CollectFacesSorted(vector< Face * > &vFacesOut, Grid &grid, Vertex *v, bool clearContainer)
Collects all Faces of a Vertex, thus, none.
Definition grid_util.cpp:401
void CollectVertices(std::vector< Vertex * > &vVertexOut, Grid &grid, Vertex *v, bool clearContainer)
Collects all vertices.
Definition grid_util.cpp:141
#define UG_ASSERT(expr, msg)
Definition assert.h:70
#define UG_CATCH_THROW(msg)
Definition error.h:64
#define UG_THROW(msg)
Definition error.h:57
#define UG_LOG(msg)
Definition log.h:367
double number
Definition types.h:124
void MatVecMult(vector_t_out &vOut, const matrix_t &m, const vector_t_in &v)
Matrix - Vector Multiplication.
Definition math_matrix_vector_functions_common_impl.hpp:49
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
vector_t::value_type VecTwoNorm(const vector_t &v)
Definition math_vector_functions_common_impl.hpp:577
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
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
static void ComputeBFMidID(const TRefElem &rRefElem, int side, MidID vMidID[], int co)
Definition fv1_geom.cpp:204
static void ComputeSCVMidID(const TRefElem &rRefElem, MidID vMidID[], int i)
Definition fv1_geom.cpp:156
void AveragePositions(TPosition &vOut, const TPosition *vCornerCoords, size_t num)
averages positions by arithmetic mean
Definition fv_util.h:59
static void ComputeSCVFMidID(const TRefElem &rRefElem, MidID vMidID[], int i)
Definition fv1_geom.cpp:118