From 9365b88d8c06242610d905a6e639d6d2f6213bff Mon Sep 17 00:00:00 2001 From: Timo Koch <timo.koch@iws.uni-stuttgart.de> Date: Mon, 22 Jun 2020 17:31:45 +0200 Subject: [PATCH] [mpfa] Pass element geometry to mpfa transmissiblity function --- .../mpfa/computetransmissibility.hh | 57 +++++++++---------- .../mpfa/omethod/localassembler.hh | 6 +- .../facet/cellcentered/mpfa/localassembler.hh | 6 +- 3 files changed, 34 insertions(+), 35 deletions(-) diff --git a/dumux/discretization/cellcentered/mpfa/computetransmissibility.hh b/dumux/discretization/cellcentered/mpfa/computetransmissibility.hh index 80e5c82235..0159500ac7 100644 --- a/dumux/discretization/cellcentered/mpfa/computetransmissibility.hh +++ b/dumux/discretization/cellcentered/mpfa/computetransmissibility.hh @@ -30,34 +30,35 @@ #include <dune/common/fvector.hh> #include <dumux/common/math.hh> +#include <dumux/discretization/extrusion.hh> -namespace Dumux -{ +namespace Dumux { /*! * \ingroup CCMpfaDiscretization * \brief Free function to evaluate the Mpfa transmissibility associated - * with the flux (in the form of flux = T*gradU) across a + * with the flux (in the form of flux = t*gradU) across a * sub-control volume face stemming from a given sub-control - * volume with corresponding tensor T. + * volume with corresponding tensor t. * * \param scv The iv-local sub-control volume * \param scvf The grid sub-control volume face - * \param T The tensor living in the scv + * \param t The tensor living in the scv * \param extrusionFactor The extrusion factor of the scv */ -template<class Scv, class Scvf, class Tensor> -Dune::FieldVector< typename Tensor::field_type, Scv::myDimension > -computeMpfaTransmissibility(const Scv& scv, - const Scvf& scvf, - const Tensor& T, - typename Scv::ctype extrusionFactor) +template<class EG, class IVSubControlVolume, class Tensor> +Dune::FieldVector<typename Tensor::field_type, IVSubControlVolume::myDimension> +computeMpfaTransmissibility(const IVSubControlVolume& scv, + const typename EG::SubControlVolumeFace& scvf, + const Tensor& t, + typename IVSubControlVolume::ctype extrusionFactor) { - Dune::FieldVector< typename Tensor::field_type, Scv::myDimension > wijk; - for (unsigned int dir = 0; dir < Scv::myDimension; ++dir) - wijk[dir] = vtmv(scvf.unitOuterNormal(), T, scv.nu(dir)); + Dune::FieldVector<typename Tensor::field_type, IVSubControlVolume::myDimension> wijk; + for (unsigned int dir = 0; dir < IVSubControlVolume::myDimension; ++dir) + wijk[dir] = vtmv(scvf.unitOuterNormal(), t, scv.nu(dir)); - wijk *= scvf.area()*extrusionFactor; + using Extrusion = Extrusion_t<typename EG::GridGeometry>; + wijk *= Extrusion::area(scvf)*extrusionFactor; wijk /= scv.detX(); return wijk; @@ -66,30 +67,28 @@ computeMpfaTransmissibility(const Scv& scv, /*! * \ingroup CCMpfaDiscretization * \brief Free function to evaluate the Mpfa transmissibility associated - * with the flux (in the form of flux = T*gradU) across a + * with the flux (in the form of flux = t*gradU) across a * sub-control volume face stemming from a given sub-control - * volume with corresponding tensor T, where T is a scalar. + * volume with corresponding tensor t, where t is a scalar. * * \param scv The iv-local sub-control volume * \param scvf The grid sub-control volume face * \param t the scalar quantity living in the scv * \param extrusionFactor The extrusion factor of the scv */ -template< class Scv, - class Scvf, - class Tensor, - std::enable_if_t< Dune::IsNumber<Tensor>::value, int > = 1 > -Dune::FieldVector<Tensor, Scv::myDimension> -computeMpfaTransmissibility(const Scv& scv, - const Scvf& scvf, - Tensor t, - typename Scv::ctype extrusionFactor) +template< class EG, class IVSubControlVolume, class Tensor, std::enable_if_t< Dune::IsNumber<Tensor>::value, int > = 1 > +Dune::FieldVector<Tensor, IVSubControlVolume::myDimension> +computeMpfaTransmissibility(const IVSubControlVolume& scv, + const typename EG::SubControlVolumeFace& scvf, + const Tensor& t, + typename IVSubControlVolume::ctype extrusionFactor) { - Dune::FieldVector<Tensor, Scv::myDimension> wijk; - for (unsigned int dir = 0; dir < Scv::myDimension; ++dir) + Dune::FieldVector<Tensor, IVSubControlVolume::myDimension> wijk; + for (unsigned int dir = 0; dir < IVSubControlVolume::myDimension; ++dir) wijk[dir] = vtmv(scvf.unitOuterNormal(), t, scv.nu(dir)); - wijk *= scvf.area()*extrusionFactor; + using Extrusion = Extrusion_t<typename EG::GridGeometry>; + wijk *= Extrusion::area(scvf)*extrusionFactor; wijk /= scv.detX(); return wijk; diff --git a/dumux/discretization/cellcentered/mpfa/omethod/localassembler.hh b/dumux/discretization/cellcentered/mpfa/omethod/localassembler.hh index 52b287045f..59bc3ea932 100644 --- a/dumux/discretization/cellcentered/mpfa/omethod/localassembler.hh +++ b/dumux/discretization/cellcentered/mpfa/omethod/localassembler.hh @@ -209,7 +209,7 @@ private: // the omega factors of the "positive" sub volume Helper::resizeVector(wijk[faceIdx], /*no outside scvs present*/1); - wijk[faceIdx][0] = computeMpfaTransmissibility(posLocalScv, curGlobalScvf, tensor, posVolVars.extrusionFactor()); + wijk[faceIdx][0] = computeMpfaTransmissibility<EG>(posLocalScv, curGlobalScvf, tensor, posVolVars.extrusionFactor()); const auto posScvLocalDofIdx = posLocalScv.localDofIndex(); for (LocalIndexType localDir = 0; localDir < dim; localDir++) @@ -245,7 +245,7 @@ private: // the omega factors of the "positive" sub volume Helper::resizeVector(wijk[faceIdx], neighborScvIndices.size()); - wijk[faceIdx][0] = computeMpfaTransmissibility(posLocalScv, curGlobalScvf, tensor, posVolVars.extrusionFactor()); + wijk[faceIdx][0] = computeMpfaTransmissibility<EG>(posLocalScv, curGlobalScvf, tensor, posVolVars.extrusionFactor()); using std::abs; bool insideZeroWij = false; @@ -307,7 +307,7 @@ private: // On surface grids, use outside face for "negative" transmissibility calculation const auto& scvf = dim < dimWorld ? this->fvGeometry().flipScvf(curGlobalScvf.index(), idxInOutside) : curGlobalScvf; - wijk[faceIdx][idxOnScvf] = computeMpfaTransmissibility(negLocalScv, scvf, negTensor, negVolVars.extrusionFactor()); + wijk[faceIdx][idxOnScvf] = computeMpfaTransmissibility<EG>(negLocalScv, scvf, negTensor, negVolVars.extrusionFactor()); // flip sign on surface grids (since we used the "outside" normal) if (dim < dimWorld) diff --git a/dumux/multidomain/facet/cellcentered/mpfa/localassembler.hh b/dumux/multidomain/facet/cellcentered/mpfa/localassembler.hh index 6fc0ea5494..217b7dea49 100644 --- a/dumux/multidomain/facet/cellcentered/mpfa/localassembler.hh +++ b/dumux/multidomain/facet/cellcentered/mpfa/localassembler.hh @@ -338,7 +338,7 @@ private: // the omega factors of the "positive" sub volume Helper::resizeVector(wijk[faceIdx], /*no outside scvs present*/1); - wijk[faceIdx][0] = computeMpfaTransmissibility(posLocalScv, curGlobalScvf, tensor, posVolVars.extrusionFactor()); + wijk[faceIdx][0] = computeMpfaTransmissibility<EG>(posLocalScv, curGlobalScvf, tensor, posVolVars.extrusionFactor()); const auto posScvLocalDofIdx = posLocalScv.localDofIndex(); for (LocalIndexType localDir = 0; localDir < dim; localDir++) @@ -377,7 +377,7 @@ private: // the omega factors of the "positive" sub volume Helper::resizeVector(wijk[faceIdx], neighborScvIndices.size()); - wijk[faceIdx][0] = computeMpfaTransmissibility(posLocalScv, curGlobalScvf, tensor, posVolVars.extrusionFactor()); + wijk[faceIdx][0] = computeMpfaTransmissibility<EG>(posLocalScv, curGlobalScvf, tensor, posVolVars.extrusionFactor()); using std::abs; bool isZeroWij = false; @@ -465,7 +465,7 @@ private: // On surface grids, use outside face for "negative" transmissibility calculation const auto& scvf = dim < dimWorld ? this->fvGeometry().flipScvf(curGlobalScvf.index(), idxInOutside) : curGlobalScvf; - wijk[faceIdx][idxOnScvf] = computeMpfaTransmissibility(negLocalScv, scvf, negTensor, negVolVars.extrusionFactor()); + wijk[faceIdx][idxOnScvf] = computeMpfaTransmissibility<EG>(negLocalScv, scvf, negTensor, negVolVars.extrusionFactor()); // flip sign on surface grids (since we used the "outside" normal) if (dim < dimWorld) -- GitLab