diff options
author | Erik Schnetter <schnetter@aei.mpg.de> | 2005-01-01 18:22:00 +0000 |
---|---|---|
committer | Erik Schnetter <schnetter@aei.mpg.de> | 2005-01-01 18:22:00 +0000 |
commit | 2a38d7eb6a6db8b150a9f6fd5f1c844b8d0ef74a (patch) | |
tree | 1c5e763b0ffee9744ee708d6016bdbd36b8f8312 /Carpet/CarpetRegrid/src/regrid.hh | |
parent | 049cec8e042a508511fdb0f0948de63f84f9b8be (diff) |
global: Turn CarpetLib templates into classes
Turn most of the templates in CarpetLib, which used to have the form
template<int D> class XXX
into classes, i.e., into something like
class XXX
by setting D to the new global integer constant dim, which in turn is set to 3.
The templates gf and data, which used to be of the form
template<typename T, int D> class XXX
are now of the form
template<typename T> class XXX
The templates vect, bbox, and bboxset remain templates.
This change simplifies the code somewhat.
darcs-hash:20050101182234-891bb-c3063528841f0d078b12cc506309ea27d8ce730d.gz
Diffstat (limited to 'Carpet/CarpetRegrid/src/regrid.hh')
-rw-r--r-- | Carpet/CarpetRegrid/src/regrid.hh | 80 |
1 files changed, 40 insertions, 40 deletions
diff --git a/Carpet/CarpetRegrid/src/regrid.hh b/Carpet/CarpetRegrid/src/regrid.hh index 460290646..fea47e073 100644 --- a/Carpet/CarpetRegrid/src/regrid.hh +++ b/Carpet/CarpetRegrid/src/regrid.hh @@ -43,25 +43,25 @@ namespace CarpetRegrid { int BaseLevel (cGH const * const cctkGH, - gh<dim> const & hh, - gh<dim>::rexts & bbsss, - gh<dim>::rbnds & obss, - gh<dim>::rprocs & pss); + gh const & hh, + gh::rexts & bbsss, + gh::rbnds & obss, + gh::rprocs & pss); int Centre (cGH const * const cctkGH, - gh<dim> const & hh, - gh<dim>::rexts & bbsss, - gh<dim>::rbnds & obss, - gh<dim>::rprocs & pss); + gh const & hh, + gh::rexts & bbsss, + gh::rbnds & obss, + gh::rprocs & pss); int ManualGridpoints (cGH const * const cctkGH, - gh<dim> const & hh, - gh<dim>::rexts & bbsss, - gh<dim>::rbnds & obss, - gh<dim>::rprocs & pss); + gh const & hh, + gh::rexts & bbsss, + gh::rbnds & obss, + gh::rprocs & pss); void ManualGridpoints_OneLevel (const cGH * const cctkGH, - const gh<dim> & hh, + const gh & hh, const int rl, const int numrl, const ivect ilower, @@ -71,13 +71,13 @@ namespace CarpetRegrid { vector<bbvect> & obs); int ManualCoordinates (cGH const * const cctkGH, - gh<dim> const & hh, - gh<dim>::rexts & bbsss, - gh<dim>::rbnds & obss, - gh<dim>::rprocs & pss); + gh const & hh, + gh::rexts & bbsss, + gh::rbnds & obss, + gh::rprocs & pss); void ManualCoordinates_OneLevel (const cGH * const cctkGH, - const gh<dim> & hh, + const gh & hh, const int rl, const int numrl, const rvect lower, @@ -87,55 +87,55 @@ namespace CarpetRegrid { vector<bbvect> & obs); ivect delta2int (const cGH * const cctkGH, - const gh<dim>& hh, + const gh& hh, const rvect & rpos, const int rl); ivect pos2int (const cGH* const cctkGH, - const gh<dim>& hh, + const gh& hh, const rvect & rpos, const int rl); int ManualGridpointList (cGH const * const cctkGH, - gh<dim> const & hh, - gh<dim>::rexts & bbsss, - gh<dim>::rbnds & obss, - gh<dim>::rprocs & pss); + gh const & hh, + gh::rexts & bbsss, + gh::rbnds & obss, + gh::rprocs & pss); int ManualCoordinateList (cGH const * const cctkGH, - gh<dim> const & hh, - gh<dim>::rexts & bbsss, - gh<dim>::rbnds & obss, - gh<dim>::rprocs & pss); + gh const & hh, + gh::rexts & bbsss, + gh::rbnds & obss, + gh::rprocs & pss); int Moving (cGH const * const cctkGH, - gh<dim> const & hh, - gh<dim>::rexts & bbsss, - gh<dim>::rbnds & obss, - gh<dim>::rprocs & pss); + gh const & hh, + gh::rexts & bbsss, + gh::rbnds & obss, + gh::rprocs & pss); int Automatic (cGH const * const cctkGH, - gh<dim> const & hh, - gh<dim>::rexts & bbsss, - gh<dim>::rbnds & obss, - gh<dim>::rprocs & pss); + gh const & hh, + gh::rexts & bbsss, + gh::rbnds & obss, + gh::rprocs & pss); void Automatic_OneLevel (const cGH * const cctkGH, - const gh<dim> & hh, + const gh & hh, const int rl, const int numrl, const int minwidth, const CCTK_REAL minfraction, const CCTK_REAL maxerror, - const gf<CCTK_REAL,dim> & errorvar, + const gf<CCTK_REAL> & errorvar, vector<ibbox> & bbs, vector<bbvect> & obs); void Automatic_Recursive (const cGH * const cctkGH, - const gh<dim> & hh, + const gh & hh, const int minwidth, const CCTK_REAL minfraction, const CCTK_REAL maxerror, - const data<CCTK_REAL,dim> & errorvar, + const data<CCTK_REAL> & errorvar, list<ibbox> & bbl, const ibbox & region); |