3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
fvmpfal3dinteractionvolumecontaineradaptive.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_FVMPFAL3D_2P2CINTERACTIONVOLUMECONTAINER_ADAPTIVE_HH
25#define DUMUX_FVMPFAL3D_2P2CINTERACTIONVOLUMECONTAINER_ADAPTIVE_HH
26
27// dumux environment
29
30namespace Dumux {
31
44template<class TypeTag>
46{
49
50 enum
51 {
52 dim = GridView::dimension, dimWorld = GridView::dimensionworld
53 };
54
56
58
61
63
64 using Element = typename GridView::Traits::template Codim<0>::Entity;
65 using ElementGeometry = typename Element::Geometry;
66 using Vertex = typename GridView::Traits::template Codim<dim>::Entity;
67
68 using IntersectionIterator = typename GridView::IntersectionIterator;
69 using Intersection = typename GridView::Intersection;
70 using IntersectionGeometry = typename Intersection::Geometry;
71
72 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
73 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
74
75 using DimVector = Dune::FieldVector<Scalar, dim>;
76
77 enum
78 {
79 pressureEqIdx = Indices::pressureEqIdx,
80 };
81 enum
82 {
83 innerEdgeFace = 2, innerSideFace = 1
84 };
85 enum
86 {
87 realFaceArea = 0, fluxFaceArea = 1
88 };
89
90public:
93
94 using GlobalInteractionVolumeVector = std::vector<InteractionVolume>;
95 using FaceAreaVector = std::vector<Dune::FieldVector<Dune::FieldVector<Scalar, 2>, 2*dim> >;
96
98
99 inline int getMpfaCase8cells(const IntersectionIterator& isIt,
100 const int localidxLarge,
102 bool& properFluxDirection);
103 inline int getMpfaCase6cells(const IntersectionIterator& isIt,
105 bool& properFluxDirection);
106 inline int getMpfaCase2or4cells(const IntersectionIterator& isIt,
108 bool& properFluxDirection);
109
111 : ParentType(problem), problem_(problem)
112 {
113 if (dim != 3)
114 {
115 DUNE_THROW(Dune::NotImplemented, "Dimension not supported!");
116 }
117 }
118private:
119 Problem& problem_;
120};
121
127template<class TypeTag>
129{
130 return;
131}
132
149template<class TypeTag>
151 const int localIdxLarge,
152 InteractionVolume& interactionVolume,
153 bool& properFluxDirection)
154{
155 int subVolumeFaceIdx = -1;
156
157 /************** determine subVolumeFace of interest */
158 int localFaceIdxInside = isIt->indexInInside();
159 switch(localIdxLarge)
160 {
161 case 0:{ //cell I = 7
162 if(localFaceIdxInside == 0) //left
163 {
164 subVolumeFaceIdx = 6;
165 properFluxDirection = false;
166 }
167 else if(localFaceIdxInside == 2) //front
168 subVolumeFaceIdx = 5;
169 else if(localFaceIdxInside == 4) //bottom
170 subVolumeFaceIdx = 10;
171 break;
172 }
173 case 1:{ //cell I = 6
174 if(localFaceIdxInside == 1) //right
175 subVolumeFaceIdx = 6;
176 else if(localFaceIdxInside == 2) //front
177 {
178 subVolumeFaceIdx = 7;
179 properFluxDirection = false;
180 }
181 else if(localFaceIdxInside == 4) //bottom
182 {
183 subVolumeFaceIdx = 11;
184 properFluxDirection = false;
185 }
186 break;
187 }
188 case 2:{ // cellI = 5
189 if(localFaceIdxInside == 0) // left
190 subVolumeFaceIdx = 4;
191 else if(localFaceIdxInside == 3) //rear
192 {
193 subVolumeFaceIdx = 5;
194 properFluxDirection = false;
195 }
196 else if(localFaceIdxInside == 4) //bottom
197 {
198 subVolumeFaceIdx = 9;
199 properFluxDirection = false;
200 }
201 break;
202 }
203 case 3:{ //cellI = 4
204 if(localFaceIdxInside == 1) //right
205 {
206 subVolumeFaceIdx = 4;
207 properFluxDirection = false;
208 }
209 else if(localFaceIdxInside == 3) //rear
210 subVolumeFaceIdx = 7;
211 else if(localFaceIdxInside == 4) //botm
212 subVolumeFaceIdx = 8;
213 break;
214 }
215 case 4:{ //cellI = 3
216 if(localFaceIdxInside == 0) //left
217 subVolumeFaceIdx = 2;
218 else if(localFaceIdxInside == 2) //front
219 {
220 subVolumeFaceIdx = 1;
221 properFluxDirection = false;
222 }
223 else if(localFaceIdxInside == 5) //top
224 {
225 subVolumeFaceIdx = 10;
226 properFluxDirection = false;
227 }
228 break;
229 }
230 case 5:{ //cellI = 2
231 if(localFaceIdxInside == 1) //right
232 {
233 subVolumeFaceIdx = 2;
234 properFluxDirection = false;
235 }
236 else if(localFaceIdxInside == 2) //front
237 subVolumeFaceIdx = 3;
238 else if(localFaceIdxInside == 5) //top
239 subVolumeFaceIdx = 11;
240 break;
241 }
242 case 6:{ //cellI = 1
243 if(localFaceIdxInside == 0) // left
244 {
245 subVolumeFaceIdx = 0;
246 properFluxDirection = false;
247 }
248 else if(localFaceIdxInside == 3) //rear
249 subVolumeFaceIdx = 1;
250 else if(localFaceIdxInside == 5) //top
251 subVolumeFaceIdx = 9;
252 break;
253 }
254 case 7:{ //cellI= 0
255 if(localFaceIdxInside == 1) //right
256 subVolumeFaceIdx = 0;
257 else if(localFaceIdxInside == 3)//rear
258 {
259 subVolumeFaceIdx = 3;
260 properFluxDirection = false;
261 }
262 else if(localFaceIdxInside == 5)//top
263 {
264 subVolumeFaceIdx = 8;
265 properFluxDirection = false;
266 }
267 break;
268 }
269 }
270 return subVolumeFaceIdx;
271}
272
285template<class TypeTag>
287 InteractionVolume& interactionVolume,
288 bool& properFluxDirection)
289{
290 int mapI(-5), mapJ(-5);
291 int mapI2(-5), mapJ2(-5);
292 // search for cell I
293 for (int i= 0; i<interactionVolume.getElementNumber(); i++)
294 {
295 if(isIt->inside() == interactionVolume.getSubVolumeElement(i))
296 {
297 if(mapI == -5)
298 mapI = i;
299 else
300 mapI2 = i;
301 }
302
303 if(isIt->outside() == interactionVolume.getSubVolumeElement(i))
304 {
305 if(mapJ == -5)
306 mapJ = i;
307 else
308 mapJ2 =i;
309 }
310 }
311
312 for(int passionfruit = 0; passionfruit <=1; passionfruit++) // loop at most twice
313 {
314 if(mapI== 0 || mapJ == 0)
315 {
316 if(mapI== 2)
317 return 3;
318 else if (mapJ== 2){
319 properFluxDirection = false;
320 return 3;
321 }
322 else if(mapJ== 1)
323 return 0;
324 else if (mapI==1){
325 properFluxDirection = false;
326 return 0;
327 }
328 }
329 // no "else if" because there is also 4 - 0 treated under I or J == 4
330 if(mapI== 1 || mapJ == 1)
331 {
332 if (mapI== 3){
333 properFluxDirection = false;
334 return 1;
335 }
336 else if(mapJ== 3)
337 return 1;
338 else if(mapI == 5){
339 properFluxDirection = false;
340 return 9;
341 }
342 else if(mapJ == 5)
343 return 9;
344 }
345 else if(mapI== 2 || mapJ == 2)
346 {
347 if(mapI == 6){
348 properFluxDirection = false;
349 return 11;
350 }
351 else if(mapJ ==6)
352 return 11;
353 }
354 else if(mapI== 3 || mapJ == 3)
355 {
356 if (mapJ== 7){
357 properFluxDirection = false;
358 return 10;
359 }
360 else if(mapI== 7)
361 return 10;
362 }
363 else if(mapI== 4 || mapJ == 4)
364 {
365 // this has to be subVolFaceIdx 8 because 7 would mean
366 // a case that should be modelled by tpfa
367 if (mapI== 0){
368 properFluxDirection = false;
369 return 8;
370 }
371 else if(mapJ== 0)
372 return 8;
373 else if(mapI == 6){
374 properFluxDirection = false;
375 return 7;
376 }
377 else if(mapJ == 6)
378 return 7;
379 }
380 else if(mapI== 5 || mapJ == 5)
381 {
382 if (mapJ== 7){
383 properFluxDirection = false;
384 return 5;
385 }
386 else if(mapI== 7)
387 return 5;
388 }
389 else if(mapI== 6 || mapJ == 6)
390 {
391 if (mapI== 7){
392 properFluxDirection = false;
393 return 6;
394 }
395 else if(mapJ== 7)
396 return 6;
397 }
398 // configuration not found: I or J are twice in the interaction volume: investigate other localIdx
399 if(mapI2 != -5)
400 mapI = mapI2;
401 else if(mapJ2 != -5)
402 mapJ = mapJ2;
403 }
404
405 Dune::dgrave << " Could not find "<< interactionVolume.getHangingNodeType() <<" case configuration for I = "
406 << problem_.variables().index(isIt->inside()) << " localIdx " << mapI << " , "
407 << problem_.variables().index(isIt->outside()) << " localIdx " << mapJ << std::endl;
408
409 return -1;
410}
411
424template<class TypeTag>
426 InteractionVolume& interactionVolume,
427 bool& properFluxDirection)
428{
429 int mapI(-5), mapJ(-5);
430 int mapI2(-5), mapJ2(-5);
431 // search for cell I
432 for (int i= 0; i<interactionVolume.getElementNumber(); i++)
433 {
434 if(isIt->inside() == interactionVolume.getSubVolumeElement(i))
435 {
436 if(mapI == -5)
437 mapI = i;
438 else
439 mapI2 = i;
440 }
441
442 if(isIt->outside() == interactionVolume.getSubVolumeElement(i))
443 {
444 if(mapJ == -5)
445 mapJ = i;
446 else
447 mapJ2 =i;
448 }
449 }
450
451 for(int passionfruit = 0; passionfruit <=1; passionfruit++) // loop at most twice
452 {
453
454 if(mapI== 0 || mapJ == 0)
455 {
456 if (mapJ== 2){
457 properFluxDirection = false;
458 return 3;
459 }
460 else if(mapI== 2)
461 return 3;
462 else if (mapJ== 4){
463 properFluxDirection = false;
464 return 8;
465 }
466 else if(mapI== 4)
467 return 8;
468 //it has to be case 0
469 else if (mapJ == 0)
470 properFluxDirection = false;
471
472 return 0;
473 }
474 else if(mapI== 1 || mapJ == 1)
475 {
476 if (mapI== 3){
477 properFluxDirection = false;
478 return 1;
479 }
480 else if(mapJ== 3)
481 return 1;
482 else if(mapI == 5){
483 properFluxDirection = false;
484 return 9;
485 }
486 else if(mapJ == 5)
487 return 9;
488 }
489 else if(mapI== 2 || mapJ == 2)
490 {
491 if (mapJ== 3){
492 properFluxDirection = false;
493 return 2;
494 }
495 else if(mapI== 3)
496 return 2;
497 else if(mapI == 6){
498 properFluxDirection = false;
499 return 11;
500 }
501 else if(mapJ ==6)
502 return 11;
503 }
504 else if(mapI== 3 || mapJ == 3)
505 {
506 if (mapJ== 7){
507 properFluxDirection = false;
508 return 10;
509 }
510 else if(mapI== 7)
511 return 10;
512 }
513 else if(mapI== 4 || mapJ == 4)
514 {
515 if (mapJ== 5){
516 properFluxDirection = false;
517 return 4;
518 }
519 else if(mapI== 5)
520 return 4;
521 else if(mapI == 6){
522 properFluxDirection = false;
523 return 7;
524 }
525 else if(mapJ == 6)
526 return 7;
527 }
528 else if(mapI== 5 || mapJ == 5)
529 {
530 if (mapJ== 7){
531 properFluxDirection = false;
532 return 5;
533 }
534 else if(mapI== 7)
535 return 5;
536 }
537 else if(mapI== 6 || mapJ == 6)
538 {
539 if (mapI== 7){
540 properFluxDirection = false;
541 return 6;
542 }
543 else if(mapJ== 7)
544 return 6;
545 }
546 // configuration not found: I or J are twice in the interaction volume: investigate other localIdx
547 if(mapI2 != -5)
548 mapI = mapI2;
549 else if(mapJ2 != -5)
550 mapJ = mapJ2;
551 }
552
553 Dune::dgrave << " Could not find "<< interactionVolume.getHangingNodeType() <<" case configuration for I = "
554 << problem_.variables().index(isIt->inside()) << " localIdx " << mapI << " and "
555 << problem_.variables().index(isIt->outside()) << " localIdx " << mapJ << std::endl;
556 return -1;
557}
558
559
560} // end namespace Dumux
561#endif
Interactionvolume container for 3-d MPFA L-method on an h-adaptive grid.
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type GetProp
get the type of a property (equivalent to old macro GET_PROP(...))
Definition: propertysystem.hh:140
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
Interactionvolume container for 3-d MPFA L-method.
Definition: 3dinteractionvolumecontainer.hh:46
GetPropType< TypeTag, Properties::MPFAInteractionVolume > InteractionVolume
Type for storing an MPFA-interaction-volume. (Usually of type FvMpfaL3dInteractionVolume or FvMpfaL3d...
Definition: 3dinteractionvolumecontainer.hh:92
InteractionVolume & interactionVolume(int vertexIdx)
Returns an interaction volume.
Definition: 3dinteractionvolumecontainer.hh:141
Interactionvolume container for 3-d MPFA L-method on an h-adaptive grid.
Definition: 3dinteractionvolumecontaineradaptive.hh:47
Interaction volume container for compositional adaptive 3-d (using MPFA L-method) Model.
Definition: fvmpfal3dinteractionvolumecontaineradaptive.hh:46
FvMpfaL3d2P2CInteractionVolumeContainerAdaptive(Problem &problem)
Definition: fvmpfal3dinteractionvolumecontaineradaptive.hh:110
int getMpfaCase6cells(const IntersectionIterator &isIt, InteractionVolume &interactionVolume, bool &properFluxDirection)
Determine the subVolumeFaceIdx for a given intersection and interactionVolume with 6 cells.
Definition: fvmpfal3dinteractionvolumecontaineradaptive.hh:286
void storeBoundaryInteractionVolume(InteractionVolume &interactionVolume, const Vertex &vertex)
Overwrites the method from the base class FvMpfaL3dInteractionVolumeContainerAdaptive On each boundar...
Definition: fvmpfal3dinteractionvolumecontaineradaptive.hh:128
int getMpfaCase2or4cells(const IntersectionIterator &isIt, InteractionVolume &interactionVolume, bool &properFluxDirection)
Determine the subVolumeFaceIdx for a given intersection and interactionVolume with 2/4 cells.
Definition: fvmpfal3dinteractionvolumecontaineradaptive.hh:425
int getMpfaCase8cells(const IntersectionIterator &isIt, const int localidxLarge, InteractionVolume &interactionVolume, bool &properFluxDirection)
Determine the subVolumeFaceIdx for a given intersection and interactionVolume with 8 cells.
Definition: fvmpfal3dinteractionvolumecontaineradaptive.hh:150