24#ifndef DUMUX_FVMPFAL3D_2P2CINTERACTIONVOLUMECONTAINER_ADAPTIVE_HH
25#define DUMUX_FVMPFAL3D_2P2CINTERACTIONVOLUMECONTAINER_ADAPTIVE_HH
44template<
class TypeTag>
52 dim = GridView::dimension, dimWorld = GridView::dimensionworld
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;
68 using IntersectionIterator =
typename GridView::IntersectionIterator;
69 using Intersection =
typename GridView::Intersection;
70 using IntersectionGeometry =
typename Intersection::Geometry;
72 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
73 using DimMatrix = Dune::FieldMatrix<Scalar, dim, dim>;
75 using DimVector = Dune::FieldVector<Scalar, dim>;
79 pressureEqIdx = Indices::pressureEqIdx,
83 innerEdgeFace = 2, innerSideFace = 1
87 realFaceArea = 0, fluxFaceArea = 1
94 using GlobalInteractionVolumeVector = std::vector<InteractionVolume>;
95 using FaceAreaVector = std::vector<Dune::FieldVector<Dune::FieldVector<Scalar, 2>, 2*dim> >;
100 const int localidxLarge,
102 bool& properFluxDirection);
105 bool& properFluxDirection);
108 bool& properFluxDirection);
115 DUNE_THROW(Dune::NotImplemented,
"Dimension not supported!");
127template<
class TypeTag>
149template<
class TypeTag>
151 const int localIdxLarge,
153 bool& properFluxDirection)
155 int subVolumeFaceIdx = -1;
158 int localFaceIdxInside = isIt->indexInInside();
159 switch(localIdxLarge)
162 if(localFaceIdxInside == 0)
164 subVolumeFaceIdx = 6;
165 properFluxDirection =
false;
167 else if(localFaceIdxInside == 2)
168 subVolumeFaceIdx = 5;
169 else if(localFaceIdxInside == 4)
170 subVolumeFaceIdx = 10;
174 if(localFaceIdxInside == 1)
175 subVolumeFaceIdx = 6;
176 else if(localFaceIdxInside == 2)
178 subVolumeFaceIdx = 7;
179 properFluxDirection =
false;
181 else if(localFaceIdxInside == 4)
183 subVolumeFaceIdx = 11;
184 properFluxDirection =
false;
189 if(localFaceIdxInside == 0)
190 subVolumeFaceIdx = 4;
191 else if(localFaceIdxInside == 3)
193 subVolumeFaceIdx = 5;
194 properFluxDirection =
false;
196 else if(localFaceIdxInside == 4)
198 subVolumeFaceIdx = 9;
199 properFluxDirection =
false;
204 if(localFaceIdxInside == 1)
206 subVolumeFaceIdx = 4;
207 properFluxDirection =
false;
209 else if(localFaceIdxInside == 3)
210 subVolumeFaceIdx = 7;
211 else if(localFaceIdxInside == 4)
212 subVolumeFaceIdx = 8;
216 if(localFaceIdxInside == 0)
217 subVolumeFaceIdx = 2;
218 else if(localFaceIdxInside == 2)
220 subVolumeFaceIdx = 1;
221 properFluxDirection =
false;
223 else if(localFaceIdxInside == 5)
225 subVolumeFaceIdx = 10;
226 properFluxDirection =
false;
231 if(localFaceIdxInside == 1)
233 subVolumeFaceIdx = 2;
234 properFluxDirection =
false;
236 else if(localFaceIdxInside == 2)
237 subVolumeFaceIdx = 3;
238 else if(localFaceIdxInside == 5)
239 subVolumeFaceIdx = 11;
243 if(localFaceIdxInside == 0)
245 subVolumeFaceIdx = 0;
246 properFluxDirection =
false;
248 else if(localFaceIdxInside == 3)
249 subVolumeFaceIdx = 1;
250 else if(localFaceIdxInside == 5)
251 subVolumeFaceIdx = 9;
255 if(localFaceIdxInside == 1)
256 subVolumeFaceIdx = 0;
257 else if(localFaceIdxInside == 3)
259 subVolumeFaceIdx = 3;
260 properFluxDirection =
false;
262 else if(localFaceIdxInside == 5)
264 subVolumeFaceIdx = 8;
265 properFluxDirection =
false;
270 return subVolumeFaceIdx;
285template<
class TypeTag>
288 bool& properFluxDirection)
290 int mapI(-5), mapJ(-5);
291 int mapI2(-5), mapJ2(-5);
293 for (
int i= 0; i<interactionVolume.getElementNumber(); i++)
295 if(isIt->inside() == interactionVolume.getSubVolumeElement(i))
303 if(isIt->outside() == interactionVolume.getSubVolumeElement(i))
312 for(
int passionfruit = 0; passionfruit <=1; passionfruit++)
314 if(mapI== 0 || mapJ == 0)
319 properFluxDirection =
false;
325 properFluxDirection =
false;
330 if(mapI== 1 || mapJ == 1)
333 properFluxDirection =
false;
339 properFluxDirection =
false;
345 else if(mapI== 2 || mapJ == 2)
348 properFluxDirection =
false;
354 else if(mapI== 3 || mapJ == 3)
357 properFluxDirection =
false;
363 else if(mapI== 4 || mapJ == 4)
368 properFluxDirection =
false;
374 properFluxDirection =
false;
380 else if(mapI== 5 || mapJ == 5)
383 properFluxDirection =
false;
389 else if(mapI== 6 || mapJ == 6)
392 properFluxDirection =
false;
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;
424template<
class TypeTag>
427 bool& properFluxDirection)
429 int mapI(-5), mapJ(-5);
430 int mapI2(-5), mapJ2(-5);
432 for (
int i= 0; i<interactionVolume.getElementNumber(); i++)
434 if(isIt->inside() == interactionVolume.getSubVolumeElement(i))
442 if(isIt->outside() == interactionVolume.getSubVolumeElement(i))
451 for(
int passionfruit = 0; passionfruit <=1; passionfruit++)
454 if(mapI== 0 || mapJ == 0)
457 properFluxDirection =
false;
463 properFluxDirection =
false;
470 properFluxDirection =
false;
474 else if(mapI== 1 || mapJ == 1)
477 properFluxDirection =
false;
483 properFluxDirection =
false;
489 else if(mapI== 2 || mapJ == 2)
492 properFluxDirection =
false;
498 properFluxDirection =
false;
504 else if(mapI== 3 || mapJ == 3)
507 properFluxDirection =
false;
513 else if(mapI== 4 || mapJ == 4)
516 properFluxDirection =
false;
522 properFluxDirection =
false;
528 else if(mapI== 5 || mapJ == 5)
531 properFluxDirection =
false;
537 else if(mapI== 6 || mapJ == 6)
540 properFluxDirection =
false;
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;
Interactionvolume container for 3-d MPFA L-method on an h-adaptive grid.
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