diff options
Adapt refinement to meta-descriptor
-rw-r--r-- | src/refinement/coupler2D.cpp | 7 | ||||
-rw-r--r-- | src/refinement/coupler2D.h | 12 | ||||
-rw-r--r-- | src/refinement/coupler2D.hh | 64 | ||||
-rw-r--r-- | src/refinement/grid2D.cpp | 5 | ||||
-rw-r--r-- | src/refinement/grid2D.h | 13 | ||||
-rw-r--r-- | src/refinement/grid2D.hh | 56 |
6 files changed, 80 insertions, 77 deletions
diff --git a/src/refinement/coupler2D.cpp b/src/refinement/coupler2D.cpp index 37997f5..733d7c9 100644 --- a/src/refinement/coupler2D.cpp +++ b/src/refinement/coupler2D.cpp @@ -25,12 +25,11 @@ #include "coupler2D.hh" #include "dynamics/latticeDescriptors.h" -#include "dynamics/latticeDescriptors.hh" namespace olb { -template class Coupler2D<double,descriptors::D2Q9Descriptor>; -template class FineCoupler2D<double,descriptors::D2Q9Descriptor>; -template class CoarseCoupler2D<double,descriptors::D2Q9Descriptor>; +template class Coupler2D<double,descriptors::D2Q9<>>; +template class FineCoupler2D<double,descriptors::D2Q9<>>; +template class CoarseCoupler2D<double,descriptors::D2Q9<>>; } diff --git a/src/refinement/coupler2D.h b/src/refinement/coupler2D.h index d43a294..9b310b7 100644 --- a/src/refinement/coupler2D.h +++ b/src/refinement/coupler2D.h @@ -29,7 +29,7 @@ namespace olb { -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> class Coupler2D { protected: Grid2D<T,DESCRIPTOR>& _coarse; @@ -57,12 +57,12 @@ public: }; -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> class FineCoupler2D : public Coupler2D<T,DESCRIPTOR> { private: - std::vector<Vector<T,1>> _c2f_rho; - std::vector<Vector<T,DESCRIPTOR<T>::d>> _c2f_u; - std::vector<Vector<T,DESCRIPTOR<T>::q>> _c2f_fneq; + std::vector<Vector<T,1>> _c2f_rho; + std::vector<Vector<T,DESCRIPTOR::d>> _c2f_u; + std::vector<Vector<T,DESCRIPTOR::q>> _c2f_fneq; public: FineCoupler2D(Grid2D<T,DESCRIPTOR>& coarse, Grid2D<T,DESCRIPTOR>& fine, @@ -74,7 +74,7 @@ public: }; -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> class CoarseCoupler2D : public Coupler2D<T,DESCRIPTOR> { public: CoarseCoupler2D(Grid2D<T,DESCRIPTOR>& coarse, Grid2D<T,DESCRIPTOR>& fine, diff --git a/src/refinement/coupler2D.hh b/src/refinement/coupler2D.hh index c2b27db..800986d 100644 --- a/src/refinement/coupler2D.hh +++ b/src/refinement/coupler2D.hh @@ -30,32 +30,32 @@ namespace olb { -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> T Coupler2D<T,DESCRIPTOR>::getScalingFactor() const { const T coarseTau = _coarse.getConverter().getLatticeRelaxationTime(); return (coarseTau - 0.25) / coarseTau; } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> T Coupler2D<T,DESCRIPTOR>::getInvScalingFactor() const { return 1./getScalingFactor(); } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> const Vector<int,3>& Coupler2D<T,DESCRIPTOR>::getFineLatticeR(int y) const { return _fineLatticeR[y]; } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> const Vector<int,3>& Coupler2D<T,DESCRIPTOR>::getCoarseLatticeR(int y) const { return _coarseLatticeR[y]; } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> Coupler2D<T,DESCRIPTOR>::Coupler2D(Grid2D<T,DESCRIPTOR>& coarse, Grid2D<T,DESCRIPTOR>& fine, Vector<T,2> origin, Vector<T,2> extend): _coarse(coarse), @@ -86,13 +86,13 @@ Coupler2D<T,DESCRIPTOR>::Coupler2D(Grid2D<T,DESCRIPTOR>& coarse, Grid2D<T,DESCRI } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> FineCoupler2D<T,DESCRIPTOR>::FineCoupler2D(Grid2D<T,DESCRIPTOR>& coarse, Grid2D<T,DESCRIPTOR>& fine, Vector<T,2> origin, Vector<T,2> extend): Coupler2D<T,DESCRIPTOR>(coarse, fine, origin, extend), _c2f_rho(this->_coarseSize), _c2f_u(this->_coarseSize, Vector<T,2>(T{})), - _c2f_fneq(this->_coarseSize, Vector<T,DESCRIPTOR<T>::q>(T{})) + _c2f_fneq(this->_coarseSize, Vector<T,DESCRIPTOR::q>(T{})) { OstreamManager clout(std::cout,"C2F"); @@ -104,7 +104,7 @@ FineCoupler2D<T,DESCRIPTOR>::FineCoupler2D(Grid2D<T,DESCRIPTOR>& coarse, Grid2D< clout << "fine size: " << this->_fineSize << std::endl; } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> void FineCoupler2D<T,DESCRIPTOR>::store() { auto& coarseLattice = this->_coarse.getSuperLattice(); @@ -114,7 +114,7 @@ void FineCoupler2D<T,DESCRIPTOR>::store() const auto pos = this->getCoarseLatticeR(y); T rho{}; T u[2] {}; - T fNeq[DESCRIPTOR<T>::q] {}; + T fNeq[DESCRIPTOR::q] {}; Cell<T,DESCRIPTOR> coarseCell; coarseLattice.get(pos, coarseCell); lbHelpers<T,DESCRIPTOR>::computeRhoU(coarseCell, rho, u); @@ -122,7 +122,7 @@ void FineCoupler2D<T,DESCRIPTOR>::store() _c2f_rho[y] = Vector<T,1>(rho); _c2f_u[y] = Vector<T,2>(u); - _c2f_fneq[y] = Vector<T,DESCRIPTOR<T>::q>(fNeq); + _c2f_fneq[y] = Vector<T,DESCRIPTOR::q>(fNeq); } } @@ -156,7 +156,7 @@ Vector<T,N> order4interpolation(const std::vector<Vector<T,N>>& data, int y) } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> void FineCoupler2D<T,DESCRIPTOR>::interpolate() { auto& coarseLattice = this->_coarse.getSuperLattice(); @@ -173,14 +173,14 @@ void FineCoupler2D<T,DESCRIPTOR>::interpolate() _c2f_rho[y] = order2interpolation(Vector<T,1>(rho), _c2f_rho[y]); _c2f_u[y] = order2interpolation(Vector<T,2>(u), _c2f_u[y]); - T fNeq[DESCRIPTOR<T>::q] {}; + T fNeq[DESCRIPTOR::q] {}; lbHelpers<T,DESCRIPTOR>::computeFneq(coarseCell, fNeq, rho, u); - _c2f_fneq[y] = order2interpolation(Vector<T,DESCRIPTOR<T>::q>(fNeq), _c2f_fneq[y]); + _c2f_fneq[y] = order2interpolation(Vector<T,DESCRIPTOR::q>(fNeq), _c2f_fneq[y]); } } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> void FineCoupler2D<T,DESCRIPTOR>::couple() { const auto& coarseLattice = this->_coarse.getSuperLattice(); @@ -191,14 +191,14 @@ void FineCoupler2D<T,DESCRIPTOR>::couple() const auto& coarsePos = this->getCoarseLatticeR(y); const auto& finePos = this->getFineLatticeR(2*y); - T fEq[DESCRIPTOR<T>::q] {}; + T fEq[DESCRIPTOR::q] {}; Cell<T,DESCRIPTOR> coarseCell; coarseLattice.get(coarsePos, coarseCell); lbHelpers<T,DESCRIPTOR>::computeFeq(coarseCell, fEq); Cell<T,DESCRIPTOR> cell; fineLattice.get(finePos, cell); - for (int iPop=0; iPop < DESCRIPTOR<T>::q; ++iPop) { + for (int iPop=0; iPop < DESCRIPTOR::q; ++iPop) { cell[iPop] = fEq[iPop] + this->getScalingFactor() * _c2f_fneq[y][iPop]; } fineLattice.set(finePos, cell); @@ -216,7 +216,7 @@ void FineCoupler2D<T,DESCRIPTOR>::couple() Cell<T,DESCRIPTOR> fineCell; fineLattice.get(finePos, fineCell); - for (int iPop=0; iPop < DESCRIPTOR<T>::q; ++iPop) { + for (int iPop=0; iPop < DESCRIPTOR::q; ++iPop) { fineCell[iPop] = lbHelpers<T,DESCRIPTOR>::equilibrium(iPop, rho[0], u.data, uSqr) + this->getScalingFactor() * fneq[iPop]; } @@ -235,7 +235,7 @@ void FineCoupler2D<T,DESCRIPTOR>::couple() Cell<T,DESCRIPTOR> fineCell; fineLattice.get(finePos, fineCell); - for (int iPop=0; iPop < DESCRIPTOR<T>::q; ++iPop) { + for (int iPop=0; iPop < DESCRIPTOR::q; ++iPop) { fineCell[iPop] = lbHelpers<T,DESCRIPTOR>::equilibrium(iPop, rho[0], u.data, uSqr) + this->getScalingFactor() * fneq[iPop]; } @@ -254,7 +254,7 @@ void FineCoupler2D<T,DESCRIPTOR>::couple() Cell<T,DESCRIPTOR> fineCell; fineLattice.get(finePos, fineCell); - for (int iPop=0; iPop < DESCRIPTOR<T>::q; ++iPop) { + for (int iPop=0; iPop < DESCRIPTOR::q; ++iPop) { fineCell[iPop] = lbHelpers<T,DESCRIPTOR>::equilibrium(iPop, rho[0], u.data, uSqr) + this->getScalingFactor() * fneq[iPop]; } @@ -264,30 +264,30 @@ void FineCoupler2D<T,DESCRIPTOR>::couple() } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> void computeRestrictedFneq(const SuperLattice2D<T,DESCRIPTOR>& lattice, Vector<int,3> latticeR, - T restrictedFneq[DESCRIPTOR<T>::q]) + T restrictedFneq[DESCRIPTOR::q]) { - for (int iPop=0; iPop < DESCRIPTOR<T>::q; ++iPop) { - const auto neighbor = latticeR + Vector<int,3> {0, DESCRIPTOR<T>::c[iPop][0], DESCRIPTOR<T>::c[iPop][1]}; + for (int iPop=0; iPop < DESCRIPTOR::q; ++iPop) { + const auto neighbor = latticeR + Vector<int,3> {0, descriptors::c<DESCRIPTOR>(iPop,0), descriptors::c<DESCRIPTOR>(iPop,1)}; Cell<T,DESCRIPTOR> cell; lattice.get(neighbor, cell); - T fNeq[DESCRIPTOR<T>::q] {}; + T fNeq[DESCRIPTOR::q] {}; lbHelpers<T,DESCRIPTOR>::computeFneq(cell, fNeq); - for (int jPop=0; jPop < DESCRIPTOR<T>::q; ++jPop) { + for (int jPop=0; jPop < DESCRIPTOR::q; ++jPop) { restrictedFneq[jPop] += fNeq[jPop]; } } - for (int iPop=0; iPop < DESCRIPTOR<T>::q; ++iPop) { - restrictedFneq[iPop] /= DESCRIPTOR<T>::q; + for (int iPop=0; iPop < DESCRIPTOR::q; ++iPop) { + restrictedFneq[iPop] /= DESCRIPTOR::q; } } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> CoarseCoupler2D<T,DESCRIPTOR>::CoarseCoupler2D( Grid2D<T,DESCRIPTOR>& coarse, Grid2D<T,DESCRIPTOR>& fine, Vector<T,2> origin, Vector<T,2> extend): @@ -303,7 +303,7 @@ CoarseCoupler2D<T,DESCRIPTOR>::CoarseCoupler2D( clout << "coarse size: " << this->_coarseSize << std::endl; } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> void CoarseCoupler2D<T,DESCRIPTOR>::couple() { const auto& fineLattice = this->_fine.getSuperLattice(); @@ -314,18 +314,18 @@ void CoarseCoupler2D<T,DESCRIPTOR>::couple() const auto& finePos = this->getFineLatticeR(2*y); const auto& coarsePos = this->getCoarseLatticeR(y); - T fEq[DESCRIPTOR<T>::q] {}; + T fEq[DESCRIPTOR::q] {}; Cell<T,DESCRIPTOR> fineCell; fineLattice.get(finePos, fineCell); lbHelpers<T,DESCRIPTOR>::computeFeq(fineCell, fEq); - T fNeq[DESCRIPTOR<T>::q] {}; + T fNeq[DESCRIPTOR::q] {}; computeRestrictedFneq(fineLattice, finePos, fNeq); Cell<T,DESCRIPTOR> coarseCell; coarseLattice.get(coarsePos, coarseCell); - for (int iPop=0; iPop < DESCRIPTOR<T>::q; ++iPop) { + for (int iPop=0; iPop < DESCRIPTOR::q; ++iPop) { coarseCell[iPop] = fEq[iPop] + this->getInvScalingFactor() * fNeq[iPop]; } diff --git a/src/refinement/grid2D.cpp b/src/refinement/grid2D.cpp index 86ca284..4472fc7 100644 --- a/src/refinement/grid2D.cpp +++ b/src/refinement/grid2D.cpp @@ -25,11 +25,10 @@ #include "grid2D.hh" #include "dynamics/latticeDescriptors.h" -#include "dynamics/latticeDescriptors.hh" namespace olb { -template class Grid2D<double,descriptors::D2Q9Descriptor>; -template class RefiningGrid2D<double,descriptors::D2Q9Descriptor>; +template class Grid2D<double,descriptors::D2Q9<>>; +template class RefiningGrid2D<double,descriptors::D2Q9<>>; } diff --git a/src/refinement/grid2D.h b/src/refinement/grid2D.h index 42937ca..51d19ce 100644 --- a/src/refinement/grid2D.h +++ b/src/refinement/grid2D.h @@ -42,9 +42,12 @@ namespace olb { -template <typename T, template<typename> class DESCRIPTOR> class FineCoupler2D; -template <typename T, template<typename> class DESCRIPTOR> class CoarseCoupler2D; -template <typename T, template<typename> class DESCRIPTOR> class RefiningGrid2D; +template <typename T, typename DESCRIPTOR> class FineCoupler2D; +template <typename T, typename DESCRIPTOR> class CoarseCoupler2D; +template <typename T, typename DESCRIPTOR> class RefiningGrid2D; + +template <typename T, typename DESCRIPTOR> class sOnLatticeBoundaryCondition2D; +template <typename T, typename DESCRIPTOR> class sOffLatticeBoundaryCondition2D; template <typename T> using RelaxationTime = utilities::NamedType<T,struct NamedRelaxationTime>; @@ -70,7 +73,7 @@ struct Characteristics { }; -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> class Grid2D { protected: FunctorPtr<IndicatorF2D<T>> _domainF; @@ -141,7 +144,7 @@ public: }; -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> class RefiningGrid2D : public Grid2D<T,DESCRIPTOR> { private: const Vector<T,2> _origin; diff --git a/src/refinement/grid2D.hh b/src/refinement/grid2D.hh index 4f76896..d4310e5 100644 --- a/src/refinement/grid2D.hh +++ b/src/refinement/grid2D.hh @@ -28,12 +28,14 @@ #include "coupler2D.hh" #include "utilities/vectorHelpers.h" +#include "boundary/superBoundaryCondition2D.h" +#include "boundary/superOffBoundaryCondition2D.h" #include "communication/heuristicLoadBalancer.h" namespace olb { -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> Grid2D<T,DESCRIPTOR>::Grid2D(FunctorPtr<IndicatorF2D<T>>&& domainF, RelaxationTime<T> tau, int resolution, @@ -68,7 +70,7 @@ Grid2D<T,DESCRIPTOR>::Grid2D(FunctorPtr<IndicatorF2D<T>>&& domainF, _converter->print(); } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> Grid2D<T,DESCRIPTOR>::Grid2D(FunctorPtr<IndicatorF2D<T>>&& domainF, LatticeVelocity<T> latticeVelocity, int resolution, @@ -103,7 +105,7 @@ Grid2D<T,DESCRIPTOR>::Grid2D(FunctorPtr<IndicatorF2D<T>>&& domainF, _converter->print(); } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> Grid2D<T,DESCRIPTOR>::Grid2D( FunctorPtr<IndicatorF2D<T>>&& domainF, RelaxationTime<T> tau, @@ -115,7 +117,7 @@ Grid2D<T,DESCRIPTOR>::Grid2D( Characteristics<T>(re)) { } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> Grid2D<T,DESCRIPTOR>::Grid2D( FunctorPtr<IndicatorF2D<T>>&& domainF, LatticeVelocity<T> latticeVelocity, @@ -127,43 +129,43 @@ Grid2D<T,DESCRIPTOR>::Grid2D( Characteristics<T>(re)) { } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> Characteristics<T> Grid2D<T,DESCRIPTOR>::getCharacteristics() const { return _characteristics; } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> UnitConverter<T,DESCRIPTOR>& Grid2D<T,DESCRIPTOR>::getConverter() { return *_converter; } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> CuboidGeometry2D<T>& Grid2D<T,DESCRIPTOR>::getCuboidGeometry() { return *_cuboids; } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> LoadBalancer<T>& Grid2D<T,DESCRIPTOR>::getLoadBalancer() { return *_balancer; } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> SuperGeometry2D<T>& Grid2D<T,DESCRIPTOR>::getSuperGeometry() { return *_geometry; } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> SuperLattice2D<T,DESCRIPTOR>& Grid2D<T,DESCRIPTOR>::getSuperLattice() { return *_lattice; } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> Dynamics<T,DESCRIPTOR>& Grid2D<T,DESCRIPTOR>::addDynamics( std::unique_ptr<Dynamics<T,DESCRIPTOR>>&& dynamics) { @@ -172,7 +174,7 @@ Dynamics<T,DESCRIPTOR>& Grid2D<T,DESCRIPTOR>::addDynamics( return ref; } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> sOnLatticeBoundaryCondition2D<T,DESCRIPTOR>& Grid2D<T,DESCRIPTOR>::getOnLatticeBoundaryCondition() { _onLatticeBoundaryConditions.emplace_back( @@ -180,7 +182,7 @@ sOnLatticeBoundaryCondition2D<T,DESCRIPTOR>& Grid2D<T,DESCRIPTOR>::getOnLatticeB return *_onLatticeBoundaryConditions.back(); } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> sOffLatticeBoundaryCondition2D<T,DESCRIPTOR>& Grid2D<T,DESCRIPTOR>::getOffLatticeBoundaryCondition() { _offLatticeBoundaryConditions.emplace_back( @@ -188,7 +190,7 @@ sOffLatticeBoundaryCondition2D<T,DESCRIPTOR>& Grid2D<T,DESCRIPTOR>::getOffLattic return *_offLatticeBoundaryConditions.back(); } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> void Grid2D<T,DESCRIPTOR>::collideAndStream() { for ( auto& fineCoupler : _fineCouplers ) { @@ -220,7 +222,7 @@ void Grid2D<T,DESCRIPTOR>::collideAndStream() } } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> FineCoupler2D<T,DESCRIPTOR>& Grid2D<T,DESCRIPTOR>::addFineCoupling( Grid2D<T,DESCRIPTOR>& fineGrid, Vector<T,2> origin, Vector<T,2> extend) { @@ -230,7 +232,7 @@ FineCoupler2D<T,DESCRIPTOR>& Grid2D<T,DESCRIPTOR>::addFineCoupling( return *_fineCouplers.back(); } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> CoarseCoupler2D<T,DESCRIPTOR>& Grid2D<T,DESCRIPTOR>::addCoarseCoupling( Grid2D<T,DESCRIPTOR>& fineGrid, Vector<T,2> origin, Vector<T,2> extend) { @@ -240,7 +242,7 @@ CoarseCoupler2D<T,DESCRIPTOR>& Grid2D<T,DESCRIPTOR>::addCoarseCoupling( return *_coarseCouplers.back(); } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> Vector<T,2> Grid2D<T,DESCRIPTOR>::alignOriginToGrid(Vector<T,2> physR) const { Vector<int,3> latticeR{}; @@ -248,7 +250,7 @@ Vector<T,2> Grid2D<T,DESCRIPTOR>::alignOriginToGrid(Vector<T,2> physR) const return _cuboids->getPhysR(latticeR.toStdVector()); } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> Vector<T,2> Grid2D<T,DESCRIPTOR>::alignExtendToGrid(Vector<T,2> extend) const { const T deltaX = _converter->getPhysDeltaX(); @@ -258,7 +260,7 @@ Vector<T,2> Grid2D<T,DESCRIPTOR>::alignExtendToGrid(Vector<T,2> extend) const }; } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> void Grid2D<T,DESCRIPTOR>::forEachGrid(std::function<void(Grid2D<T,DESCRIPTOR>&)>&& f) { f(*this); @@ -267,7 +269,7 @@ void Grid2D<T,DESCRIPTOR>::forEachGrid(std::function<void(Grid2D<T,DESCRIPTOR>&) } } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> void Grid2D<T,DESCRIPTOR>::forEachGrid( const std::string& id, std::function<void(Grid2D<T,DESCRIPTOR>&,const std::string&)>&& f) @@ -280,7 +282,7 @@ void Grid2D<T,DESCRIPTOR>::forEachGrid( } } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> Grid2D<T,DESCRIPTOR>& Grid2D<T,DESCRIPTOR>::locate(Vector<T,2> pos) { int iC; @@ -292,7 +294,7 @@ Grid2D<T,DESCRIPTOR>& Grid2D<T,DESCRIPTOR>::locate(Vector<T,2> pos) return *this; } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> std::size_t Grid2D<T,DESCRIPTOR>::getActiveVoxelN() const { std::size_t n = _geometry->getStatistics().getNvoxel(); @@ -302,7 +304,7 @@ std::size_t Grid2D<T,DESCRIPTOR>::getActiveVoxelN() const return n; } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> RefiningGrid2D<T,DESCRIPTOR>& Grid2D<T,DESCRIPTOR>::refine( Vector<T,2> wantedOrigin, Vector<T,2> wantedExtend, bool addCouplers) { @@ -348,7 +350,7 @@ RefiningGrid2D<T,DESCRIPTOR>& Grid2D<T,DESCRIPTOR>::refine( } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> RefiningGrid2D<T,DESCRIPTOR>::RefiningGrid2D( Grid2D<T,DESCRIPTOR>& parentGrid, Vector<T,2> origin, Vector<T,2> extend): Grid2D<T,DESCRIPTOR>( @@ -360,19 +362,19 @@ RefiningGrid2D<T,DESCRIPTOR>::RefiningGrid2D( _extend(extend), _parentGrid(parentGrid) { } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> Vector<T,2> RefiningGrid2D<T,DESCRIPTOR>::getOrigin() const { return _origin; } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> Vector<T,2> RefiningGrid2D<T,DESCRIPTOR>::getExtend() const { return _extend; } -template <typename T, template<typename> class DESCRIPTOR> +template <typename T, typename DESCRIPTOR> std::unique_ptr<IndicatorF2D<T>> RefiningGrid2D<T,DESCRIPTOR>::getRefinedOverlap() const { const T coarseDeltaX = _parentGrid.getConverter().getPhysDeltaX(); |