aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetRegrid/src/regrid.hh
diff options
context:
space:
mode:
authorschnetter <>2004-01-25 13:57:00 +0000
committerschnetter <>2004-01-25 13:57:00 +0000
commitcfcb3f61bcd9a88abbb97d1c49325cc67dbe70a8 (patch)
tree7f2c49463717ab281fec6bb139c91533466e25dd /Carpet/CarpetRegrid/src/regrid.hh
parent55f76be95c2272b35b6941d6ec38a77bfb23a101 (diff)
Import the recently announced changes:
Import the recently announced changes: 1. Carpet has now an infrastructure for multiple maps (aka "grid patches"). Instead of a single grid hierarchy there can now be several. This is largely untested, because the remainder of Cactus cannot handle multiple coordinate systems. 2. The order in which the schedule bins are called has changed. As Ian Hawke pointed out, the previous order during time evolution was inconsistent. The initial data ordering did not allow for recovering and was not usable for progressively solving elliptic equations for initial data. 3. Carpet now supports convergence levels. The convergence level specifies by how many factors of two the resolution in the parameter file should be coarsened (or refined, if negative). This should make convergence tests and test runs much easier. It is, in principle, also possible to run several convergence levels at once. This has not been tested because the remainder of Cactus cannot handle multiple resolutions. This will be necessary for a multigrid solver, and also for having a shadow hierarchy to determine where to refine adaptively. 4. Carpet works together with the new CoordBase domain specification parameters. Without these, using convergence levels will lead to very strange results. 5. The "modes" have changed. There are now: meta mode: the whole simulation global mode: one convergence level level mode: one refinement level singlemap mode: one map on one refinement level local mode: as previously The whole mode handling has been cleaned up. 6. The regridding thorn has been cleaned up. 7. The kind of prolongation stencil is now determined in Carpet, i.e. at a fairly hight level, instead of in CarpetLib. 8. The low-order prolongation operators have been made much more efficient (as have previously the higher-order ones). 9. Assorted smaller changes. For Carpet users, there should be no major incompatibilities. The major improvements are 3 and 4 combined. Here is an example: CoordBase::domainsize = extent CoordBase::spacing = gridspacing CoordBase::zero_origin_x = yes CoordBase::zero_origin_y = yes CoordBase::zero_origin_z = yes CoordBase::xextent = 20.0 CoordBase::yextent = 20.0 CoordBase::zextent = 20.0 CoordBase::dx = 1.0 CoordBase::dy = 1.0 CoordBase::dz = 1.0 CoordBase::boundary_shiftout_x_lower = 1 CoordBase::boundary_shiftout_y_lower = 1 CoordBase::boundary_shiftout_z_lower = 1 Carpet::domain_from_coordbase = yes Carpet::convergence_level = 0 grid::type = coordbase grid::domain = octant grid::avoid_origin = no This gives you a grid that extends from the origin ("zero_origin") up to 20.0 with a grid spacing of 1.0. Symmetry zones and boundary zones are added automatically. The "shiftout" says that there is no boundary point on the origin. The staggering parameters (not shown) default to "no". In order to change the resolution, only the convergence level has to be adjusted. Note that the old way of specifying the domain extent still works. For Carpet developers, one major change is the new mode handling. As described in 5, the looping macros (that loop over all refinement levels, or all components) have changed. darcs-hash:20040125135727-07bb3-51c9647c1b5080e7e180b52a1b81fa155cfd19e9.gz
Diffstat (limited to 'Carpet/CarpetRegrid/src/regrid.hh')
-rw-r--r--Carpet/CarpetRegrid/src/regrid.hh219
1 files changed, 169 insertions, 50 deletions
diff --git a/Carpet/CarpetRegrid/src/regrid.hh b/Carpet/CarpetRegrid/src/regrid.hh
index 45a390a97..27edc99c9 100644
--- a/Carpet/CarpetRegrid/src/regrid.hh
+++ b/Carpet/CarpetRegrid/src/regrid.hh
@@ -1,4 +1,4 @@
-// $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetRegrid/src/regrid.hh,v 1.10 2003/11/13 16:04:37 schnetter Exp $
+// $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetRegrid/src/regrid.hh,v 1.11 2004/01/25 14:57:30 schnetter Exp $
#ifndef CARPETREGRID_HH
#define CARPETREGRID_HH
@@ -6,14 +6,15 @@
#include <list>
#include "cctk.h"
+#include "cctk_Arguments.h"
+#include "bbox.hh"
#include "gf.hh"
#include "gh.hh"
+#include "vect.hh"
#include "carpet.hh"
-#include "regrid.h"
-
namespace CarpetRegrid {
@@ -23,53 +24,171 @@ namespace CarpetRegrid {
- typedef vect<int,dim> ivect;
- typedef bbox<int,dim> ibbox;
-
- typedef vect<vect<bool,2>,dim> bvect;
-
- typedef vect<CCTK_REAL,dim> rvect;
- typedef bbox<CCTK_REAL,dim> rbbox;
-
-
-
- int CarpetRegridRegrid (const cGH * const cctkGH,
- gh<dim>::rexts& bbsss,
- gh<dim>::rbnds& obss,
- gh<dim>::rprocs& pss);
-
-
-
- void MakeRegions_BaseLevel (const cGH* cctkGH,
- list<ibbox>& bbl, list<bvect>& obl);
-
- void MakeRegions_RefineCentre (const cGH* cctkGH, const int reflevels,
- list<ibbox>& bbl, list<bvect>& obl);
-
- void MakeRegions_AsSpecified (const cGH* cctkGH, const int reflevels,
- const vector<ivect> lower,
- const vector<ivect> upper,
- list<ibbox>& bbl, list<bvect>& obl);
- void MakeRegions_AsSpecified (const cGH* cctkGH, const int reflevels,
- const vector<rvect> lower,
- const vector<rvect> upper,
- list<ibbox>& bbl, list<bvect>& obl);
-
- void MakeRegions_AsSpecified (const cGH* cctkGH, const int reflevels,
- const vector<vector<ibbox> > bbss,
- const vector<vector<bvect> > obss,
- list<ibbox>& bbl, list<bvect>& obl);
- void MakeRegions_AsSpecified (const cGH* cctkGH, const int reflevels,
- const vector<vector<rbbox> > bbss,
- const vector<vector<bvect> > obss,
- list<ibbox>& bbl, list<bvect>& obl);
-
- void MakeRegions_Adaptively (const cGH* cctkGH,
- const int minwidth,
- const CCTK_REAL minfraction,
- const CCTK_REAL maxerror,
- const gf<CCTK_REAL,dim>& error,
- list<ibbox>& bbl, list<bvect>& obl);
+ extern "C" {
+
+ /* Scheduled functions */
+ int CarpetRegridParamcheck (CCTK_ARGUMENTS);
+
+ /* Aliased functions */
+// CCTK_INT CarpetRegrid_Regrid (const cGH * const cctkGH,
+// gh<dim>::rexts * bbsss,
+// gh<dim>::rbnds * obss,
+// gh<dim>::rprocs * pss);
+ CCTK_INT CarpetRegrid_Regrid (CCTK_POINTER_TO_CONST const cctkGH_,
+ CCTK_INT const reflevel,
+ CCTK_INT const map,
+ CCTK_INT const size,
+ CCTK_INT const * const nboundaryzones,
+ CCTK_INT const * const is_internal,
+ CCTK_INT const * const is_staggered,
+ CCTK_INT const * const shiftout,
+ CCTK_POINTER const bbsss_,
+ CCTK_POINTER const obss_,
+ CCTK_POINTER const pss_);
+ }
+
+
+
+ int BaseLevel (cGH const * const cctkGH,
+ gh<dim> const & hh,
+ int const reflevel,
+ int const map,
+ int const size,
+ jjvect const & nboundaryzones,
+ jjvect const & is_internal,
+ jjvect const & is_staggered,
+ jjvect const & shiftout,
+ gh<dim>::rexts & bbsss,
+ gh<dim>::rbnds & obss,
+ gh<dim>::rprocs & pss);
+
+ int Centre (cGH const * const cctkGH,
+ gh<dim> const & hh,
+ int const reflevel,
+ int const map,
+ int const size,
+ jjvect const & nboundaryzones,
+ jjvect const & is_internal,
+ jjvect const & is_staggered,
+ jjvect const & shiftout,
+ gh<dim>::rexts & bbsss,
+ gh<dim>::rbnds & obss,
+ gh<dim>::rprocs & pss);
+
+ int ManualGridpoints (cGH const * const cctkGH,
+ gh<dim> const & hh,
+ int const reflevel,
+ int const map,
+ int const size,
+ jjvect const & nboundaryzones,
+ jjvect const & is_internal,
+ jjvect const & is_staggered,
+ jjvect const & shiftout,
+ gh<dim>::rexts & bbsss,
+ gh<dim>::rbnds & obss,
+ gh<dim>::rprocs & pss);
+
+ void ManualGridpoints_OneLevel (const cGH * const cctkGH,
+ const gh<dim> & hh,
+ const int reflevel,
+ const int reflevels,
+ const ivect ilower,
+ const ivect iupper,
+ const bbvect obound,
+ vector<ibbox> & bbs,
+ vector<bbvect> & obs);
+
+ int ManualCoordinates (cGH const * const cctkGH,
+ gh<dim> const & hh,
+ int const reflevel,
+ int const map,
+ int const size,
+ jjvect const & nboundaryzones,
+ jjvect const & is_internal,
+ jjvect const & is_staggered,
+ jjvect const & shiftout,
+ gh<dim>::rexts & bbsss,
+ gh<dim>::rbnds & obss,
+ gh<dim>::rprocs & pss);
+
+ void ManualCoordinates_OneLevel (const cGH * const cctkGH,
+ const gh<dim> & hh,
+ const int reflevel,
+ const int reflevels,
+ const rvect ilower,
+ const rvect iupper,
+ const bbvect obound,
+ vector<ibbox> & bbs,
+ vector<bbvect> & obs);
+
+ ivect delta2int (const cGH * const cctkGH, const gh<dim>& hh,
+ const rvect & rpos, const int reflevel);
+ ivect pos2int (const cGH* const cctkGH, const gh<dim>& hh,
+ const rvect & rpos, const int reflevel);
+
+ int ManualGridpointList (cGH const * const cctkGH,
+ gh<dim> const & hh,
+ int const reflevel,
+ int const map,
+ int const size,
+ jjvect const & nboundaryzones,
+ jjvect const & is_internal,
+ jjvect const & is_staggered,
+ jjvect const & shiftout,
+ gh<dim>::rexts & bbsss,
+ gh<dim>::rbnds & obss,
+ gh<dim>::rprocs & pss);
+
+ int ManualCoordinateList (cGH const * const cctkGH,
+ gh<dim> const & hh,
+ int const reflevel,
+ int const map,
+ int const size,
+ jjvect const & nboundaryzones,
+ jjvect const & is_internal,
+ jjvect const & is_staggered,
+ jjvect const & shiftout,
+ gh<dim>::rexts & bbsss,
+ gh<dim>::rbnds & obss,
+ gh<dim>::rprocs & pss);
+
+ int Automatic (cGH const * const cctkGH,
+ gh<dim> const & hh,
+ int const reflevel,
+ int const map,
+ int const size,
+ jjvect const & nboundaryzones,
+ jjvect const & is_internal,
+ jjvect const & is_staggered,
+ jjvect const & shiftout,
+ gh<dim>::rexts & bbsss,
+ gh<dim>::rbnds & obss,
+ gh<dim>::rprocs & pss);
+
+ void Automatic_OneLevel (const cGH * const cctkGH,
+ const gh<dim> & hh,
+ const int reflevel,
+ const int minwidth,
+ const CCTK_REAL minfraction,
+ const CCTK_REAL maxerror,
+ const gf<CCTK_REAL,dim> & errorvar,
+ vector<ibbox> & bbs,
+ vector<bbvect> & obs);
+
+ void Automatic_Recursive (const cGH * const cctkGH,
+ const gh<dim> & hh,
+ const int minwidth,
+ const CCTK_REAL minfraction,
+ const CCTK_REAL maxerror,
+ const data<CCTK_REAL,dim> & errorvar,
+ list<ibbox> & bbl,
+ const ibbox & region);
+
+ void Automatic_Recombine (list<ibbox> & bbl1,
+ list<ibbox> & bbl2,
+ list<ibbox> & bbl,
+ const ibbox & iface,
+ const int dir);
} // namespace CarpetRegrid