24#ifndef DUMUX_VARIABLECLASS2P2C_ADAPTIVE_HH
25#define DUMUX_VARIABLECLASS2P2C_ADAPTIVE_HH
27#include <dune/common/power.hh>
45template<
class TypeTag>
56 using AdaptedValues =
typename CellData::AdaptedValues;
58 using Grid =
typename GridView::Grid;
59 using IdType =
typename Grid::LocalIdSet::IdType;
60 using IntersectionIterator =
typename GridView::IntersectionIterator;
61 using Intersection =
typename GridView::Intersection;
64 dim = GridView::dimension
74 using GlobalPosition = Dune::FieldVector<Scalar, GridView::dimensionworld>;
75 using TransmissivityMatrix = Dune::FieldVector<Scalar,dim+1>;
93 TransmissivityMatrix
T1_[(dim-1)*2];
162 problem.pressureModel().adaptPressure();
180 const TransmissivityMatrix& T1,
181 const GlobalPosition& globalPos3,
182 const int& globalIdx3)
184 IdType intersectionID =
grid_.localIdSet().subId(
185 irregularIs.inside(), irregularIs.indexInInside(), 1);
187 if (irregularIs.inside().level() < irregularIs.outside().level())
190 intersectionID =
grid_.localIdSet().subId(
191 irregularIs.outside(), irregularIs.indexInOutside(), 1);
223 IntersectionIterator& secondHalfEdgeIntersectionIt,
224 const TransmissivityMatrix& T1,
225 const TransmissivityMatrix& T1_secondHalfEdge,
226 const GlobalPosition& globalPos3,
227 const int& globalIdx3)
229 IdType intersectionID
230 =
grid_.localIdSet().subId(irregularIs.inside(),
231 irregularIs.indexInInside(), 1);
234 if (irregularIs.inside().level() < irregularIs.outside().level())
238 =
grid_.localIdSet().subId(irregularIs.outside(),
239 irregularIs.indexInOutside(), 1);
282 const TransmissivityMatrix& T1,
283 const GlobalPosition& globalPos3,
284 const int& globalIdx3,
285 const GlobalPosition& globalPos4,
286 const int& globalIdx4,
291 const int offset = subFaceIdx * 2;
293 IdType intersectionID
294 =
grid_.localIdSet().subId(irregularIs.inside(),
295 irregularIs.indexInInside(), 1);
298 if (irregularIs.inside().level() < irregularIs.outside().level())
302 =
grid_.localIdSet().subId(irregularIs.outside(),
303 irregularIs.indexInOutside(), 1);
343 Scalar weight,
int subFaceIdx = -1)
345 IdType intersectionID
346 =
grid_.localIdSet().subId(irregularIs.inside(),
347 irregularIs.indexInInside(), 1);
350 if (irregularIs.inside().level() < irregularIs.outside().level())
354 =
grid_.localIdSet().subId(irregularIs.outside(),
355 irregularIs.indexInOutside(), 1);
385 IntersectionIterator& secondHalfEdgeIntersectionIt,
386 TransmissivityMatrix& T1,
387 TransmissivityMatrix& T1_secondHalfEdge,
388 GlobalPosition& globalPos3,
391 IdType intersectionID
392 =
grid_.localIdSet().subId(irregularIs.inside(),
393 irregularIs.indexInInside(), 1);
395 if (irregularIs.inside().level() < irregularIs.outside().level())
399 =
grid_.localIdSet().subId(irregularIs.outside(),
400 irregularIs.indexInOutside(), 1);
462 TransmissivityMatrix& T1,
463 GlobalPosition& globalPos3,
465 GlobalPosition& globalPos4,
471 const int offset = subFaceIdx * 2;
473 IdType intersectionID
474 =
grid_.localIdSet().subId(irregularIs.inside(),
475 irregularIs.indexInInside(), 1);
477 if (irregularIs.inside().level() < irregularIs.outside().level())
481 =
grid_.localIdSet().subId(irregularIs.outside(),
482 irregularIs.indexInOutside(), 1);
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
Base class holding the variables and discretized data for sequential models.
Definition: variableclass.hh:49
const GridView & gridView() const
Return gridView.
Definition: variableclass.hh:156
void adaptVariableSize(const int size)
Resizes sequential variable vectors.
Definition: variableclass.hh:102
Class holding additionally mpfa data of adaptive compositional models.
Definition: 2p2c/sequential/variableclassadaptive.hh:47
VariableClass2P2CAdaptive(const GridView &gridView)
Constructs a VariableClass object.
Definition: 2p2c/sequential/variableclassadaptive.hh:132
const Grid & grid_
The grid.
Definition: 2p2c/sequential/variableclassadaptive.hh:124
int getMpfaData3D(const Intersection &irregularIs, TransmissivityMatrix &T1, GlobalPosition &globalPos3, int &globalIdx3, GlobalPosition &globalPos4, int &globalIdx4, int subFaceIdx=0)
Provides access to stored 3D Mpfa Data on an intersection for up to 4 subfaces.
Definition: 2p2c/sequential/variableclassadaptive.hh:461
int getMpfaData(const Intersection &irregularIs, IntersectionIterator &secondHalfEdgeIntersectionIt, TransmissivityMatrix &T1, TransmissivityMatrix &T1_secondHalfEdge, GlobalPosition &globalPos3, int &globalIdx3)
Provides access to stored Mpfa Data on an intersection for both half-edges.
Definition: 2p2c/sequential/variableclassadaptive.hh:384
std::map< IdType, mpfaData > irregularInterfaceMap_
Storage container for mpfa information.
Definition: 2p2c/sequential/variableclassadaptive.hh:123
void reconstructPrimVars(Problem &problem)
Reconstruct missing primary variables (where elements are created/deleted)
Definition: 2p2c/sequential/variableclassadaptive.hh:158
void storeMpfaData(const typename GridView::Intersection &irregularIs, IntersectionIterator &secondHalfEdgeIntersectionIt, const TransmissivityMatrix &T1, const TransmissivityMatrix &T1_secondHalfEdge, const GlobalPosition &globalPos3, const int &globalIdx3)
Stores Mpfa Data on an intersection for both half-edges.
Definition: 2p2c/sequential/variableclassadaptive.hh:222
void adaptVariableSize(int size)
Resizes sequential variable vectors.
Definition: 2p2c/sequential/variableclassadaptive.hh:144
void storeMpfaData(const typename GridView::Intersection &irregularIs, const TransmissivityMatrix &T1, const GlobalPosition &globalPos3, const int &globalIdx3)
Stores Mpfa Data of one interaction region on an intersection.
Definition: 2p2c/sequential/variableclassadaptive.hh:179
void performTransmissitivityWeighting(const typename GridView::Intersection &irregularIs, Scalar weight, int subFaceIdx=-1)
Weigths the transmissivity coefficient by the flux area (3D)
Definition: 2p2c/sequential/variableclassadaptive.hh:342
static const int storageRequirement
Definition: 2p2c/sequential/variableclassadaptive.hh:82
void storeMpfaData3D(const typename GridView::Intersection &irregularIs, const TransmissivityMatrix &T1, const GlobalPosition &globalPos3, const int &globalIdx3, const GlobalPosition &globalPos4, const int &globalIdx4, int subFaceIdx=0)
Stores 3D Mpfa Data on an intersection.
Definition: 2p2c/sequential/variableclassadaptive.hh:281
Storage object for data related to the MPFA method.
Definition: 2p2c/sequential/variableclassadaptive.hh:92
mpfaData()
Constructor for the local storage object of mpfa data.
Definition: 2p2c/sequential/variableclassadaptive.hh:100
void setIntersection(IntersectionIterator &is23)
Stores an intersection for the 2D implementation.
Definition: 2p2c/sequential/variableclassadaptive.hh:112
int interactionRegionsStored
Definition: 2p2c/sequential/variableclassadaptive.hh:97
const IntersectionIterator & getIntersection()
Acess method to the stored intersection (only 2D)
Definition: 2p2c/sequential/variableclassadaptive.hh:118
TransmissivityMatrix T1_[(dim-1) *2]
Definition: 2p2c/sequential/variableclassadaptive.hh:93
GlobalPosition globalPos3_[storageRequirement]
Definition: 2p2c/sequential/variableclassadaptive.hh:95
int globalIdx3_[storageRequirement]
Definition: 2p2c/sequential/variableclassadaptive.hh:94
std::vector< IntersectionIterator > secondHalfEdgeIntersection_
Definition: 2p2c/sequential/variableclassadaptive.hh:96
Base class holding the variables and discretized data for sequential models.
Definition: sequential/variableclassadaptive.hh:47
void reconstructPrimVars(const Problem &problem)
Definition: sequential/variableclassadaptive.hh:140
Base class holding the variables for sequential models.