20#ifndef DUMUX_FVMPFAL3DINTERACTIONVOLUME_HH
21#define DUMUX_FVMPFAL3DINTERACTIONVOLUME_HH
41 subVolumeTotalNum = 8,
42 fluxFacesTotalNum = 12,
43 fluxFacesNumOnSubVolume = 3,
44 fluxEdgesTotalNum = 6,
45 edgesNumOnFluxFace = 2
48 static int getFaceIndexFromSubVolume(
int subVolumeIdx,
int subVolumeFaceIdx)
50 return faceIndexOnSubVolume_[subVolumeIdx][subVolumeFaceIdx];
53 static int getEdgeIndexFromSubVolumeFace(
int subVolumeIdx,
int subVolumeFaceIdx,
int subVolumeEdgeIdx)
55 return edgeIndexOnSubVolumeFace_[subVolumeIdx][subVolumeFaceIdx][subVolumeEdgeIdx];
58 static int getFaceIndexFromElements(
int elem1Idx,
int elem2Idx)
60 return faceIndexFromElements_[elem1Idx][elem2Idx];
65 static const int faceIndexOnSubVolume_[subVolumeTotalNum][fluxFacesNumOnSubVolume];
66 static const int edgeIndexOnSubVolumeFace_[subVolumeTotalNum][fluxFacesNumOnSubVolume][edgesNumOnFluxFace];
67 static const int faceIndexFromElements_[subVolumeTotalNum][subVolumeTotalNum];
70const int IndexTranslator::faceIndexOnSubVolume_[subVolumeTotalNum][fluxFacesNumOnSubVolume] =
82const int IndexTranslator::edgeIndexOnSubVolumeFace_[subVolumeTotalNum][fluxFacesNumOnSubVolume][edgesNumOnFluxFace] =
84 {{2, 0},{0, 5},{5, 2}},
85 {{3, 0},{0, 2},{2, 3}},
86 {{5, 0},{0, 4},{4, 5}},
87 {{4, 0},{0, 3},{3, 4}},
88 {{2, 5},{1, 2},{5, 1}},
89 {{3, 2},{1, 3},{2, 1}},
90 {{5, 4},{1, 5},{4, 1}},
91 {{4, 3},{1, 4},{3, 1}}
94const int IndexTranslator::faceIndexFromElements_[subVolumeTotalNum][subVolumeTotalNum] =
96 {-1, 0, 3, -1, 8, -1, -1, -1},
97 {0, -1, -1, 1, -1, 9, -1, -1},
98 {3, -1, -1, 2, -1, -1, 11, -1},
99 {-1, 1, 2, -1, -1, -1, -1, 10},
100 {8, -1, -1, -1, -1, 4, 7, -1},
101 {-1, 9, -1, -1, 4, -1, -1, 5},
102 {-1, -1, 11, -1, 7, -1, -1, 6},
103 {-1, -1, -1, 10, -1, 5, 6, -1}
113template<
class TypeTag>
124 dim = GridView::dimension,
125 dimWorld = GridView::dimensionworld,
128 using Element =
typename GridView::template Codim<0>::Entity;
129 using ElementSeed =
typename Grid::template Codim<0>::EntitySeed;
131 using BoundaryTypes =
typename GET_PROP_TYPE(TypeTag, BoundaryTypes);
132 using SolutionTypes =
typename GET_PROP(TypeTag, SolutionTypes);
133 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
135 using DimVector = Dune::FieldVector<Scalar, dim>;
136 using FieldVectorVector = Dune::FieldVector<DimVector, dim>;
137 using FieldVectorVector2 = Dune::FieldVector<DimVector, 2>;
138 using FieldVectorVectorVector = Dune::FieldVector<FieldVectorVector2, dim>;
139 using IndexVector = Dune::FieldVector<int, dim>;
140 using BCTypeVector = std::vector<BoundaryTypes>;
141 using BCVector = std::vector<PrimaryVariables>;
165 , normal_(FieldVectorVector(DimVector(0.0)))
166 , facePos_(DimVector(0.0))
167 , edgePos_((DimVector(0.0)))
170 , indexOnElement_(IndexVector(0.0))
172 , centerVertexPos_(0)
181 centerVertexPos_ = 0;
182 indexOnElement_ = IndexVector(0.0);
186 facePos_ = DimVector(0.0);
187 edgePos_ = DimVector(0.0);
188 normal_ = FieldVectorVector(DimVector(0.0));
190 boundaryTypes_.clear();
191 neumannValues_.clear();
192 dirichletValues_.clear();
201 centerVertexPos_ = centerVertexPos;
212 elements_[subVolumeIdx].push_back(element.seed());
217 elements_[subVolumeIdx].insert(elements_[subVolumeIdx].begin(), element.seed());
228 facePos_[fIdx] = pos;
238 edgePos_[edgeIdx] = pos;
248 faceArea_[fIdx] = faceArea;
257 void setNormal(DimVector& normal,
int subVolumeIdx,
int subVolumeFaceIdx)
259 normal_[subVolumeIdx][subVolumeFaceIdx] = normal;
267 void setBoundary(BoundaryTypes& boundaryTypes,
int subVolumeFaceIdx)
269 if (boundaryTypes_.size() == 0)
273 boundaryTypes_[subVolumeFaceIdx] = boundaryTypes;
274 faceType_[subVolumeFaceIdx] =
boundary;
283 faceType_[subVolumeFaceIdx] =
outside;
293 if (dirichletValues_.size() == 0)
297 dirichletValues_[subVolumeFaceIdx] = condition;
307 if (neumannValues_.size() == 0)
311 neumannValues_[subVolumeFaceIdx] = condition;
322 indexOnElement_[subVolumeIdx][subVolumeFaceIdx] = indexInInside;
328 return centerVertexPos_;
340 return IndexTranslator::getFaceIndexFromSubVolume(subVolumeIdx, subVolumeFaceIdx);
353 return IndexTranslator::getEdgeIndexFromSubVolumeFace(subVolumeIdx, subVolumeFaceIdx, subVolumeEdgeIdx);
371 return indexOnElement_[subVolumeIdx][subVolumeFaceIdx];
383 return grid_->entity(elements_[subVolumeIdx][0]);
386 std::cout<<
"Problems when calling getSubVolumeElement("<<subVolumeIdx<<
")\n";
388 DUNE_THROW(Dune::RangeError,
"element not in interaction volume!");
400 return (elements_[subVolumeIdx].size() > 0);
411 return boundaryTypes_[subVolumeFaceIdx];
422 return (faceType_[subVolumeFaceIdx] ==
outside);
433 return (faceType_[subVolumeFaceIdx] ==
inside);
444 return (faceType_[subVolumeFaceIdx] ==
boundary);
453 for (
unsigned int i = 0; i < faceType_.size(); i++)
467 return (boundaryTypes_.size() > 0);
478 return dirichletValues_[subVolumeFaceIdx];
489 return neumannValues_[subVolumeFaceIdx];
499 DimVector&
getNormal(
int subVolumeIdx,
int subVolumeFaceIdx)
501 return normal_[subVolumeIdx][subVolumeFaceIdx];
513 return facePos_[IndexTranslator::getFaceIndexFromSubVolume(subVolumeIdx, subVolumeFaceIdx)];
524 DimVector&
getEdgePosition(
int subVolumeIdx,
int subVolumeFaceIdx,
int subVolumeEdgeIdx)
526 return edgePos_[IndexTranslator::getEdgeIndexFromSubVolumeFace(subVolumeIdx, subVolumeFaceIdx, subVolumeEdgeIdx)];
538 return faceArea_[IndexTranslator::getFaceIndexFromSubVolume(subVolumeIdx, subVolumeFaceIdx)];
549 return facePos_[fIdx];
560 return edgePos_[edgeIdx];
571 return faceArea_[fIdx];
577 std::cout<<
"\nNumber of stored elements: "<<elementNum_<<
"\n";
578 std::cout<<
" center position: "<<centerVertexPos_<<
"\n";
581 if (elements_[i].size() > 0)
583 std::cout<<
"element "<<i<<
":\n";
584 std::cout<<
"element level: "<<grid_->entity(elements_[i][0]).level()<<
"\n";
585 std::cout<<
"element position: "<<grid_->entity(elements_[i][0]).geometry().center()<<
"\n";
586 std::cout<<
"element volume: "<<grid_->entity(elements_[i][0]).geometry().volume()<<
"\n";
587 std::cout<<
"face indices on element: "<<indexOnElement_[i]<<
"\n";
588 std::cout<<
"face normals on element: "<<normal_[i]<<
"\n";
589 std::cout<<
"face areas on element: ";
590 for (
int j = 0; j < 3; j++)
595 std::cout<<
"face position on element: ";
596 for (
int j = 0; j < 3; j++)
601 std::cout<<
"edgePos:";
602 for (
int j = 0; j < 3; j++)
605 for (
int k = 0; k < 2; k++)
690 Dune::FieldVector<FieldVectorVector, subVolumeTotalNum> normal_;
691 Dune::FieldVector<DimVector, fluxFacesTotalNum> facePos_;
692 Dune::FieldVector<DimVector, fluxEdgesTotalNum> edgePos_;
693 Dune::FieldVector<Scalar, fluxFacesTotalNum> faceArea_;
694 BCTypeVector boundaryTypes_;
695 std::vector<int> faceType_;
696 Dune::FieldVector<IndexVector, subVolumeTotalNum> indexOnElement_;
697 std::vector<std::vector<ElementSeed> > elements_;
698 BCVector neumannValues_;
699 BCVector dirichletValues_;
700 DimVector centerVertexPos_;
#define GET_PROP(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:281
#define GET_PROP_TYPE(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:283
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Class including the information of a 3d interaction volume of a MPFA L-method that does not change wi...
Definition: linteractionvolume3d.hh:115
int getIndexOnElement(int subVolumeIdx, int subVolumeFaceIdx)
The local dune face index dependent on the local interaction volume indices.
Definition: linteractionvolume3d.hh:369
DimVector & getEdgePosition(int edgeIdx)
The position of the edge center.
Definition: linteractionvolume3d.hh:558
bool isInsideFace(int subVolumeFaceIdx)
Check if an interaction volume face is inside the model domain (for boundary vertices)
Definition: linteractionvolume3d.hh:431
void setEdgePosition(const DimVector &pos, int edgeIdx)
Store the center of the edges.
Definition: linteractionvolume3d.hh:236
void setSubVolumeElement(const Element &element, int subVolumeIdx)
Store a dune element as a sub volume element.
Definition: linteractionvolume3d.hh:208
bool isBoundaryFace(int subVolumeFaceIdx)
Check if an interaction volume face is a model domain boundary (for boundary vertices)
Definition: linteractionvolume3d.hh:442
BoundaryTypes & getBoundaryType(int subVolumeFaceIdx)
The boundary types.
Definition: linteractionvolume3d.hh:409
void setNormal(DimVector &normal, int subVolumeIdx, int subVolumeFaceIdx)
Store the normals.
Definition: linteractionvolume3d.hh:257
PrimaryVariables & getDirichletValues(int subVolumeFaceIdx)
The Dirichlet boundary values.
Definition: linteractionvolume3d.hh:476
Element getSubVolumeElement(int subVolumeIdx)
The dune element of the interaction volume sub-volume.
Definition: linteractionvolume3d.hh:380
void setBoundary(BoundaryTypes &boundaryTypes, int subVolumeFaceIdx)
Store the types of boundary conditions.
Definition: linteractionvolume3d.hh:267
PrimaryVariables & getNeumannValues(int subVolumeFaceIdx)
The Neumann boundary values.
Definition: linteractionvolume3d.hh:487
void setIndexOnElement(int indexInInside, int subVolumeIdx, int subVolumeFaceIdx)
Store the local dune face index dependent on the local interaction volume indices.
Definition: linteractionvolume3d.hh:320
@ fluxFacesTotalNum
Number of flux faces in the interaction volume.
Definition: linteractionvolume3d.hh:156
@ subVolumeTotalNum
Number of sub-volumes in the interaction volume.
Definition: linteractionvolume3d.hh:155
@ fluxEdgesTotalNum
Number of edges in the interaction volume.
Definition: linteractionvolume3d.hh:157
DimVector & getEdgePosition(int subVolumeIdx, int subVolumeFaceIdx, int subVolumeEdgeIdx)
The position of the edge center.
Definition: linteractionvolume3d.hh:524
void setFaceArea(Scalar faceArea, int fIdx)
Store the flux face areas.
Definition: linteractionvolume3d.hh:246
void setNeumannCondition(PrimaryVariables &condition, int subVolumeFaceIdx)
Store Neumann boundary conditions.
Definition: linteractionvolume3d.hh:305
bool isInnerVolume()
Check if an interaction volume is not located around a boundary vertex.
Definition: linteractionvolume3d.hh:451
DimVector & getNormal(int subVolumeIdx, int subVolumeFaceIdx)
The face normal.
Definition: linteractionvolume3d.hh:499
int getFaceIndexFromSubVolume(int subVolumeIdx, int subVolumeFaceIdx)
The local interaction volume face index dependent on the local sub volume indices.
Definition: linteractionvolume3d.hh:338
void printInteractionVolumeInfo()
Print the stored interaction volume data.
Definition: linteractionvolume3d.hh:575
DimVector & getCenterPosition()
The position of the central vertex.
Definition: linteractionvolume3d.hh:326
Scalar & getFaceArea(int fIdx)
The interaction volume flux face area.
Definition: linteractionvolume3d.hh:569
bool hasSubVolumeElement(int subVolumeIdx)
Check if a dune element is stored for an interaction volume sub-volume.
Definition: linteractionvolume3d.hh:398
bool isOutsideFace(int subVolumeFaceIdx)
Check if an interaction volume face is outside the model domain (for boundary vertices)
Definition: linteractionvolume3d.hh:420
int getEdgeIndexFromSubVolumeFace(int subVolumeIdx, int subVolumeFaceIdx, int subVolumeEdgeIdx)
The local interaction volume edge index dependent on the local sub volume indices.
Definition: linteractionvolume3d.hh:351
DimVector & getFacePosition(int fIdx)
The position of the face center.
Definition: linteractionvolume3d.hh:547
void setDirichletCondition(PrimaryVariables &condition, int subVolumeFaceIdx)
Store Dirichlet boundary conditions.
Definition: linteractionvolume3d.hh:291
FaceTypes
Definition: linteractionvolume3d.hh:147
@ outside
Flux face is outside the model domain.
Definition: linteractionvolume3d.hh:150
@ inside
Flux face is inside the model domain.
Definition: linteractionvolume3d.hh:148
@ boundary
Flux face is a boundary face.
Definition: linteractionvolume3d.hh:149
Scalar & getFaceArea(int subVolumeIdx, int subVolumeFaceIdx)
The interaction volume flux face area.
Definition: linteractionvolume3d.hh:536
int getElementNumber()
Number of stored dune elements.
Definition: linteractionvolume3d.hh:357
FvMpfaL3dInteractionVolume(const Grid &grid)
Constructs a FvMpfaL3dInteractionVolume object.
Definition: linteractionvolume3d.hh:163
bool isBoundaryInteractionVolume()
Check if an interaction volume is located around a boundary vertex.
Definition: linteractionvolume3d.hh:465
DimVector & getFacePosition(int subVolumeIdx, int subVolumeFaceIdx)
The position of the face center.
Definition: linteractionvolume3d.hh:511
void setCenterPosition(const DimVector ¢erVertexPos)
Store position of the central vertex.
Definition: linteractionvolume3d.hh:199
void reset()
Reset the interaction volume (deletes stored data)
Definition: linteractionvolume3d.hh:177
void setOutsideFace(int subVolumeFaceIdx)
Define a flux face to be outside the model domain (for boundary vertices)
Definition: linteractionvolume3d.hh:281
void setFacePosition(const DimVector &pos, int fIdx)
Store the position of a flux face.
Definition: linteractionvolume3d.hh:226
Base file for properties related to sequential models.