Skip to content
Snippets Groups Projects
Commit cd916a72 authored by Bernd Flemisch's avatar Bernd Flemisch
Browse files

2pfluxvariables, boxfvelementgeometry: add possibility for generalizing

the flux approximation.
2plocalresidual: derive from the already existing property
BaseLocalResidual instead of explicitly deriving from BoxLocalResidual.

Reviewed by Christoph.


git-svn-id: svn://svn.iws.uni-stuttgart.de/DUMUX/dumux/trunk@8062 2fb0f335-1f38-0410-981e-8018bf24f1b0
parent 4267b202
No related branches found
No related tags found
No related merge requests found
...@@ -157,18 +157,21 @@ private: ...@@ -157,18 +157,21 @@ private:
{ {
// calculate gradients // calculate gradients
for (int idx = 0; for (int idx = 0;
idx < fvElemGeom_.numVertices; idx < fvElemGeom_.numFAP;
idx++) // loop over adjacent vertices idx++) // loop over adjacent vertices
{ {
// FE gradient at vertex idx // FE gradient at vertex idx
const Vector &feGrad = face().grad[idx]; const Vector &feGrad = face().grad[idx];
// index for the element volume variables
int volVarsIdx = face().fapIndices[idx];
// compute sum of pressure gradients for each phase // compute sum of pressure gradients for each phase
for (int phase = 0; phase < numPhases; phase++) for (int phase = 0; phase < numPhases; phase++)
{ {
// the pressure gradient // the pressure gradient
Vector tmp(feGrad); Vector tmp(feGrad);
tmp *= elemVolVars[idx].pressure(phase); tmp *= elemVolVars[volVarsIdx].pressure(phase);
potentialGrad_[phase] += tmp; potentialGrad_[phase] += tmp;
} }
} }
......
...@@ -44,7 +44,7 @@ namespace Dumux ...@@ -44,7 +44,7 @@ namespace Dumux
* that it uses static polymorphism. * that it uses static polymorphism.
*/ */
template<class TypeTag> template<class TypeTag>
class TwoPLocalResidual : public BoxLocalResidual<TypeTag> class TwoPLocalResidual : public GET_PROP_TYPE(TypeTag, BaseLocalResidual)
{ {
protected: protected:
typedef typename GET_PROP_TYPE(TypeTag, LocalResidual) Implementation; typedef typename GET_PROP_TYPE(TypeTag, LocalResidual) Implementation;
......
...@@ -325,6 +325,7 @@ class BoxFVElementGeometry ...@@ -325,6 +325,7 @@ class BoxFVElementGeometry
enum{maxNF = (dim < 3 ? 1 : 6)}; enum{maxNF = (dim < 3 ? 1 : 6)};
enum{maxCOS = (dim < 3 ? 2 : 4)}; enum{maxCOS = (dim < 3 ? 2 : 4)};
enum{maxBF = (dim < 3 ? 8 : 24)}; enum{maxBF = (dim < 3 ? 8 : 24)};
enum{maxNFAP = (dim < 3 ? 4 : 8)};
typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar; typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
typedef typename GridView::ctype CoordScalar; typedef typename GridView::ctype CoordScalar;
typedef typename GridView::Traits::template Codim<0>::Entity Element; typedef typename GridView::Traits::template Codim<0>::Entity Element;
...@@ -571,6 +572,7 @@ public: ...@@ -571,6 +572,7 @@ public:
Scalar area; //!< area of face Scalar area; //!< area of face
Dune::FieldVector<Vector, maxNC> grad; //!< derivatives of shape functions at ip Dune::FieldVector<Vector, maxNC> grad; //!< derivatives of shape functions at ip
Dune::FieldVector<Scalar, maxNC> shapeValue; //!< value of shape functions at ip Dune::FieldVector<Scalar, maxNC> shapeValue; //!< value of shape functions at ip
Dune::FieldVector<int, maxNFAP> fapIndices; //!< indices w.r.t.neighbors of the flux approximation points
}; };
typedef SubControlVolumeFace BoundaryFace; //!< compatibility typedef typedef SubControlVolumeFace BoundaryFace; //!< compatibility typedef
...@@ -586,6 +588,7 @@ public: ...@@ -586,6 +588,7 @@ public:
int numVertices; //!< number of verts int numVertices; //!< number of verts
int numEdges; //!< number of edges int numEdges; //!< number of edges
int numFaces; //!< number of faces (0 in < 3D) int numFaces; //!< number of faces (0 in < 3D)
int numFAP; //!< number of flux approximation points
const LocalFiniteElementCache feCache_; const LocalFiniteElementCache feCache_;
bool computeGradientAtScvCenters; bool computeGradientAtScvCenters;
...@@ -618,6 +621,7 @@ public: ...@@ -618,6 +621,7 @@ public:
numVertices = referenceElement.size(dim); numVertices = referenceElement.size(dim);
numEdges = referenceElement.size(dim-1); numEdges = referenceElement.size(dim-1);
numFaces = (dim<3)?0:referenceElement.size(1); numFaces = (dim<3)?0:referenceElement.size(1);
numFAP = numVertices;
// corners: // corners:
for (int vert = 0; vert < numVertices; vert++) { for (int vert = 0; vert < numVertices; vert++) {
...@@ -713,6 +717,7 @@ public: ...@@ -713,6 +717,7 @@ public:
for (int vert = 0; vert < numVertices; vert++) { for (int vert = 0; vert < numVertices; vert++) {
jacInvT.mv(localJac[vert][0], subContVolFace[k].grad[vert]); jacInvT.mv(localJac[vert][0], subContVolFace[k].grad[vert]);
subContVolFace[k].shapeValue[vert] = Scalar(shapeVal[vert]); subContVolFace[k].shapeValue[vert] = Scalar(shapeVal[vert]);
subContVolFace[k].fapIndices[vert] = vert;
} }
} // end loop over edges / sub control volume faces } // end loop over edges / sub control volume faces
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment