diff --git a/dune/solvers/transferoperators/genericmultigridtransfer.hh b/dune/solvers/transferoperators/genericmultigridtransfer.hh index 60313b54ef929c02b3a1e0b2dcddce0e3030439f..64352907d7f8edab35a3507727e3b5ee78be1a8e 100644 --- a/dune/solvers/transferoperators/genericmultigridtransfer.hh +++ b/dune/solvers/transferoperators/genericmultigridtransfer.hh @@ -88,46 +88,34 @@ public: matrix.setSize(rows,cols); matrix.setBuildMode(TransferOperatorType::random); - typedef typename GridType::template Codim<0>::LevelIterator ElementIterator; - - typename GridType::LevelGridView levelView = grid.levelGridView(cL); - ElementIterator cIt = levelView.template begin<0>(); - ElementIterator cEndIt = levelView.template end<0>(); - - // /////////////////////////////////////////// // Determine the indices present in the matrix // ///////////////////////////////////////////////// Dune::MatrixIndexSet indices(rows, cols); - for (; cIt != cEndIt; ++cIt) { + for (const auto& cIt : elements(grid.levelGridView(cL))) { typedef typename GridType::template Codim<0>::Entity EntityType; - typedef typename EntityType::HierarchicIterator HierarchicIterator; // Get local finite element - const FEType& coarseBaseSet = p1FECache.get(cIt->type()); + const FEType& coarseBaseSet = p1FECache.get(cIt.type()); const size_t numCoarseBaseFct = coarseBaseSet.localBasis().size(); // preallocate vector for function evaluations std::vector<Dune::FieldVector<field_type,1> > values(numCoarseBaseFct); - HierarchicIterator fIt = cIt->hbegin(fL); - HierarchicIterator fEndIt = cIt->hend(fL); - - for (; fIt != fEndIt; ++fIt) { + for (const auto& fIt : descendantElements(cIt,fL)) { - if (fIt->level()==cIt->level()) + if (fIt.level()==cIt.level()) continue; - const Dune::ReferenceElement<ctype,dim>& fineRefElement - = Dune::ReferenceElements<ctype, dim>::general(fIt->type()); + const auto& fineRefElement = Dune::ReferenceElements<ctype, dim>::general(fIt.type()); - const typename EntityType::LocalGeometry& fGeometryInFather = fIt->geometryInFather(); + const typename EntityType::LocalGeometry& fGeometryInFather = fIt.geometryInFather(); // Get local finite element - const FEType& fineBaseSet = p1FECache.get(fIt->type()); + const FEType& fineBaseSet = p1FECache.get(fIt.type()); const size_t numFineBaseFct = fineBaseSet.localBasis().size(); @@ -136,7 +124,7 @@ public: { const Dune::LocalKey& jLocalKey = fineBaseSet.localCoefficients().localKey(j); - int globalFine = fineIndexSet.subIndex(*fIt, jLocalKey.subEntity(), jLocalKey.codim()); + int globalFine = fineIndexSet.subIndex(fIt, jLocalKey.subEntity(), jLocalKey.codim()); Dune::FieldVector<ctype, dim> fineBasePosition = fineRefElement.position(jLocalKey.subEntity(), jLocalKey.codim()); Dune::FieldVector<ctype, dim> local = fGeometryInFather.global(fineBasePosition); @@ -149,7 +137,7 @@ public: if (values[i] > 0.001) { const Dune::LocalKey& iLocalKey = coarseBaseSet.localCoefficients().localKey(i); - int globalCoarse = coarseIndexSet.subIndex(*cIt, iLocalKey.subEntity(), iLocalKey.codim()); + int globalCoarse = coarseIndexSet.subIndex(cIt, iLocalKey.subEntity(), iLocalKey.codim()); indices.add(globalFine, globalCoarse); } } @@ -164,35 +152,30 @@ public: // ///////////////////////////////////////////// // Compute the matrix // ///////////////////////////////////////////// - cIt = levelView.template begin<0>(); - for (; cIt != cEndIt; ++cIt) { + + for (const auto& cIt : elements(grid.levelGridView(cL))) { // Get local finite element - const FEType& coarseBaseSet = p1FECache.get(cIt->type()); + const FEType& coarseBaseSet = p1FECache.get(cIt.type()); const size_t numCoarseBaseFct = coarseBaseSet.localBasis().size(); typedef typename GridType::template Codim<0>::Entity EntityType; - typedef typename EntityType::HierarchicIterator HierarchicIterator; // preallocate vector for function evaluations std::vector<Dune::FieldVector<field_type,1> > values(numCoarseBaseFct); - HierarchicIterator fIt = cIt->hbegin(fL); - HierarchicIterator fEndIt = cIt->hend(fL); - - for (; fIt != fEndIt; ++fIt) { + for (const auto& fIt : descendantElements(cIt,fL)) { - if (fIt->level()==cIt->level()) + if (fIt.level()==cIt.level()) continue; - const Dune::ReferenceElement<ctype,dim>& fineRefElement - = Dune::ReferenceElements<ctype, dim>::general(fIt->type()); + const auto& fineRefElement = Dune::ReferenceElements<ctype, dim>::general(fIt.type()); - const typename EntityType::LocalGeometry& fGeometryInFather = fIt->geometryInFather(); + const typename EntityType::LocalGeometry& fGeometryInFather = fIt.geometryInFather(); // Get local finite element - const FEType& fineBaseSet = p1FECache.get(fIt->type()); + const FEType& fineBaseSet = p1FECache.get(fIt.type()); const size_t numFineBaseFct = fineBaseSet.localBasis().size(); @@ -200,7 +183,7 @@ public: { const Dune::LocalKey& jLocalKey = fineBaseSet.localCoefficients().localKey(j); - int globalFine = fineIndexSet.subIndex(*fIt, jLocalKey.subEntity(), jLocalKey.codim()); + int globalFine = fineIndexSet.subIndex(fIt, jLocalKey.subEntity(), jLocalKey.codim()); Dune::FieldVector<ctype, dim> fineBasePosition = fineRefElement.position(jLocalKey.subEntity(), jLocalKey.codim()); Dune::FieldVector<ctype, dim> local = fGeometryInFather.global(fineBasePosition); @@ -213,16 +196,14 @@ public: if (values[i] > 0.001) { const Dune::LocalKey& iLocalKey = coarseBaseSet.localCoefficients().localKey(i); - int globalCoarse = coarseIndexSet.subIndex(*cIt, iLocalKey.subEntity(), iLocalKey.codim()); + int globalCoarse = coarseIndexSet.subIndex(cIt, iLocalKey.subEntity(), iLocalKey.codim()); matrix[globalFine][globalCoarse] = Dune::ScaledIdentityMatrix<ctype,TransferMatrixBlock::rows>(values[i]); } } } } - } - } @@ -237,7 +218,6 @@ public: { typedef typename TransferOperatorType::block_type TransferMatrixBlock; typedef typename GridViewCoarse::ctype ctype; - typedef typename GridViewCoarse::Grid GridType; const int dim = GridViewCoarse::dimension; @@ -255,24 +235,18 @@ public: matrix.setSize(rows,cols); matrix.setBuildMode(TransferOperatorType::random); - typedef typename GridViewFine::template Codim<0>::Iterator ElementIterator; - - ElementIterator fIt = gvFine.template begin<0>(); - ElementIterator fEndIt = gvFine.template end<0>(); - - // /////////////////////////////////////////// // Determine the indices present in the matrix // ///////////////////////////////////////////////// + Dune::MatrixIndexSet indices(rows, cols); - for (; fIt != fEndIt; ++fIt) { + for (const auto& fIt : elements(gvFine)) { - const Dune::ReferenceElement<ctype,dim>& fineRefElement - = Dune::ReferenceElements<ctype, dim>::general(fIt->type()); + const auto& fineRefElement = Dune::ReferenceElements<ctype, dim>::general(fIt.type()); // Get local finite element - const FEType& fineBaseSet = p1FECache.get(fIt->type());; + const FEType& fineBaseSet = p1FECache.get(fIt.type());; const size_t numFineBaseFct = fineBaseSet.localBasis().size(); @@ -285,15 +259,15 @@ public: } // Get ancestor element in the coarse grid view, and the local position there. - typename GridViewFine::template Codim<0>::EntityPointer ancestor(fIt); + auto ancestor = fIt; - while (not gvCoarse.contains(ancestor) && ancestor->level() != 0 ) { + while (not gvCoarse.contains(ancestor) && ancestor.level() != 0 ) { - typename GridType::template Codim<0>::LocalGeometry geometryInFather = ancestor->geometryInFather(); + auto geometryInFather = ancestor.geometryInFather(); for (size_t i=0; i<numFineBaseFct; i++) local[i] = geometryInFather.global(local[i]); - ancestor = ancestor->father(); + ancestor = ancestor.father(); } @@ -302,7 +276,7 @@ public: for (size_t j=0; j<numFineBaseFct; j++) { // Get local finite element - const FEType& coarseBaseSet = p1FECache.get(ancestor->type()); + const FEType& coarseBaseSet = p1FECache.get(ancestor.type()); const size_t numCoarseBaseFct = coarseBaseSet.localBasis().size(); @@ -313,14 +287,14 @@ public: coarseBaseSet.localBasis().evaluateFunction(local[j], values); const Dune::LocalKey& jLocalKey = fineBaseSet.localCoefficients().localKey(j); - int globalFine = fineIndexSet.subIndex(*fIt, jLocalKey.subEntity(), jLocalKey.codim()); + int globalFine = fineIndexSet.subIndex(fIt, jLocalKey.subEntity(), jLocalKey.codim()); for (size_t i=0; i<numCoarseBaseFct; i++) { if (values[i] > 0.001) { const Dune::LocalKey& iLocalKey = coarseBaseSet.localCoefficients().localKey(i); - int globalCoarse = coarseIndexSet.subIndex(*ancestor, iLocalKey.subEntity(), iLocalKey.codim()); + int globalCoarse = coarseIndexSet.subIndex(ancestor, iLocalKey.subEntity(), iLocalKey.codim()); indices.add(globalFine, globalCoarse); } } @@ -333,14 +307,13 @@ public: // ///////////////////////////////////////////// // Compute the matrix // ///////////////////////////////////////////// - for (fIt = gvFine.template begin<0>(); fIt != fEndIt; ++fIt) { + for (const auto& fIt : elements(gvFine)) { - const Dune::ReferenceElement<ctype,dim>& fineRefElement - = Dune::ReferenceElements<ctype, dim>::general(fIt->type()); + const auto& fineRefElement = Dune::ReferenceElements<ctype, dim>::general(fIt.type()); // Get local finite element - const FEType& fineBaseSet = p1FECache.get(fIt->type());; + const FEType& fineBaseSet = p1FECache.get(fIt.type());; const size_t numFineBaseFct = fineBaseSet.localBasis().size(); @@ -353,15 +326,15 @@ public: } // Get ancestor element in the coarse grid view, and the local position there. - typename GridViewFine::template Codim<0>::EntityPointer ancestor(fIt); + auto ancestor = fIt; - while (not gvCoarse.contains(ancestor) && ancestor->level() != 0 ) { + while (not gvCoarse.contains(ancestor) && ancestor.level() != 0 ) { - typename GridType::template Codim<0>::LocalGeometry geometryInFather = ancestor->geometryInFather(); + auto geometryInFather = ancestor.geometryInFather(); for (size_t i=0; i<numFineBaseFct; i++) local[i] = geometryInFather.global(local[i]); - ancestor = ancestor->father(); + ancestor = ancestor.father(); } @@ -370,7 +343,7 @@ public: for (size_t j=0; j<numFineBaseFct; j++) { // Get local finite element - const FEType& coarseBaseSet = p1FECache.get(ancestor->type()); + const FEType& coarseBaseSet = p1FECache.get(ancestor.type()); const size_t numCoarseBaseFct = coarseBaseSet.localBasis().size(); @@ -381,22 +354,20 @@ public: coarseBaseSet.localBasis().evaluateFunction(local[j], values); const Dune::LocalKey& jLocalKey = fineBaseSet.localCoefficients().localKey(j); - int globalFine = fineIndexSet.subIndex(*fIt, jLocalKey.subEntity(), jLocalKey.codim()); + int globalFine = fineIndexSet.subIndex(fIt, jLocalKey.subEntity(), jLocalKey.codim()); for (size_t i=0; i<numCoarseBaseFct; i++) { if (values[i] > 0.001) { const Dune::LocalKey& iLocalKey = coarseBaseSet.localCoefficients().localKey(i); - int globalCoarse = coarseIndexSet.subIndex(*ancestor, iLocalKey.subEntity(), iLocalKey.codim()); + int globalCoarse = coarseIndexSet.subIndex(ancestor, iLocalKey.subEntity(), iLocalKey.codim()); matrix[globalFine][globalCoarse] = Dune::ScaledIdentityMatrix<ctype,TransferMatrixBlock::rows>(values[i]); } } } - } - } //! Multiply the vector f from the right to the prolongation matrix