From 52dfa0e2db7dbe2a1b27080e08da6ad2020d61a8 Mon Sep 17 00:00:00 2001 From: Kilian Weishaupt <kilian.weishaupt@iws.uni-stuttgart.de> Date: Sat, 9 Dec 2017 11:59:02 +0100 Subject: [PATCH] [navierstokes] Make indices independent of discretization * clean up --- dumux/assembly/staggeredlocalassembler.hh | 5 +++-- dumux/freeflow/navierstokes/indices.hh | 22 +++++++++++--------- dumux/implicit/staggered/primaryvariables.hh | 10 +++++---- 3 files changed, 21 insertions(+), 16 deletions(-) diff --git a/dumux/assembly/staggeredlocalassembler.hh b/dumux/assembly/staggeredlocalassembler.hh index f2aeea6500..9de574c3a7 100644 --- a/dumux/assembly/staggeredlocalassembler.hh +++ b/dumux/assembly/staggeredlocalassembler.hh @@ -84,6 +84,7 @@ class StaggeredLocalAssembler<TypeTag, using SubControlVolumeFace = typename GET_PROP_TYPE(TypeTag, SubControlVolumeFace); static constexpr bool enableGlobalFluxVarsCache = GET_PROP_VALUE(TypeTag, EnableGlobalFluxVariablesCache); + static constexpr auto faceOffset = GET_PROP_VALUE(TypeTag, NumEqCellCenter); public: @@ -402,7 +403,7 @@ static void dCCdFace_(Assembler& assembler, partialDeriv /= eps; // update the global jacobian matrix with the current partial derivatives - updateGlobalJacobian_(matrix[cellCenterIdx][faceIdx], cellCenterGlobalI, globalJ, pvIdx - Indices::faceOffset, partialDeriv); + updateGlobalJacobian_(matrix[cellCenterIdx][faceIdx], cellCenterGlobalI, globalJ, pvIdx - faceOffset, partialDeriv); // restore the original faceVars faceVars = origFaceVars; @@ -532,7 +533,7 @@ static void dFacedFace_(Assembler& assembler, partialDeriv /= eps; // update the global jacobian matrix with the current partial derivatives - updateGlobalJacobian_(matrix[faceIdx][faceIdx], faceGlobalI, globalJ, pvIdx - Indices::faceOffset, partialDeriv); + updateGlobalJacobian_(matrix[faceIdx][faceIdx], faceGlobalI, globalJ, pvIdx - faceOffset, partialDeriv); // restore the original faceVars faceVars = origFaceVars; diff --git a/dumux/freeflow/navierstokes/indices.hh b/dumux/freeflow/navierstokes/indices.hh index 39fec61a10..165ec691c4 100644 --- a/dumux/freeflow/navierstokes/indices.hh +++ b/dumux/freeflow/navierstokes/indices.hh @@ -18,7 +18,7 @@ *****************************************************************************/ /*! * \file - * \brief Defines the indices for the one-phase fully implicit model. + * \brief Defines the indices for the one-phase isothermal Navier-Stoke model. */ #ifndef DUMUX_NAVIERSTOKES_COMMON_INDICES_HH #define DUMUX_NAVIERSTOKES_COMMON_INDICES_HH @@ -30,8 +30,7 @@ namespace Dumux // \{ /*! * \ingroup NavierStokesModel - * \ingroup ImplicitIndices - * \brief The common indices for the isothermal stokes model. + * \brief The common indices for the isothermal Navier-Stoke model. * * \tparam PVOffset The first index in a primary variable vector. */ @@ -47,17 +46,20 @@ struct NavierStokesCommonIndices static constexpr int conti0EqIdx = massBalanceIdx; //!< Index of first (for C-guys: 0th) mass conservation equation static constexpr int pressureIdx = massBalanceIdx; //!< Index of the pressure in a solution vector - static constexpr int faceOffset = GET_PROP_VALUE(TypeTag, NumEqCellCenter); //!< Index of the momentum balance equation - static constexpr int momentumBalanceIdx = PVOffset + faceOffset; //!< Index of the momentum balance equation + using GridView = typename GET_PROP_TYPE(TypeTag, GridView); + static constexpr auto dim = GridView::dimension; + static constexpr auto numEq = GET_PROP_VALUE(TypeTag, NumEq); + static constexpr auto momentumBalanceOffset = GET_PROP_VALUE(TypeTag, NumEq) - dim; + + static constexpr int momentumBalanceIdx = PVOffset + momentumBalanceOffset; //!< Index of the momentum balance equation static constexpr int momentumXBalanceIdx = momentumBalanceIdx; //!< Index of the momentum balance equation static constexpr int momentumYBalanceIdx = momentumBalanceIdx + 1; //!< Index of the momentum balance equation static constexpr int momentumZBalanceIdx = momentumBalanceIdx + 2; //!< Index of the momentum balance equation - static constexpr int velocityIdx = momentumBalanceIdx; //!< Index of the velocity in a solution vector (NOTE: This PV lives in a different vector than the pressure) - static constexpr int velocityXIdx = velocityIdx; //!< Index of the velocity in a solution vector (NOTE: This PV lives in a different vector than the pressure) - static constexpr int velocityYIdx = velocityIdx + 1; //!< Index of the velocity in a solution vector (NOTE: This PV lives in a different vector than the pressure) - static constexpr int velocityZIdx = velocityIdx + 2; //!< Index of the velocity in a solution vector (NOTE: This PV lives in a different vector than the pressure) + static constexpr int velocityIdx = momentumBalanceIdx; //!< Index of the velocity in a solution vector + static constexpr int velocityXIdx = velocityIdx; //!< Index of the velocity in a solution vector + static constexpr int velocityYIdx = velocityIdx + 1; //!< Index of the velocity in a solution vector + static constexpr int velocityZIdx = velocityIdx + 2; //!< Index of the velocity in a solution vector - static constexpr int phaseIdx = 0; //!< Index of the fluid phase (required to use the same fluid system in coupled models) }; // \} diff --git a/dumux/implicit/staggered/primaryvariables.hh b/dumux/implicit/staggered/primaryvariables.hh index 06c1b66ecb..64114a8387 100644 --- a/dumux/implicit/staggered/primaryvariables.hh +++ b/dumux/implicit/staggered/primaryvariables.hh @@ -46,6 +46,8 @@ class StaggeredPrimaryVariables : public Dune::MultiTypeBlockVector<CellCenterPr using ParentType = Dune::MultiTypeBlockVector<CellCenterPrimaryVariables, FacePrimaryVariables>; + static constexpr auto faceOffset = GET_PROP_VALUE(TypeTag, NumEqCellCenter); + public: StaggeredPrimaryVariables() = default; @@ -87,10 +89,10 @@ public: using ParentType::operator []; const Scalar& operator [](const unsigned int pvIdx) const { - if(pvIdx < Indices::faceOffset) + if(pvIdx < faceOffset) return ParentType::operator[](cellCenterIdx)[pvIdx]; else - return ParentType::operator[](faceIdx)[pvIdx - Indices::faceOffset]; + return ParentType::operator[](faceIdx)[pvIdx - faceOffset]; } /*! @@ -101,10 +103,10 @@ public: */ Scalar& operator [](const unsigned int pvIdx) { - if(pvIdx < Indices::faceOffset) + if(pvIdx < faceOffset) return ParentType::operator[](cellCenterIdx)[pvIdx]; else - return ParentType::operator[](faceIdx)[pvIdx - Indices::faceOffset]; + return ParentType::operator[](faceIdx)[pvIdx - faceOffset]; } }; -- GitLab