From ca320997e050cdcd6718e5c8c41de62d76abda71 Mon Sep 17 00:00:00 2001 From: Bernd Flemisch <bernd@iws.uni-stuttgart.de> Date: Wed, 14 Oct 2015 13:12:09 +0200 Subject: [PATCH] [range-for] adapt the test folder Some tests still use the GridView's begin method explicitly in order to obtain the first element with which some calculations are performed. It should be evaluated if this can be circumvented. --- test/decoupled/1p/resultevaluation.hh | 20 ++--- test/decoupled/1p/resultevaluation3d.hh | 73 +++++++------------ test/decoupled/1p/test_1pspatialparams.hh | 7 +- test/decoupled/1p/test_diffusionproblem.hh | 7 +- test/decoupled/1p/test_diffusionproblem3d.hh | 7 +- .../1p/test_diffusionspatialparams.hh | 13 +--- .../2p/buckleyleverettanalyticsolution.hh | 27 +++---- .../decoupled/2p/mcwhorteranalyticsolution.hh | 31 ++++---- test/decoupled/2p/test_3d2pproblem.hh | 11 +-- test/decoupled/2p/test_transportproblem.hh | 14 ++-- test/geomechanics/el2p/el2pproblem.hh | 33 +++------ test/implicit/1p/1pniconductionproblem.hh | 26 +++---- test/implicit/1p/1pniconvectionproblem.hh | 24 +++--- test/implicit/1p2c/1p2cniconductionproblem.hh | 26 +++---- test/implicit/1p2c/1p2cniconvectionproblem.hh | 24 +++--- test/implicit/2p/cc2pcornerpointproblem.hh | 11 +-- test/implicit/2pnc/fuelcellproblem.hh | 10 +-- test/implicit/3p/3pniconductionproblem.hh | 26 +++---- test/implicit/3p/3pniconvectionproblem.hh | 24 +++--- test/implicit/3p3c/columnxylolproblem.hh | 11 +-- test/implicit/3p3c/infiltrationproblem.hh | 11 +-- test/implicit/3p3c/kuevetteproblem.hh | 11 +-- test/implicit/co2/heterogeneousproblem.hh | 17 ++--- test/implicit/co2/heterogeneousproblemni.hh | 17 ++--- .../co2/heterogeneousspatialparameters.hh | 9 +-- .../richards/richardsanalyticalproblem.hh | 7 +- .../richards/richardsniconductionproblem.hh | 26 +++---- .../richards/richardsniconvectionproblem.hh | 24 +++--- 28 files changed, 221 insertions(+), 326 deletions(-) diff --git a/test/decoupled/1p/resultevaluation.hh b/test/decoupled/1p/resultevaluation.hh index 60298b128b..aff224927b 100644 --- a/test/decoupled/1p/resultevaluation.hh +++ b/test/decoupled/1p/resultevaluation.hh @@ -76,8 +76,6 @@ public: enum {dim=Grid::dimension}; typedef typename Grid::template Codim<0>::Entity Element; typedef typename Element::Geometry Geometry; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; - typedef typename GridView::IntersectionIterator IntersectionIterator; typedef typename Geometry::JacobianInverseTransposed JacobianInverseTransposed; uMin = 1e100; @@ -100,11 +98,8 @@ public: double denominatorGrad = 0; double numeratorFlux = 0; double denominatorFlux = 0; - ElementIterator endEIt = gridView.template end<0>(); - for (ElementIterator eIt = gridView.template begin<0>(); eIt != endEIt; ++eIt) + for (const auto& element : Dune::elements(gridView)) { - const Element& element = *eIt; - // element geometry const Geometry& geometry = element.geometry(); Dune::GeometryType geomType = geometry.type(); @@ -136,22 +131,21 @@ public: int isIdx = -1; Dune::FieldVector<Scalar,2*dim> fluxVector; Dune::FieldVector<Scalar,dim> exactGradient; - IntersectionIterator endis = gridView.iend(element); - for (IntersectionIterator is = gridView.ibegin(element); is!=endis; ++is) + for (const auto& intersection : Dune::intersections(gridView, element)) { // local number of facet - int fIdx = is->indexInInside(); + int fIdx = intersection.indexInInside(); if (consecutiveNumbering) isIdx++; else isIdx = fIdx; - Dune::FieldVector<double,dim> faceGlobal = is->geometry().center(); - double faceVol = is->geometry().volume(); + Dune::FieldVector<double,dim> faceGlobal = intersection.geometry().center(); + double faceVol = intersection.geometry().volume(); // get normal vector - Dune::FieldVector<double,dim> unitOuterNormal = is->centerUnitOuterNormal(); + Dune::FieldVector<double,dim> unitOuterNormal = intersection.centerUnitOuterNormal(); // get the exact gradient exactGradient = problem.exactGrad(faceGlobal); @@ -180,7 +174,7 @@ public: approximateFlux *= faceVol; fluxVector[fIdx] = approximateFlux; - if (!is->neighbor()) { + if (!intersection.neighbor()) { if (fabs(faceGlobal[1]) < 1e-6) { fluy0 += approximateFlux; exactfluy0 += exactFlux; diff --git a/test/decoupled/1p/resultevaluation3d.hh b/test/decoupled/1p/resultevaluation3d.hh index 2a98b6a866..e80d92d3ca 100644 --- a/test/decoupled/1p/resultevaluation3d.hh +++ b/test/decoupled/1p/resultevaluation3d.hh @@ -86,8 +86,6 @@ public: typedef typename Entity::Geometry Geometry; typedef typename Grid::LevelGridView GV; typedef typename GV::IndexSet IS; - typedef typename GV::template Codim<0>::Iterator Iterator; - typedef typename GV::IntersectionIterator IntersectionIterator; typedef Dune::MultipleCodimMultipleGeomTypeMapper<GV,Dune::MCMGElementLayout> EM; typedef Dune::MultipleCodimMultipleGeomTypeMapper<GV,FaceLayout> FM; typedef typename Grid::ctype ct; @@ -104,12 +102,8 @@ public: uMean = 0; double domainVolume = 0; - Iterator eendit = gridview.template end<0>(); - for (Iterator it = gridview.template begin<0>(); it != eendit; ++it) + for (const auto& element : Dune::elements(gridview)) { - // get entity - const Entity& element = *it; - // get volume double volume = element.geometry().volume(); @@ -161,11 +155,8 @@ public: double numeratorFlux = 0; double denominatorFlux = 0; bool exactOutput = true; - for (Iterator it = gridview.template begin<0>(); it != eendit; ++it) + for (const auto& element : Dune::elements(gridview)) { - // get entity - const Entity& element = *it; - // element geometry const Geometry& geometry = element.geometry(); @@ -214,20 +205,19 @@ public: int i = -1; Dune::FieldVector<ct,2*dim> fluxVector; Dune::FieldVector<ct,dim> exactGradient; - IntersectionIterator endis = gridview.iend(element); - for (IntersectionIterator is = gridview.ibegin(element); is!=endis; ++is) + for (const auto& intersection : Dune::intersections(gridview, element)) { // local number of facet - i = is->indexInInside(); + i = intersection.indexInInside(); // global number of face int faceIndex = facemapper.template map<1>(element, i); - Dune::FieldVector<double,dim> faceGlobal = is->geometry().center(); - double faceVol = is->geometry().volume(); + Dune::FieldVector<double,dim> faceGlobal = intersection.geometry().center(); + double faceVol = intersection.geometry().volume(); // get normal vector - Dune::FieldVector<double,dim> unitOuterNormal = is->centerUnitOuterNormal(); + Dune::FieldVector<double,dim> unitOuterNormal = intersection.centerUnitOuterNormal(); // get the approximate solution on the face double approximateFace = (*solution.pressTrace)[faceIndex]; @@ -260,7 +250,7 @@ public: fluxVector[i] = approximateFlux; //if (is.boundary()) { - if (!is->neighbor()) { + if (!intersection.neighbor()) { if (fabs(faceGlobal[2]) < 1e-6) { fluz0 += approximateFlux; exactfluz0 += exactFlux; @@ -427,8 +417,6 @@ public: enum {dim=Grid::dimension}; typedef typename Grid::template Codim<0>::Entity Element; typedef typename Element::Geometry Geometry; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; - typedef typename GridView::IntersectionIterator IntersectionIterator; typedef Dune::MultipleCodimMultipleGeomTypeMapper<GridView, Dune::MCMGElementLayout> ElementMapper; typedef Dune::BlockVector<Dune::FieldVector<Scalar, 1> > SolVector; typedef typename Geometry::JacobianInverseTransposed JacobianInverseTransposed; @@ -472,11 +460,8 @@ public: double numeratorFluxIn = 0; double denominatorFluxIn = 0; bool exactOutput = true; - ElementIterator endEIt = gridView.template end<0>(); - for (ElementIterator eIt = gridView.template begin<0>(); eIt != endEIt; ++eIt) + for (const auto& element : Dune::elements(gridView)) { - const Element& element = *eIt; - // element geometry const Geometry& geometry = element.geometry(); @@ -505,11 +490,9 @@ public: // calculate the relative error for inner cells bool bndCell = false; - IntersectionIterator beginis = gridView.ibegin(element); - IntersectionIterator endis = gridView.iend(element); - for (IntersectionIterator is = beginis; is!=endis; ++is) + for (const auto& intersection : Dune::intersections(gridView, element)) { - if (is->boundary()) + if (intersection.boundary()) { bndCell = true; break; @@ -545,21 +528,21 @@ public: int isIdx = -1; Dune::FieldVector<Scalar,2*dim> fluxVector; Dune::FieldVector<Scalar,dim> exactGradient; - for (IntersectionIterator is = beginis; is!=endis; ++is) + for (const auto& intersection : Dune::intersections(gridView, element)) { // local number of facet - int fIdx = is->indexInInside(); + int fIdx = intersection.indexInInside(); if (consecutiveNumbering) isIdx++; else isIdx = fIdx; - Dune::FieldVector<double,dim> faceGlobal = is->geometry().center(); - double faceVol = is->geometry().volume(); + Dune::FieldVector<double,dim> faceGlobal = intersection.geometry().center(); + double faceVol = intersection.geometry().volume(); // get normal vector - Dune::FieldVector<double,dim> unitOuterNormal = is->centerUnitOuterNormal(); + Dune::FieldVector<double,dim> unitOuterNormal = intersection.centerUnitOuterNormal(); // get the exact gradient exactGradient = problem.exactGrad(faceGlobal); @@ -594,7 +577,7 @@ public: approximateFlux *= faceVol; fluxVector[fIdx] = approximateFlux; - if (!is->neighbor()) { + if (!intersection.neighbor()) { if (fabs(faceGlobal[2]) < 1e-6) { fluz0 += approximateFlux; exactfluz0 += exactFlux; @@ -661,14 +644,14 @@ public: denominatorGrad += volume*(exactGradient*exactGradient); // calculate the maximum of the diagonal length of all elements on leaf grid - for (int i = 0; i < eIt->geometry().corners(); ++i) + for (int i = 0; i < element.geometry().corners(); ++i) { - Dune::FieldVector<Scalar,dim> corner1 = eIt->geometry().corner(i); + Dune::FieldVector<Scalar,dim> corner1 = element.geometry().corner(i); - for (int j = 0; j < eIt->geometry().corners(); ++j) + for (int j = 0; j < element.geometry().corners(); ++j) { // get all corners of current element and compare the distances between them - Dune::FieldVector<Scalar,dim> corner2 = eIt->geometry().corner(j); + Dune::FieldVector<Scalar,dim> corner2 = element.geometry().corner(j); // distance vector between corners Dune::FieldVector<Scalar,dim> distVec = corner1 - corner2; @@ -733,11 +716,8 @@ public: Scalar maxFlux = -1; Scalar numerator = 0; Scalar denominator = 0; - ElementIterator endEIt = gridView.template end<0>(); - for (ElementIterator eIt = gridView.template begin<0>(); eIt != endEIt; ++eIt) + for (const auto& element : Dune::elements(gridView)) { - const Element& element = *eIt; - // element geometry const Geometry& geometry = element.geometry(); @@ -772,11 +752,10 @@ public: int i = -1; Dune::FieldVector<Scalar,dim> exactGradient; - IntersectionIterator endis = gridView.iend(element); - for (IntersectionIterator is = gridView.ibegin(element); is!=endis; ++is) + for (const auto& intersection : Dune::intersections(gridView, element)) { // get geometry type of face - Dune::GeometryType gtf = is->geometryInInside().type(); + Dune::GeometryType gtf = intersection.geometryInInside().type(); // local number of facet i++; @@ -785,10 +764,10 @@ public: const Dune::FieldVector<Scalar,dim-1>& faceLocalNm1 = ReferenceFaces::general(gtf).position(0,0); // center of face in global coordinates - Dune::FieldVector<Scalar,dim> faceGlobal = is->geometry().global(faceLocalNm1); + Dune::FieldVector<Scalar,dim> faceGlobal = intersection.geometry().global(faceLocalNm1); // get normal vector - Dune::FieldVector<Scalar,dim> unitOuterNormal = is->unitOuterNormal(faceLocalNm1); + Dune::FieldVector<Scalar,dim> unitOuterNormal = intersection.unitOuterNormal(faceLocalNm1); if (switchNormals) unitOuterNormal *= -1.0; diff --git a/test/decoupled/1p/test_1pspatialparams.hh b/test/decoupled/1p/test_1pspatialparams.hh index aa30c43948..e8c93b060a 100644 --- a/test/decoupled/1p/test_1pspatialparams.hh +++ b/test/decoupled/1p/test_1pspatialparams.hh @@ -46,7 +46,6 @@ class TestOnePSpatialParams: public FVSpatialParamsOneP<TypeTag> enum {dim=Grid::dimension, dimWorld=Grid::dimensionworld}; typedef typename Grid::Traits::template Codim<0>::Entity Element; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; typedef Dune::FieldVector<CoordScalar, dimWorld> GlobalPosition; typedef Dune::FieldMatrix<Scalar,dim,dim> FieldMatrix; @@ -68,11 +67,9 @@ public: delta_ = delta; permeability_.resize(gridView_.size(0)); - ElementIterator eIt = gridView_.template begin<0>(); - ElementIterator eEndIt = gridView_.template end<0>(); - for(;eIt != eEndIt; ++eIt) + for(const auto& element : Dune::elements(gridView_)) { - setPermeability_(permeability_[indexSet_.index(*eIt)], eIt->geometry().center()); + setPermeability_(permeability_[indexSet_.index(element)], element.geometry().center()); } } diff --git a/test/decoupled/1p/test_diffusionproblem.hh b/test/decoupled/1p/test_diffusionproblem.hh index 2d5d864a5c..5079079225 100644 --- a/test/decoupled/1p/test_diffusionproblem.hh +++ b/test/decoupled/1p/test_diffusionproblem.hh @@ -173,7 +173,6 @@ class TestDiffusionProblem: public DiffusionProblem2P<TypeTag> typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; typedef typename GridView::Traits::template Codim<0>::Entity Element; typedef Dune::FieldVector<Scalar, dimWorld> GlobalPosition; typedef Dune::FieldVector<Scalar, dim> LocalPosition; @@ -220,11 +219,9 @@ public: { ScalarSolution *exactPressure = this->resultWriter().allocateManagedBuffer(this->gridView().size(0)); - ElementIterator eIt = this->gridView().template begin<0>(); - ElementIterator eEndIt = this->gridView().template end<0>(); - for(;eIt != eEndIt; ++eIt) + for(const auto& element : Dune::elements(this->gridView())) { - (*exactPressure)[this->elementMapper().index(*eIt)][0] = exact(eIt->geometry().center()); + (*exactPressure)[this->elementMapper().index(element)][0] = exact(element.geometry().center()); } this->resultWriter().attachCellData(*exactPressure, "exact pressure"); diff --git a/test/decoupled/1p/test_diffusionproblem3d.hh b/test/decoupled/1p/test_diffusionproblem3d.hh index c92fcec158..8cee46df6f 100644 --- a/test/decoupled/1p/test_diffusionproblem3d.hh +++ b/test/decoupled/1p/test_diffusionproblem3d.hh @@ -137,7 +137,6 @@ class TestDiffusion3DProblem: public DiffusionProblem2P<TypeTag> typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; typedef typename GridView::Traits::template Codim<0>::Entity Element; typedef typename GridView::Intersection Intersection; typedef Dune::FieldVector<Scalar, dimWorld> GlobalPosition; @@ -174,11 +173,9 @@ public: { ScalarSolution *exactPressure = this->resultWriter().allocateManagedBuffer(this->gridView().size(0)); - ElementIterator eIt = this->gridView().template begin<0>(); - ElementIterator eEndIt = this->gridView().template end<0>(); - for(;eIt != eEndIt; ++eIt) + for(const auto& element : Dune::elements(this->gridView())) { - (*exactPressure)[this->elementMapper().index(*eIt)][0] = exact(eIt->geometry().center()); + (*exactPressure)[this->elementMapper().index(element)][0] = exact(element.geometry().center()); } this->resultWriter().attachCellData(*exactPressure, "exact pressure"); diff --git a/test/decoupled/1p/test_diffusionspatialparams.hh b/test/decoupled/1p/test_diffusionspatialparams.hh index 31e5f9031c..b021516826 100644 --- a/test/decoupled/1p/test_diffusionspatialparams.hh +++ b/test/decoupled/1p/test_diffusionspatialparams.hh @@ -74,7 +74,6 @@ class TestDiffusionSpatialParams: public FVSpatialParams<TypeTag> enum {dim=Grid::dimension, dimWorld=Grid::dimensionworld}; typedef typename Grid::Traits::template Codim<0>::Entity Element; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; typedef Dune::FieldVector<CoordScalar, dimWorld> GlobalPosition; typedef Dune::FieldMatrix<Scalar,dim,dim> FieldMatrix; @@ -105,11 +104,9 @@ public: delta_ = delta; permeability_.resize(gridView_.size(0)); - ElementIterator eIt = gridView_.template begin<0>(); - ElementIterator eEndIt = gridView_.template end<0>(); - for(;eIt != eEndIt; ++eIt) + for(const auto& element : Dune::elements(gridView_)) { - perm(permeability_[indexSet_.index(*eIt)], eIt->geometry().center()); + perm(permeability_[indexSet_.index(element)], element.geometry().center()); } } @@ -121,11 +118,9 @@ public: ScalarSolution *permXY = writer.allocateManagedBuffer(gridView_.size(0)); ScalarSolution *permYY = writer.allocateManagedBuffer(gridView_.size(0)); - ElementIterator eIt = gridView_.template begin<0>(); - ElementIterator eEndIt = gridView_.template end<0>(); - for(;eIt != eEndIt; ++eIt) + for(const auto& element : Dune::elements(gridView_)) { - int eIdxGlobal = indexSet_.index(*eIt); + int eIdxGlobal = indexSet_.index(element); (*permXX)[eIdxGlobal][0] = permeability_[eIdxGlobal][0][0]; (*permXY)[eIdxGlobal][0] = permeability_[eIdxGlobal][0][1]; (*permYY)[eIdxGlobal][0] = permeability_[eIdxGlobal][1][1]; diff --git a/test/decoupled/2p/buckleyleverettanalyticsolution.hh b/test/decoupled/2p/buckleyleverettanalyticsolution.hh index 5ec62412f9..e65b4f5aca 100644 --- a/test/decoupled/2p/buckleyleverettanalyticsolution.hh +++ b/test/decoupled/2p/buckleyleverettanalyticsolution.hh @@ -88,7 +88,6 @@ template<class TypeTag> class BuckleyLeverettAnalytic }; typedef Dune::BlockVector<Dune::FieldVector<Scalar, 1> > BlockVector; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; typedef Dune::FieldVector<Scalar, dimworld> GlobalPosition; private: @@ -113,17 +112,17 @@ private: */ void prepareAnalytic() { - ElementIterator dummyElement = problem_.gridView().template begin<0>(); - const MaterialLawParams& materialLawParams(problem_.spatialParams().materialLawParams(*dummyElement)); + const auto& dummyElement = *problem_.gridView().template begin<0>(); + const MaterialLawParams& materialLawParams(problem_.spatialParams().materialLawParams(dummyElement)); swr_ = materialLawParams.swr(); snr_ = materialLawParams.snr(); - Scalar porosity = problem_.spatialParams().porosity(*dummyElement); + Scalar porosity = problem_.spatialParams().porosity(dummyElement); FluidState fluidState; - fluidState.setTemperature(problem_.temperature(*dummyElement)); - fluidState.setPressure(wPhaseIdx, problem_.referencePressure(*dummyElement)); - fluidState.setPressure(nPhaseIdx, problem_.referencePressure(*dummyElement)); + fluidState.setTemperature(problem_.temperature(dummyElement)); + fluidState.setPressure(wPhaseIdx, problem_.referencePressure(dummyElement)); + fluidState.setPressure(nPhaseIdx, problem_.referencePressure(dummyElement)); Scalar viscosityW = FluidSystem::viscosity(fluidState, wPhaseIdx); Scalar viscosityNW = FluidSystem::viscosity(fluidState, nPhaseIdx); @@ -196,15 +195,14 @@ private: Scalar globalVolume = 0; Scalar errorNorm = 0.0; - ElementIterator eEndIt = problem_.gridView().template end<0> (); - for (ElementIterator eIt = problem_.gridView().template begin<0> (); eIt != eEndIt; ++eIt) + for (const auto& element : Dune::elements(problem_.gridView())) { // get entity - int index = problem_.variables().index(*eIt); + int index = problem_.variables().index(element); Scalar sat = problem_.variables().cellData(index).saturation(wPhaseIdx); - Scalar volume = eIt->geometry().volume(); + Scalar volume = element.geometry().volume(); Scalar error = analyticSolution_[index] - sat; @@ -239,13 +237,12 @@ private: Scalar xMax = frontParams_[0].second * time; // iterate over vertices and get analytic saturation solution - ElementIterator eEndIt = problem_.gridView().template end<0> (); - for (ElementIterator eIt = problem_.gridView().template begin<0> (); eIt != eEndIt; ++eIt) + for (const auto& element : Dune::elements(problem_.gridView())) { // get global coordinate of cell center - GlobalPosition globalPos = eIt->geometry().center(); + GlobalPosition globalPos = element.geometry().center(); - int index = problem_.variables().index(*eIt); + int index = problem_.variables().index(element); if (globalPos[0] > xMax) analyticSolution_[index] = swr_; diff --git a/test/decoupled/2p/mcwhorteranalyticsolution.hh b/test/decoupled/2p/mcwhorteranalyticsolution.hh index c82ca7984b..6a5523bcc3 100644 --- a/test/decoupled/2p/mcwhorteranalyticsolution.hh +++ b/test/decoupled/2p/mcwhorteranalyticsolution.hh @@ -70,7 +70,6 @@ class McWhorterAnalytic }; typedef Dune::BlockVector<Dune::FieldVector<Scalar, 1> > BlockVector; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; typedef Dune::FieldVector<Scalar, dimworld> GlobalPosition; private: @@ -95,15 +94,14 @@ private: Scalar globalVolume = 0; Scalar errorNorm = 0.0; - ElementIterator eEndIt = problem_.gridView().template end<0> (); - for (ElementIterator eIt = problem_.gridView().template begin<0> (); eIt != eEndIt; ++eIt) + for (const auto& element : Dune::elements(problem_.gridView())) { // get entity - int index = problem_.variables().index(*eIt); + int index = problem_.variables().index(element); Scalar sat = problem_.variables().cellData(index).saturation(wPhaseIdx); - Scalar volume = eIt->geometry().volume(); + Scalar volume = element.geometry().volume(); Scalar error = analyticSolution_[index] - sat; @@ -131,15 +129,15 @@ private: void prepareAnalytic() { - ElementIterator dummyElement = problem_.gridView().template begin<0>(); - const MaterialLawParams& materialLawParams(problem_.spatialParams().materialLawParams(*dummyElement)); + const auto& dummyElement = *problem_.gridView().template begin<0>(); + const MaterialLawParams& materialLawParams(problem_.spatialParams().materialLawParams(dummyElement)); swr_ = materialLawParams.swr(); snr_ = materialLawParams.snr(); - porosity_ = problem_.spatialParams().porosity(*dummyElement); - permeability_ = problem_.spatialParams().intrinsicPermeability(*dummyElement)[0][0]; + porosity_ = problem_.spatialParams().porosity(dummyElement); + permeability_ = problem_.spatialParams().intrinsicPermeability(dummyElement)[0][0]; PrimaryVariables initVec; - problem_.initial(initVec, *dummyElement); + problem_.initial(initVec, dummyElement); sInit_ = initVec[saturationIdx]; Scalar s0 =(1 - snr_ - swr_); // std::cerr<<"\n swr, snr: "<< swr_ << snr_ << "\n"; @@ -156,9 +154,9 @@ private: satVec_[i]=satVec_[i-1]+h_; } FluidState fluidState; - fluidState.setTemperature(problem_.temperature(*dummyElement)); - fluidState.setPressure(wPhaseIdx, problem_.referencePressure(*dummyElement)); - fluidState.setPressure(nPhaseIdx, problem_.referencePressure(*dummyElement)); + fluidState.setTemperature(problem_.temperature(dummyElement)); + fluidState.setPressure(wPhaseIdx, problem_.referencePressure(dummyElement)); + fluidState.setPressure(nPhaseIdx, problem_.referencePressure(dummyElement)); Scalar viscosityW = FluidSystem::viscosity(fluidState, wPhaseIdx); Scalar viscosityNW = FluidSystem::viscosity(fluidState, nPhaseIdx); @@ -301,16 +299,15 @@ private: // std::cout<<" xf_ = "<<xf_<<std::endl; // iterate over vertices and get analytical saturation solution - ElementIterator eEndIt = problem_.gridView().template end<0> (); - for (ElementIterator eIt = problem_.gridView().template begin<0> (); eIt!= eEndIt; ++eIt) + for (const auto& element : Dune::elements(problem_.gridView())) { // get global coordinate of cell center - const GlobalPosition& globalPos = eIt->geometry().center(); + const GlobalPosition& globalPos = element.geometry().center(); // std::cout<<"globalPos = "<<globalPos[0]<<", x0 = "<<xf_[0]<<"\n"; // find index of current vertex - int index = problem_.variables().index(*eIt); + int index = problem_.variables().index(element); // find x_f next to global coordinate of the vertex int xnext = 0; diff --git a/test/decoupled/2p/test_3d2pproblem.hh b/test/decoupled/2p/test_3d2pproblem.hh index dfd8a692e2..b27bcac32b 100644 --- a/test/decoupled/2p/test_3d2pproblem.hh +++ b/test/decoupled/2p/test_3d2pproblem.hh @@ -174,7 +174,6 @@ typedef typename GridView::Traits::template Codim<0>::Entity Element; typedef typename GridView::Intersection Intersection; typedef Dune::FieldVector<Scalar, dimWorld> GlobalPosition; typedef Dune::FieldVector<Scalar, dim> LocalPosition; -typedef typename GridView::template Codim<dim>::Iterator VertexIterator; typedef typename GET_PROP_TYPE(TypeTag, GridCreator) GridCreator; typedef typename GET_PROP(TypeTag, ParameterTree) ParameterTree; @@ -197,21 +196,19 @@ ParentType(timeManager, gridView), inflowEdge_(0), outflowEdge_(0) Scalar maxDist = this->bBoxMin().two_norm(); // calculate the bounding box of the grid view - VertexIterator vIt = gridView.template begin<dim>(); - const VertexIterator vEndIt = gridView.template end<dim>(); - for (; vIt!=vEndIt; ++vIt) { - GlobalPosition vertexCoord(vIt->geometry().center()); + for (const auto& vertex : Dune::vertices(gridView)) { + GlobalPosition vertexCoord(vertex.geometry().center()); Scalar dist = vertexCoord.two_norm(); if (dist > maxDist) { maxDist = dist; - outflowEdge_ = vIt->geometry().center(); + outflowEdge_ = vertex.geometry().center(); } if (dist < minDist) { minDist = dist; - inflowEdge_ = vIt->geometry().center(); + inflowEdge_ = vertex.geometry().center(); } } diff --git a/test/decoupled/2p/test_transportproblem.hh b/test/decoupled/2p/test_transportproblem.hh index 0358af2e44..84b4a9e7af 100644 --- a/test/decoupled/2p/test_transportproblem.hh +++ b/test/decoupled/2p/test_transportproblem.hh @@ -106,8 +106,6 @@ class TestTransportProblem: public TransportProblem2P<TypeTag> typedef typename GET_PROP(TypeTag, SolutionTypes) SolutionTypes; typedef typename SolutionTypes::PrimaryVariables PrimaryVariables; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; - typedef typename GridView::IntersectionIterator IntersectionIterator; typedef typename GET_PROP_TYPE(TypeTag, CellData) CellData; enum @@ -139,24 +137,22 @@ public: vel[0] = 1e-5; // compute update vector - ElementIterator eEndIt = this->gridView().template end<0> (); - for (ElementIterator eIt = this->gridView().template begin<0> (); eIt != eEndIt; ++eIt) + for (const auto& element : Dune::elements(this->gridView())) { // cell index - int eIdxGlobal = this->elementMapper().index(*eIt); + int eIdxGlobal = this->elementMapper().index(element); CellData& cellData = this->variables().cellData(eIdxGlobal); // run through all intersections with neighbors and boundary - IntersectionIterator isEndIt = this->gridView().iend(*eIt); - for (IntersectionIterator isIt = this->gridView().ibegin(*eIt); isIt != isEndIt; ++isIt) + for (const auto& intersection : Dune::intersections(this->gridView(), element)) { // local number of facet - int indexInInside = isIt->indexInInside(); + int indexInInside = intersection.indexInInside(); cellData.fluxData().setVelocity(wPhaseIdx, indexInInside, vel); - const GlobalPosition& unitOuterNormal = isIt->centerUnitOuterNormal(); + const GlobalPosition& unitOuterNormal = intersection.centerUnitOuterNormal(); Scalar pot = vel * unitOuterNormal; diff --git a/test/geomechanics/el2p/el2pproblem.hh b/test/geomechanics/el2p/el2pproblem.hh index 6d914cf5b4..ba961971a8 100644 --- a/test/geomechanics/el2p/el2pproblem.hh +++ b/test/geomechanics/el2p/el2pproblem.hh @@ -186,7 +186,6 @@ class El2P_TestProblem : public ImplicitPorousMediaProblem<TypeTag> typedef typename GridView::template Codim<0>::Entity Element; typedef typename GridView::template Codim<dim>::Entity Vertex; typedef typename GridView::Intersection Intersection; - typedef typename GridView::template Codim<dim>::Iterator VertexIterator; typedef typename GET_PROP_TYPE(TypeTag, VertexMapper) VertexMapper; typedef typename GET_PROP_TYPE(TypeTag, FVElementGeometry) FVElementGeometry; @@ -275,12 +274,10 @@ public: // initialization of the actual model run and for the pressure Dirichlet boundary values. void initializePressure() { - VertexIterator vIt = gridView_.template begin<dim>(); - VertexIterator vEndIt = gridView_.template end<dim>(); - for(; vIt != vEndIt; ++vIt) + for(const auto& vertex : Dune::vertices(gridView_)) { - int vIdxGlobal = this->vertexMapper().index(*vIt); - GlobalPosition globalPos = (*vIt).geometry().corner(0); + int vIdxGlobal = this->vertexMapper().index(vertex); + GlobalPosition globalPos = vertex.geometry().corner(0); // initial approximate pressure distribution at start of initialization run pInit_[vIdxGlobal] = -(1.013e5 + (depthBOR_ - globalPos[2]) * brineDensity_ * 9.81); @@ -317,11 +314,9 @@ public: this->setInitializationRun(initializationRun_); std::cout<<"El2P_TestProblem: initialized pressure field copied to pInit_"<<std::endl; - VertexIterator vIt = gridView_.template begin<dim>(); - VertexIterator vEndIt = gridView_.template end<dim>(); - for(; vIt != vEndIt; ++vIt) + for(const auto& vertex : Dune::vertices(gridView_)) { - int vIdxGlobal = this->vertexMapper().index(*vIt); + int vIdxGlobal = this->vertexMapper().index(vertex); pInit_[vIdxGlobal] = -this->model().curSol().base()[vIdxGlobal*2][0]; } } @@ -837,19 +832,15 @@ public: inline void evaluateGlobal(const DomainType & position, RangeType & values) const { bool valueSet; - VertexIterator vIt = - gridView_.template begin<GridView::dimension> (); - VertexIterator vEndIt = - gridView_.template end<GridView::dimension> (); valueSet = false; // loop over all vertices - for (; vIt != vEndIt; ++vIt) + for (const auto& vertex : Dune::vertices(gridView_)) { // get global index of current vertex - int vIdxGlobal = vertexMapper_.index(*vIt); + int vIdxGlobal = vertexMapper_.index(vertex); Dune::FieldVector<double, 3> globalPos = - (*vIt).geometry().corner(0); + (vertex).geometry().corner(0); // compare coordinates of current vertex with position coordinates if (globalPos[0] >= position[0] - eps_ && globalPos[0] <= position[0] + eps_ @@ -886,13 +877,9 @@ public: void setPressure(std::vector<Scalar> pInit) { std::cout << "InitialPressSat: setPressure function called" << std::endl; - VertexIterator vIt = - gridView_.template begin<GridView::dimension> (); - VertexIterator vEndIt = - gridView_.template end<GridView::dimension> (); - for (; vIt != vEndIt; ++vIt) + for (const auto& vertex : Dune::vertices(gridView_)) { - int vIdxGlobal = vertexMapper_.index(*vIt); + int vIdxGlobal = vertexMapper_.index(vertex); pInit_[vIdxGlobal] = -pInit[vIdxGlobal]; } } diff --git a/test/implicit/1p/1pniconductionproblem.hh b/test/implicit/1p/1pniconductionproblem.hh index 0743431ae0..ec516b544e 100644 --- a/test/implicit/1p/1pniconductionproblem.hh +++ b/test/implicit/1p/1pniconductionproblem.hh @@ -127,7 +127,6 @@ class OnePNIConductionProblem : public ImplicitPorousMediaProblem<TypeTag> typedef typename GridView::template Codim<0>::Entity Element; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; typedef typename GridView::Intersection Intersection; typedef Dune::FieldVector<Scalar, dimWorld> GlobalPosition; @@ -177,9 +176,8 @@ public: FVElementGeometry fvGeometry; VolumeVariables volVars; - ElementIterator eIt = this->gridView().template begin<0>(); - ElementIterator eEndIt = this->gridView().template end<0>(); - fvGeometry.update(this->gridView(), *eIt); + const auto& element = *this->gridView().template begin<0>(); + fvGeometry.update(this->gridView(), element); PrimaryVariables initialPriVars(0); GlobalPosition globalPos(0); initial_(initialPriVars, globalPos); @@ -187,33 +185,33 @@ public: //update the constant volume variables volVars.update(initialPriVars, *this, - *eIt, + element, fvGeometry, 0, false); - Scalar porosity = this->spatialParams().porosity(*eIt, fvGeometry, 0); + Scalar porosity = this->spatialParams().porosity(element, fvGeometry, 0); Scalar densityW = volVars.density(); Scalar heatCapacityW = FluidSystem::heatCapacity(volVars.fluidState(), 0); - Scalar densityS = this->spatialParams().solidDensity(*eIt, fvGeometry, 0); - Scalar heatCapacityS = this->spatialParams().solidHeatCapacity(*eIt, fvGeometry, 0); + Scalar densityS = this->spatialParams().solidDensity(element, fvGeometry, 0); + Scalar heatCapacityS = this->spatialParams().solidHeatCapacity(element, fvGeometry, 0); Scalar storage = densityW*heatCapacityW*porosity + densityS*heatCapacityS*(1 - porosity); Scalar effectiveThermalConductivity = ThermalConductivityModel::effectiveThermalConductivity(volVars, this->spatialParams(), - *eIt, fvGeometry, 0); + element, fvGeometry, 0); Scalar time = std::max(this->timeManager().time() + this->timeManager().timeStepSize(), 1e-10); - for (; eIt != eEndIt; ++eIt) + for (const auto& element : Dune::elements(this->gridView())) { - fvGeometry.update(this->gridView(), *eIt); + fvGeometry.update(this->gridView(), element); for (int scvIdx = 0; scvIdx < fvGeometry.numScv; ++scvIdx) { - int globalIdx = this->model().dofMapper().subIndex(*eIt, scvIdx, dofCodim); + int globalIdx = this->model().dofMapper().subIndex(element, scvIdx, dofCodim); if (isBox) - globalPos = eIt->geometry().corner(scvIdx); + globalPos = element.geometry().corner(scvIdx); else - globalPos = eIt->geometry().center(); + globalPos = element.geometry().center(); (*temperatureExact)[globalIdx] = temperatureHigh_ + (initialPriVars[temperatureIdx] - temperatureHigh_) *std::erf(0.5*std::sqrt(globalPos[0]*globalPos[0]*storage/time/effectiveThermalConductivity)); diff --git a/test/implicit/1p/1pniconvectionproblem.hh b/test/implicit/1p/1pniconvectionproblem.hh index a0110fbcce..1a7b03781a 100644 --- a/test/implicit/1p/1pniconvectionproblem.hh +++ b/test/implicit/1p/1pniconvectionproblem.hh @@ -131,7 +131,6 @@ class OnePNIConvectionProblem : public ImplicitPorousMediaProblem<TypeTag> typedef typename GridView::template Codim<0>::Entity Element; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; typedef typename GridView::Intersection Intersection; typedef Dune::FieldVector<Scalar, dimWorld> GlobalPosition; @@ -185,9 +184,8 @@ public: FVElementGeometry fvGeometry; VolumeVariables volVars; - ElementIterator eIt = this->gridView().template begin<0>(); - ElementIterator eEndIt = this->gridView().template end<0>(); - fvGeometry.update(this->gridView(), *eIt); + const auto& element = *this->gridView().template begin<0>(); + fvGeometry.update(this->gridView(), element); PrimaryVariables initialPriVars(0); GlobalPosition globalPos(0); initial_(initialPriVars, globalPos); @@ -195,34 +193,34 @@ public: //update the constant volume variables volVars.update(initialPriVars, *this, - *eIt, + element, fvGeometry, 0, false); - Scalar porosity = this->spatialParams().porosity(*eIt, fvGeometry, 0); + Scalar porosity = this->spatialParams().porosity(element, fvGeometry, 0); Scalar densityW = volVars.density(); Scalar heatCapacityW = FluidSystem::heatCapacity(volVars.fluidState(), 0); Scalar storageW = densityW*heatCapacityW*porosity; - Scalar densityS = this->spatialParams().solidDensity(*eIt, fvGeometry, 0); - Scalar heatCapacityS = this->spatialParams().solidHeatCapacity(*eIt, fvGeometry, 0); + Scalar densityS = this->spatialParams().solidDensity(element, fvGeometry, 0); + Scalar heatCapacityS = this->spatialParams().solidHeatCapacity(element, fvGeometry, 0); Scalar storageTotal = storageW + densityS*heatCapacityS*(1 - porosity); std::cout<<"storage: "<<storageTotal<<std::endl; Scalar time = std::max(this->timeManager().time() + this->timeManager().timeStepSize(), 1e-10); Scalar retardedFrontVelocity = darcyVelocity_*storageW/storageTotal/porosity; std::cout<<"retarded velocity: "<<retardedFrontVelocity<<std::endl; - for (; eIt != eEndIt; ++eIt) + for (const auto& element : Dune::elements(this->gridView())) { - fvGeometry.update(this->gridView(), *eIt); + fvGeometry.update(this->gridView(), element); for (int scvIdx = 0; scvIdx < fvGeometry.numScv; ++scvIdx) { - int globalIdx = this->model().dofMapper().subIndex(*eIt, scvIdx, dofCodim); + int globalIdx = this->model().dofMapper().subIndex(element, scvIdx, dofCodim); if (isBox) - globalPos = eIt->geometry().corner(scvIdx); + globalPos = element.geometry().corner(scvIdx); else - globalPos = eIt->geometry().center(); + globalPos = element.geometry().center(); (*temperatureExact)[globalIdx] = globalPos[0] < retardedFrontVelocity*time ? temperatureHigh_ : temperatureLow_; } diff --git a/test/implicit/1p2c/1p2cniconductionproblem.hh b/test/implicit/1p2c/1p2cniconductionproblem.hh index 8f60c9ee40..7bafbbeda6 100644 --- a/test/implicit/1p2c/1p2cniconductionproblem.hh +++ b/test/implicit/1p2c/1p2cniconductionproblem.hh @@ -138,7 +138,6 @@ class OnePTwoCNIConductionProblem : public ImplicitPorousMediaProblem<TypeTag> typedef typename GridView::template Codim<0>::Entity Element; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; typedef typename GridView::Intersection Intersection; typedef Dune::FieldVector<Scalar, dimWorld> GlobalPosition; @@ -191,9 +190,8 @@ public: FVElementGeometry fvGeometry; VolumeVariables volVars; - ElementIterator eIt = this->gridView().template begin<0>(); - ElementIterator eEndIt = this->gridView().template end<0>(); - fvGeometry.update(this->gridView(), *eIt); + const auto& element = *this->gridView().template begin<0>(); + fvGeometry.update(this->gridView(), element); PrimaryVariables initialPriVars(0); GlobalPosition globalPos(0); initial_(initialPriVars, globalPos); @@ -201,33 +199,33 @@ public: //update the constant volume variables volVars.update(initialPriVars, *this, - *eIt, + element, fvGeometry, 0, false); - Scalar porosity = this->spatialParams().porosity(*eIt, fvGeometry, 0); + Scalar porosity = this->spatialParams().porosity(element, fvGeometry, 0); Scalar densityW = volVars.density(); Scalar heatCapacityW = FluidSystem::heatCapacity(volVars.fluidState(), 0); - Scalar densityS = this->spatialParams().solidDensity(*eIt, fvGeometry, 0); - Scalar heatCapacityS = this->spatialParams().solidHeatCapacity(*eIt, fvGeometry, 0); + Scalar densityS = this->spatialParams().solidDensity(element, fvGeometry, 0); + Scalar heatCapacityS = this->spatialParams().solidHeatCapacity(element, fvGeometry, 0); Scalar storage = densityW*heatCapacityW*porosity + densityS*heatCapacityS*(1 - porosity); Scalar effectiveThermalConductivity = ThermalConductivityModel::effectiveThermalConductivity(volVars, this->spatialParams(), - *eIt, fvGeometry, 0); + element, fvGeometry, 0); Scalar time = std::max(this->timeManager().time() + this->timeManager().timeStepSize(), 1e-10); - for (; eIt != eEndIt; ++eIt) + for (const auto& element : Dune::elements(this->gridView())) { - fvGeometry.update(this->gridView(), *eIt); + fvGeometry.update(this->gridView(), element); for (int scvIdx = 0; scvIdx < fvGeometry.numScv; ++scvIdx) { - int dofIdxGlobal = this->model().dofMapper().subIndex(*eIt, scvIdx, dofCodim); + int dofIdxGlobal = this->model().dofMapper().subIndex(element, scvIdx, dofCodim); if (isBox) - globalPos = eIt->geometry().corner(scvIdx); + globalPos = element.geometry().corner(scvIdx); else - globalPos = eIt->geometry().center(); + globalPos = element.geometry().center(); (*temperatureExact)[dofIdxGlobal] = temperatureHigh_ + (initialPriVars[temperatureIdx] - temperatureHigh_) *std::erf(0.5*std::sqrt(globalPos[0]*globalPos[0]*storage/time/effectiveThermalConductivity)); diff --git a/test/implicit/1p2c/1p2cniconvectionproblem.hh b/test/implicit/1p2c/1p2cniconvectionproblem.hh index 19f6f97e31..a35665579b 100644 --- a/test/implicit/1p2c/1p2cniconvectionproblem.hh +++ b/test/implicit/1p2c/1p2cniconvectionproblem.hh @@ -137,7 +137,6 @@ class OnePTwoCNIConvectionProblem : public ImplicitPorousMediaProblem<TypeTag> typedef typename GridView::template Codim<0>::Entity Element; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; typedef typename GridView::Intersection Intersection; typedef Dune::FieldVector<Scalar, dimWorld> GlobalPosition; @@ -194,9 +193,8 @@ public: FVElementGeometry fvGeometry; VolumeVariables volVars; - ElementIterator eIt = this->gridView().template begin<0>(); - ElementIterator eEndIt = this->gridView().template end<0>(); - fvGeometry.update(this->gridView(), *eIt); + const auto& element = *this->gridView().template begin<0>(); + fvGeometry.update(this->gridView(), element); PrimaryVariables initialPriVars(0); GlobalPosition globalPos(0); initial_(initialPriVars, globalPos); @@ -204,34 +202,34 @@ public: //update the constant volume variables volVars.update(initialPriVars, *this, - *eIt, + element, fvGeometry, 0, false); - Scalar porosity = this->spatialParams().porosity(*eIt, fvGeometry, 0); + Scalar porosity = this->spatialParams().porosity(element, fvGeometry, 0); Scalar densityW = volVars.density(); Scalar heatCapacityW = FluidSystem::heatCapacity(volVars.fluidState(), 0); Scalar storageW = densityW*heatCapacityW*porosity; - Scalar densityS = this->spatialParams().solidDensity(*eIt, fvGeometry, 0); - Scalar heatCapacityS = this->spatialParams().solidHeatCapacity(*eIt, fvGeometry, 0); + Scalar densityS = this->spatialParams().solidDensity(element, fvGeometry, 0); + Scalar heatCapacityS = this->spatialParams().solidHeatCapacity(element, fvGeometry, 0); Scalar storageTotal = storageW + densityS*heatCapacityS*(1 - porosity); std::cout<<"storage: "<<storageTotal<<std::endl; Scalar time = std::max(this->timeManager().time() + this->timeManager().timeStepSize(), 1e-10); Scalar retardedFrontVelocity = darcyVelocity_*storageW/storageTotal/porosity; std::cout<<"retarded velocity: "<<retardedFrontVelocity<<std::endl; - for (; eIt != eEndIt; ++eIt) + for (const auto& element : Dune::elements(this->gridView())) { - fvGeometry.update(this->gridView(), *eIt); + fvGeometry.update(this->gridView(), element); for (int scvIdx = 0; scvIdx < fvGeometry.numScv; ++scvIdx) { - int dofIdxGlobal = this->model().dofMapper().subIndex(*eIt, scvIdx, dofCodim); + int dofIdxGlobal = this->model().dofMapper().subIndex(element, scvIdx, dofCodim); if (isBox) - globalPos = eIt->geometry().corner(scvIdx); + globalPos = element.geometry().corner(scvIdx); else - globalPos = eIt->geometry().center(); + globalPos = element.geometry().center(); (*temperatureExact)[dofIdxGlobal] = globalPos[0] < retardedFrontVelocity*time ? temperatureHigh_ : temperatureLow_; } diff --git a/test/implicit/2p/cc2pcornerpointproblem.hh b/test/implicit/2p/cc2pcornerpointproblem.hh index aead6d6379..3b06af54f1 100644 --- a/test/implicit/2p/cc2pcornerpointproblem.hh +++ b/test/implicit/2p/cc2pcornerpointproblem.hh @@ -86,7 +86,6 @@ class CC2PCornerPointProblem : public ImplicitPorousMediaProblem<TypeTag> typedef ImplicitPorousMediaProblem<TypeTag> ParentType; typedef typename GET_PROP_TYPE(TypeTag, GridCreator) GridCreator; typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; typedef typename GridView::template Codim<0>::Entity Element; typedef typename GridView::Intersection Intersection; typedef typename GET_PROP_TYPE(TypeTag, FVElementGeometry) FVElementGeometry; @@ -314,15 +313,13 @@ public: ScalarField *permX = this->resultWriter().allocateManagedBuffer(numElements); ScalarField *permZ = this->resultWriter().allocateManagedBuffer(numElements); - ElementIterator eIt = this->gridView().template begin<0>(); - ElementIterator eEndIt = this->gridView().template end<0>(); - for (; eIt != eEndIt; ++eIt) + for (const auto& element : Dune::elements(this->gridView())) { FVElementGeometry fvGeometry; - fvGeometry.update(this->gridView(), *eIt); - int eIdx = this->elementMapper().map(*eIt); + fvGeometry.update(this->gridView(), element); + int eIdx = this->elementMapper().map(element); - const DimWorldMatrix K = this->spatialParams().intrinsicPermeability(*eIt, fvGeometry, /*element data*/ 0); + const DimWorldMatrix K = this->spatialParams().intrinsicPermeability(element, fvGeometry, /*element data*/ 0); // transfer output to mD = 9.86923e-16 m^2 (*permX)[eIdx] = K[0][0]/9.86923e-16; (*permZ)[eIdx] = K[2][2]/9.86923e-16; diff --git a/test/implicit/2pnc/fuelcellproblem.hh b/test/implicit/2pnc/fuelcellproblem.hh index e735751b49..f61c35938f 100644 --- a/test/implicit/2pnc/fuelcellproblem.hh +++ b/test/implicit/2pnc/fuelcellproblem.hh @@ -322,24 +322,24 @@ public: ScalarField *reactionSourceH2O = this->resultWriter().allocateManagedBuffer (numDofs); ScalarField *reactionSourceO2 = this->resultWriter().allocateManagedBuffer (numDofs); - for (auto eIt = this->gridView().template begin<0>(); eIt != this->gridView().template end<0>(); ++eIt) + for (const auto& element : Dune::elements(this->gridView())) { FVElementGeometry fvGeometry; - fvGeometry.update(this->gridView(), *eIt); + fvGeometry.update(this->gridView(), element); ElementVolumeVariables elemVolVars; elemVolVars.update(*this, - *eIt, + element, fvGeometry, false /* oldSol? */); for (int scvIdx = 0; scvIdx < fvGeometry.numScv; ++scvIdx) { - int dofIdxGlobal = this->model().dofMapper().subIndex(*eIt, scvIdx, dofCodim); + int dofIdxGlobal = this->model().dofMapper().subIndex(element, scvIdx, dofCodim); //reactionSource Output PrimaryVariables source; - this->solDependentSource(source, *eIt, fvGeometry, scvIdx, elemVolVars); + this->solDependentSource(source, element, fvGeometry, scvIdx, elemVolVars); (*reactionSourceH2O)[dofIdxGlobal] = source[wPhaseIdx]; (*reactionSourceO2)[dofIdxGlobal] = source[numComponents-1]; diff --git a/test/implicit/3p/3pniconductionproblem.hh b/test/implicit/3p/3pniconductionproblem.hh index ac8fc04070..63d1feee53 100644 --- a/test/implicit/3p/3pniconductionproblem.hh +++ b/test/implicit/3p/3pniconductionproblem.hh @@ -128,7 +128,6 @@ class ThreePNIConductionProblem : public ImplicitPorousMediaProblem<TypeTag> typedef typename GridView::template Codim<0>::Entity Element; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; typedef typename GridView::Intersection Intersection; typedef Dune::FieldVector<Scalar, dimWorld> GlobalPosition; @@ -180,9 +179,8 @@ public: FVElementGeometry fvGeometry; VolumeVariables volVars; - ElementIterator eIt = this->gridView().template begin<0>(); - ElementIterator eEndIt = this->gridView().template end<0>(); - fvGeometry.update(this->gridView(), *eIt); + const auto& element = this->gridView().template begin<0>(); + fvGeometry.update(this->gridView(), element); PrimaryVariables initialPriVars(0); GlobalPosition globalPos(0); initial_(initialPriVars, globalPos); @@ -190,33 +188,33 @@ public: //update the constant volume variables volVars.update(initialPriVars, *this, - *eIt, + element, fvGeometry, 0, false); - Scalar porosity = this->spatialParams().porosity(*eIt, fvGeometry, 0); + Scalar porosity = this->spatialParams().porosity(element, fvGeometry, 0); Scalar densityW = volVars.density(swIdx); Scalar heatCapacityW = IapwsH2O::liquidHeatCapacity(initialPriVars[temperatureIdx], initialPriVars[pressureIdx]); - Scalar densityS = this->spatialParams().solidDensity(*eIt, fvGeometry, 0); - Scalar heatCapacityS = this->spatialParams().solidHeatCapacity(*eIt, fvGeometry, 0); + Scalar densityS = this->spatialParams().solidDensity(element, fvGeometry, 0); + Scalar heatCapacityS = this->spatialParams().solidHeatCapacity(element, fvGeometry, 0); Scalar storage = densityW*heatCapacityW*porosity + densityS*heatCapacityS*(1 - porosity); Scalar effectiveThermalConductivity = ThermalConductivityModel::effectiveThermalConductivity(volVars, this->spatialParams(), - *eIt, fvGeometry, 0); + element, fvGeometry, 0); Scalar time = std::max(this->timeManager().time() + this->timeManager().timeStepSize(), 1e-10); - for (; eIt != eEndIt; ++eIt) + for (const auto& element : Dune::elements(this->gridView())) { - fvGeometry.update(this->gridView(), *eIt); + fvGeometry.update(this->gridView(), element); for (int scvIdx = 0; scvIdx < fvGeometry.numScv; ++scvIdx) { - int globalIdx = this->model().dofMapper().subIndex(*eIt, scvIdx, dofCodim); + int globalIdx = this->model().dofMapper().subIndex(element, scvIdx, dofCodim); if (isBox) - globalPos = eIt->geometry().corner(scvIdx); + globalPos = element.geometry().corner(scvIdx); else - globalPos = eIt->geometry().center(); + globalPos = element.geometry().center(); (*temperatureExact)[globalIdx] = temperatureHigh_ + (initialPriVars[temperatureIdx] - temperatureHigh_) *std::erf(0.5*std::sqrt(globalPos[0]*globalPos[0]*storage/time/effectiveThermalConductivity)); diff --git a/test/implicit/3p/3pniconvectionproblem.hh b/test/implicit/3p/3pniconvectionproblem.hh index 3cc8aaebbf..57631bc7e1 100644 --- a/test/implicit/3p/3pniconvectionproblem.hh +++ b/test/implicit/3p/3pniconvectionproblem.hh @@ -128,7 +128,6 @@ class ThreePNIConvectionProblem : public ImplicitPorousMediaProblem<TypeTag> typedef typename GridView::template Codim<0>::Entity Element; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; typedef typename GridView::Intersection Intersection; typedef Dune::FieldVector<Scalar, dimWorld> GlobalPosition; @@ -185,9 +184,8 @@ public: FVElementGeometry fvGeometry; VolumeVariables volVars; - ElementIterator eIt = this->gridView().template begin<0>(); - ElementIterator eEndIt = this->gridView().template end<0>(); - fvGeometry.update(this->gridView(), *eIt); + const auto& element = *this->gridView().template begin<0>(); + fvGeometry.update(this->gridView(), element); PrimaryVariables initialPriVars(0); GlobalPosition globalPos(0); initial_(initialPriVars, globalPos); @@ -195,34 +193,34 @@ public: //update the constant volume variables volVars.update(initialPriVars, *this, - *eIt, + element, fvGeometry, 0, false); - Scalar porosity = this->spatialParams().porosity(*eIt, fvGeometry, 0); + Scalar porosity = this->spatialParams().porosity(element, fvGeometry, 0); Scalar densityW = volVars.density(swIdx); Scalar heatCapacityW = IapwsH2O::liquidHeatCapacity(initialPriVars[temperatureIdx], initialPriVars[pressureIdx]); Scalar storageW = densityW*heatCapacityW*porosity; - Scalar densityS = this->spatialParams().solidDensity(*eIt, fvGeometry, 0); - Scalar heatCapacityS = this->spatialParams().solidHeatCapacity(*eIt, fvGeometry, 0); + Scalar densityS = this->spatialParams().solidDensity(element, fvGeometry, 0); + Scalar heatCapacityS = this->spatialParams().solidHeatCapacity(element, fvGeometry, 0); Scalar storageTotal = storageW + densityS*heatCapacityS*(1 - porosity); std::cout<<"storage: "<<storageTotal<<std::endl; Scalar time = std::max(this->timeManager().time() + this->timeManager().timeStepSize(), 1e-10); Scalar retardedFrontVelocity = darcyVelocity_*storageW/storageTotal/porosity; std::cout<<"retarded velocity: "<<retardedFrontVelocity<<std::endl; - for (; eIt != eEndIt; ++eIt) + for (const auto& element : Dune::elements(this->gridView())) { - fvGeometry.update(this->gridView(), *eIt); + fvGeometry.update(this->gridView(), element); for (int scvIdx = 0; scvIdx < fvGeometry.numScv; ++scvIdx) { - int globalIdx = this->model().dofMapper().subIndex(*eIt, scvIdx, dofCodim); + int globalIdx = this->model().dofMapper().subIndex(element, scvIdx, dofCodim); if (isBox) - globalPos = eIt->geometry().corner(scvIdx); + globalPos = element.geometry().corner(scvIdx); else - globalPos = eIt->geometry().center(); + globalPos = element.geometry().center(); (*temperatureExact)[globalIdx] = globalPos[0] < retardedFrontVelocity*time ? temperatureHigh_ : temperatureLow_; } diff --git a/test/implicit/3p3c/columnxylolproblem.hh b/test/implicit/3p3c/columnxylolproblem.hh index 36e48f81e6..7460490480 100644 --- a/test/implicit/3p3c/columnxylolproblem.hh +++ b/test/implicit/3p3c/columnxylolproblem.hh @@ -119,7 +119,6 @@ class ColumnProblem : public ImplicitPorousMediaProblem<TypeTag> typedef typename GET_PROP_TYPE(TypeTag, TimeManager) TimeManager; typedef typename GridView::template Codim<0>::Entity Element; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; typedef typename GridView::template Codim<dim>::Entity Vertex; typedef typename GridView::Intersection Intersection; @@ -293,16 +292,14 @@ public: FVElementGeometry fvGeometry; - ElementIterator eIt = this->gridView().template begin<0>(); - ElementIterator eEndIt = this->gridView().template end<0>(); - for (; eIt != eEndIt; ++eIt) + for (const auto& element : Dune::elements(this->gridView())) { - fvGeometry.update(this->gridView(), *eIt); + fvGeometry.update(this->gridView(), element); for (int scvIdx = 0; scvIdx < fvGeometry.numScv; ++scvIdx) { - int dofIdxGlobal = this->model().dofMapper().subIndex(*eIt, scvIdx, dofCodim); - (*Kxx)[dofIdxGlobal] = this->spatialParams().intrinsicPermeability(*eIt, fvGeometry, scvIdx); + int dofIdxGlobal = this->model().dofMapper().subIndex(element, scvIdx, dofCodim); + (*Kxx)[dofIdxGlobal] = this->spatialParams().intrinsicPermeability(element, fvGeometry, scvIdx); } } diff --git a/test/implicit/3p3c/infiltrationproblem.hh b/test/implicit/3p3c/infiltrationproblem.hh index af05330f23..a2a07d57fa 100644 --- a/test/implicit/3p3c/infiltrationproblem.hh +++ b/test/implicit/3p3c/infiltrationproblem.hh @@ -132,7 +132,6 @@ class InfiltrationProblem : public ImplicitPorousMediaProblem<TypeTag> typedef typename GET_PROP_TYPE(TypeTag, TimeManager) TimeManager; typedef typename GridView::template Codim<0>::Entity Element; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; typedef typename GridView::template Codim<dim>::Entity Vertex; typedef typename GridView::Intersection Intersection; @@ -352,16 +351,14 @@ public: FVElementGeometry fvGeometry; - ElementIterator eIt = this->gridView().template begin<0>(); - ElementIterator eEndIt = this->gridView().template end<0>(); - for (; eIt != eEndIt; ++eIt) + for (const auto& element : Dune::elements(this->gridView())) { - fvGeometry.update(this->gridView(), *eIt); + fvGeometry.update(this->gridView(), element); for (int scvIdx = 0; scvIdx < fvGeometry.numScv; ++scvIdx) { - int dofIdxGlobal = this->model().dofMapper().subIndex(*eIt, scvIdx, dofCodim); - (*Kxx)[dofIdxGlobal] = this->spatialParams().intrinsicPermeability(*eIt, fvGeometry, scvIdx); + int dofIdxGlobal = this->model().dofMapper().subIndex(element, scvIdx, dofCodim); + (*Kxx)[dofIdxGlobal] = this->spatialParams().intrinsicPermeability(element, fvGeometry, scvIdx); } } diff --git a/test/implicit/3p3c/kuevetteproblem.hh b/test/implicit/3p3c/kuevetteproblem.hh index 5708c8a330..5da7d3aab7 100644 --- a/test/implicit/3p3c/kuevetteproblem.hh +++ b/test/implicit/3p3c/kuevetteproblem.hh @@ -132,7 +132,6 @@ class KuevetteProblem : public ImplicitPorousMediaProblem<TypeTag> typedef typename GET_PROP_TYPE(TypeTag, TimeManager) TimeManager; typedef typename GridView::template Codim<0>::Entity Element; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; typedef typename GridView::template Codim<dim>::Entity Vertex; typedef typename GridView::Intersection Intersection; @@ -300,16 +299,14 @@ public: FVElementGeometry fvGeometry; - ElementIterator eIt = this->gridView().template begin<0>(); - ElementIterator eEndIt = this->gridView().template end<0>(); - for (; eIt != eEndIt; ++eIt) + for (const auto& element : Dune::elements(this->gridView())) { - fvGeometry.update(this->gridView(), *eIt); + fvGeometry.update(this->gridView(), element); for (int scvIdx = 0; scvIdx < fvGeometry.numScv; ++scvIdx) { - int dofIdxGlobal = this->model().dofMapper().subIndex(*eIt, scvIdx, dofCodim); - (*Kxx)[dofIdxGlobal] = this->spatialParams().intrinsicPermeability(*eIt, fvGeometry, scvIdx); + int dofIdxGlobal = this->model().dofMapper().subIndex(element, scvIdx, dofCodim); + (*Kxx)[dofIdxGlobal] = this->spatialParams().intrinsicPermeability(element, fvGeometry, scvIdx); } } diff --git a/test/implicit/co2/heterogeneousproblem.hh b/test/implicit/co2/heterogeneousproblem.hh index 013b1bc744..1e29643e21 100644 --- a/test/implicit/co2/heterogeneousproblem.hh +++ b/test/implicit/co2/heterogeneousproblem.hh @@ -147,7 +147,6 @@ class HeterogeneousProblem : public ImplicitPorousMediaProblem<TypeTag> typedef typename GET_PROP_TYPE(TypeTag, TimeManager) TimeManager; typedef typename GridView::template Codim<0>::Entity Element; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; typedef typename GridView::template Codim<dim>::Entity Vertex; typedef typename GridView::Intersection Intersection; @@ -267,28 +266,26 @@ public: FVElementGeometry fvGeometry; VolumeVariables volVars; - ElementIterator eIt = this->gridView().template begin<0>(); - ElementIterator eEndIt = this->gridView().template end<0>(); - for (; eIt != eEndIt; ++eIt) + for (const auto& element : Dune::elements(this->gridView())) { - int eIdx = this->elementMapper().map(*eIt); + int eIdx = this->elementMapper().map(element); (*rank)[eIdx] = this->gridView().comm().rank(); - fvGeometry.update(this->gridView(), *eIt); + fvGeometry.update(this->gridView(), element); for (int scvIdx = 0; scvIdx < fvGeometry.numScv; ++scvIdx) { - int dofIdxGlobal = this->model().dofMapper().map(*eIt, scvIdx, dofCodim); + int dofIdxGlobal = this->model().dofMapper().map(element, scvIdx, dofCodim); volVars.update(this->model().curSol()[dofIdxGlobal], *this, - *eIt, + element, fvGeometry, scvIdx, false); (*boxVolume)[dofIdxGlobal] += fvGeometry.subContVol[scvIdx].volume; } - (*Kxx)[eIdx] = this->spatialParams().intrinsicPermeability(*eIt, fvGeometry, /*element data*/ 0); - (*cellPorosity)[eIdx] = this->spatialParams().porosity(*eIt, fvGeometry, /*element data*/ 0); + (*Kxx)[eIdx] = this->spatialParams().intrinsicPermeability(element, fvGeometry, /*element data*/ 0); + (*cellPorosity)[eIdx] = this->spatialParams().porosity(element, fvGeometry, /*element data*/ 0); } //pass the scalar fields to the vtkwriter diff --git a/test/implicit/co2/heterogeneousproblemni.hh b/test/implicit/co2/heterogeneousproblemni.hh index 89fa68194d..2ecfa6d60e 100644 --- a/test/implicit/co2/heterogeneousproblemni.hh +++ b/test/implicit/co2/heterogeneousproblemni.hh @@ -151,7 +151,6 @@ class HeterogeneousNIProblem : public ImplicitPorousMediaProblem<TypeTag> typedef typename GET_PROP_TYPE(TypeTag, TimeManager) TimeManager; typedef typename GridView::template Codim<0>::Entity Element; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; typedef typename GridView::template Codim<dim>::Entity Vertex; typedef typename GridView::Intersection Intersection; @@ -274,21 +273,19 @@ public: FVElementGeometry fvGeometry; VolumeVariables volVars; - ElementIterator eIt = this->gridView().template begin<0>(); - ElementIterator eEndIt = this->gridView().template end<0>(); - for (; eIt != eEndIt; ++eIt) + for (const auto& element : Dune::elements(this->gridView())) { - int eIdx = this->elementMapper().map(*eIt); + int eIdx = this->elementMapper().map(element); (*rank)[eIdx] = this->gridView().comm().rank(); - fvGeometry.update(this->gridView(), *eIt); + fvGeometry.update(this->gridView(), element); for (int scvIdx = 0; scvIdx < fvGeometry.numScv; ++scvIdx) { - int dofIdxGlobal = this->model().dofMapper().map(*eIt, scvIdx, dofCodim); + int dofIdxGlobal = this->model().dofMapper().map(element, scvIdx, dofCodim); volVars.update(this->model().curSol()[dofIdxGlobal], *this, - *eIt, + element, fvGeometry, scvIdx, false); @@ -296,8 +293,8 @@ public: (*enthalpyW)[dofIdxGlobal] = volVars.enthalpy(lPhaseIdx); (*enthalpyN)[dofIdxGlobal] = volVars.enthalpy(gPhaseIdx); } - (*Kxx)[eIdx] = this->spatialParams().intrinsicPermeability(*eIt, fvGeometry, /*element data*/ 0); - (*cellPorosity)[eIdx] = this->spatialParams().porosity(*eIt, fvGeometry, /*element data*/ 0); + (*Kxx)[eIdx] = this->spatialParams().intrinsicPermeability(element, fvGeometry, /*element data*/ 0); + (*cellPorosity)[eIdx] = this->spatialParams().porosity(element, fvGeometry, /*element data*/ 0); } //pass the scalar fields to the vtkwriter diff --git a/test/implicit/co2/heterogeneousspatialparameters.hh b/test/implicit/co2/heterogeneousspatialparameters.hh index ed7a33ee28..91a2ef264c 100644 --- a/test/implicit/co2/heterogeneousspatialparameters.hh +++ b/test/implicit/co2/heterogeneousspatialparameters.hh @@ -83,7 +83,6 @@ class HeterogeneousSpatialParams : public ImplicitSpatialParams<TypeTag> typedef typename GET_PROP_TYPE(TypeTag, FVElementGeometry) FVElementGeometry; typedef typename GridView::template Codim<0>::Entity Element; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; public: typedef typename GET_PROP_TYPE(TypeTag, MaterialLaw) MaterialLaw; @@ -136,12 +135,10 @@ public: int numElements = gridView_.size(0); paramIdx_.resize(numElements); - ElementIterator eIt = gridView_.template begin<0>(); - const ElementIterator eEndIt = gridView_.template end<0>(); - for (; eIt != eEndIt; ++eIt) + for (const auto& element : Dune::elements(gridView_)) { - int eIdx = gridView_.indexSet().index(*eIt); - int param = GridCreator::parameters(*eIt)[0]; + int eIdx = gridView_.indexSet().index(element); + int param = GridCreator::parameters(element)[0]; paramIdx_[eIdx] = param; } } diff --git a/test/implicit/richards/richardsanalyticalproblem.hh b/test/implicit/richards/richardsanalyticalproblem.hh index 312c5a495a..5f50ebd023 100644 --- a/test/implicit/richards/richardsanalyticalproblem.hh +++ b/test/implicit/richards/richardsanalyticalproblem.hh @@ -342,15 +342,14 @@ public: Scalar l2analytic = 0.0; const Scalar time = this->timeManager().time() + this->timeManager().timeStepSize(); - for (auto elementIt = this->gridView().template begin<0>(); - elementIt != this->gridView().template end<0>(); ++elementIt) + for (const auto& element : Dune::elements(this->gridView())) { // value from numerical approximation Scalar numericalSolution = - this->model().curSol()[this->model().dofMapper().subIndex(*elementIt, 0, 0)]; + this->model().curSol()[this->model().dofMapper().subIndex(element, 0, 0)]; // integrate over element using a quadrature rule - Geometry geometry = elementIt->geometry(); + Geometry geometry = element.geometry(); Dune::GeometryType gt = geometry.type(); Dune::QuadratureRule<Scalar, dim> rule = Dune::QuadratureRules<Scalar, dim>::rule(gt, qOrder); diff --git a/test/implicit/richards/richardsniconductionproblem.hh b/test/implicit/richards/richardsniconductionproblem.hh index 111a834901..dee59b0297 100644 --- a/test/implicit/richards/richardsniconductionproblem.hh +++ b/test/implicit/richards/richardsniconductionproblem.hh @@ -124,7 +124,6 @@ class RichardsNIConductionProblem : public RichardsProblem<TypeTag> }; typedef typename GridView::template Codim<0>::Entity Element; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; typedef typename GridView::Intersection Intersection; typedef Dune::FieldVector<Scalar, dimWorld> GlobalPosition; @@ -169,37 +168,36 @@ public: FVElementGeometry fvGeometry; VolumeVariables volVars; - ElementIterator eIt = this->gridView().template begin<0>(); - ElementIterator eEndIt = this->gridView().template end<0>(); - fvGeometry.update(this->gridView(), *eIt); + const auto& element = *this->gridView().template begin<0>(); + fvGeometry.update(this->gridView(), element); PrimaryVariables initialPriVars(0); GlobalPosition globalPos(0); initial_(initialPriVars, globalPos); //update the constant volume variables - volVars.update(initialPriVars, *this, *eIt, fvGeometry, 0, false); + volVars.update(initialPriVars, *this, element, fvGeometry, 0, false); - Scalar porosity = this->spatialParams().porosity(*eIt, fvGeometry, 0); + Scalar porosity = this->spatialParams().porosity(element, fvGeometry, 0); Scalar densityW = volVars.density(wPhaseIdx); Scalar heatCapacityW = FluidSystem::heatCapacity(volVars.fluidState(), 0); - Scalar densityS = this->spatialParams().solidDensity(*eIt, fvGeometry, 0); - Scalar heatCapacityS = this->spatialParams().solidHeatCapacity(*eIt, fvGeometry, 0); + Scalar densityS = this->spatialParams().solidDensity(element, fvGeometry, 0); + Scalar heatCapacityS = this->spatialParams().solidHeatCapacity(element, fvGeometry, 0); Scalar storage = densityW*heatCapacityW*porosity + densityS*heatCapacityS*(1 - porosity); Scalar effectiveThermalConductivity = ThermalConductivityModel::effectiveThermalConductivity(volVars, this->spatialParams(), - *eIt, fvGeometry, 0); + element, fvGeometry, 0); Scalar time = std::max(this->timeManager().time() + this->timeManager().timeStepSize(), 1e-10); - for (; eIt != eEndIt; ++eIt) + for (const auto& element : Dune::elements(this->gridView())) { - fvGeometry.update(this->gridView(), *eIt); + fvGeometry.update(this->gridView(), element); for (int scvIdx = 0; scvIdx < fvGeometry.numScv; ++scvIdx) { - int globalIdx = this->model().dofMapper().subIndex(*eIt, scvIdx, dofCodim); + int globalIdx = this->model().dofMapper().subIndex(element, scvIdx, dofCodim); if (isBox) - globalPos = eIt->geometry().corner(scvIdx); + globalPos = element.geometry().corner(scvIdx); else - globalPos = eIt->geometry().center(); + globalPos = element.geometry().center(); (*temperatureExact)[globalIdx] = temperatureHigh_ + (initialPriVars[temperatureIdx] - temperatureHigh_) *std::erf(0.5*std::sqrt(globalPos[0]*globalPos[0]*storage/time/effectiveThermalConductivity)); diff --git a/test/implicit/richards/richardsniconvectionproblem.hh b/test/implicit/richards/richardsniconvectionproblem.hh index 5c42afdb1f..08df81a1a4 100644 --- a/test/implicit/richards/richardsniconvectionproblem.hh +++ b/test/implicit/richards/richardsniconvectionproblem.hh @@ -129,7 +129,6 @@ class RichardsNIConvectionProblem : public RichardsProblem<TypeTag> }; typedef typename GridView::template Codim<0>::Entity Element; - typedef typename GridView::template Codim<0>::Iterator ElementIterator; typedef typename GridView::Intersection Intersection; typedef Dune::FieldVector<Scalar, dimWorld> GlobalPosition; @@ -178,39 +177,38 @@ public: FVElementGeometry fvGeometry; VolumeVariables volVars; - ElementIterator eIt = this->gridView().template begin<0>(); - ElementIterator eEndIt = this->gridView().template end<0>(); - fvGeometry.update(this->gridView(), *eIt); + const auto& element = *this->gridView().template begin<0>(); + fvGeometry.update(this->gridView(), element); PrimaryVariables initialPriVars(0); GlobalPosition globalPos(0); initial_(initialPriVars, globalPos); //update the constant volume variables - volVars.update(initialPriVars, *this, *eIt, fvGeometry, 0, false); + volVars.update(initialPriVars, *this, element, fvGeometry, 0, false); - Scalar porosity = this->spatialParams().porosity(*eIt, fvGeometry, 0); + Scalar porosity = this->spatialParams().porosity(element, fvGeometry, 0); Scalar densityW = volVars.density(wPhaseIdx); Scalar heatCapacityW = FluidSystem::heatCapacity(volVars.fluidState(), 0); Scalar storageW = densityW*heatCapacityW*porosity; - Scalar densityS = this->spatialParams().solidDensity(*eIt, fvGeometry, 0); - Scalar heatCapacityS = this->spatialParams().solidHeatCapacity(*eIt, fvGeometry, 0); + Scalar densityS = this->spatialParams().solidDensity(element, fvGeometry, 0); + Scalar heatCapacityS = this->spatialParams().solidHeatCapacity(element, fvGeometry, 0); Scalar storageTotal = storageW + densityS*heatCapacityS*(1 - porosity); std::cout<<"storage: "<<storageTotal<<std::endl; Scalar time = std::max(this->timeManager().time() + this->timeManager().timeStepSize(), 1e-10); Scalar retardedFrontVelocity = darcyVelocity_*storageW/storageTotal/porosity; std::cout<<"retarded velocity: "<<retardedFrontVelocity<<std::endl; - for (; eIt != eEndIt; ++eIt) + for (const auto& element : Dune::elements(this->gridView())) { - fvGeometry.update(this->gridView(), *eIt); + fvGeometry.update(this->gridView(), element); for (int scvIdx = 0; scvIdx < fvGeometry.numScv; ++scvIdx) { - int globalIdx = this->model().dofMapper().subIndex(*eIt, scvIdx, dofCodim); + int globalIdx = this->model().dofMapper().subIndex(element, scvIdx, dofCodim); if (isBox) - globalPos = eIt->geometry().corner(scvIdx); + globalPos = element.geometry().corner(scvIdx); else - globalPos = eIt->geometry().center(); + globalPos = element.geometry().center(); (*temperatureExact)[globalIdx] = globalPos[0] < retardedFrontVelocity*time ? temperatureHigh_ : temperatureLow_; } -- GitLab