aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorschnetter <schnetter@c78560ca-4b45-4335-b268-5f3340f3cb52>2004-06-17 11:42:30 +0000
committerschnetter <schnetter@c78560ca-4b45-4335-b268-5f3340f3cb52>2004-06-17 11:42:30 +0000
commit4c1fd0811071d8e9178a9ba1850fd6eda640fa4e (patch)
tree5562be94005543a68d3b73e6493518d4352eb83f
parent703e3a02858aad19477488f4eee3ea43af100c18 (diff)
Recalculate the coordinate grid functions in the postregrid bin.
git-svn-id: http://svn.cactuscode.org/arrangements/CactusBase/CartGrid3D/trunk@206 c78560ca-4b45-4335-b268-5f3340f3cb52
-rw-r--r--schedule.ccl5
-rw-r--r--src/CartGrid3D.c40
2 files changed, 45 insertions, 0 deletions
diff --git a/schedule.ccl b/schedule.ccl
index c0d5421..87a34ce 100644
--- a/schedule.ccl
+++ b/schedule.ccl
@@ -29,6 +29,11 @@ schedule CartGrid3D as SpatialCoordinates at CCTK_BASEGRID
LANG:C
} "Set up spatial 3D Cartesian coordinates on the GH"
+schedule CartGrid3D_SetCoordinates at CCTK_POSTREGRID
+{
+ LANG: C
+} "Set Coordinates after regridding"
+
schedule CartGrid3D_ApplyBC in BoundaryConditions
{
LANG: C
diff --git a/src/CartGrid3D.c b/src/CartGrid3D.c
index bcd2223..a9a5398 100644
--- a/src/CartGrid3D.c
+++ b/src/CartGrid3D.c
@@ -450,3 +450,43 @@ void CartGrid3D(CCTK_ARGUMENTS)
"GAINDEX");
ierr += Util_TableSetReal (coord_handle, cctk_delta_space[0], "DELTA");
}
+
+
+ /*@@
+ @routine CartGrid3D_SetCoordinates
+ @date 2004-06-17
+ @author Christian Ott
+ @desc
+ Sets up Cartesian coordinates.
+ @enddesc
+ @var CCTK_ARGUMENTS
+ @vdesc Cactus argument list
+ @vtype
+ @vio in/out
+ @endvar
+@@*/
+void CartGrid3D_SetCoordinates(CCTK_ARGUMENTS)
+{
+ int i, j, k, idx;
+
+ DECLARE_CCTK_ARGUMENTS;
+ DECLARE_CCTK_PARAMETERS;
+
+ CCTK_VInfo(CCTK_THORNSTRING,"Resetting coordinates after regridding.");
+
+ for(k=0; k < cctk_lsh[2]; k++)
+ {
+ for(j=0; j < cctk_lsh[1]; j++)
+ {
+ for(i=0; i < cctk_lsh[0]; i++)
+ {
+ idx = CCTK_GFINDEX3D (cctkGH, i, j, k);
+ x[idx] = CCTK_DELTA_SPACE(0)*(i+cctk_lbnd[0]) + CCTK_ORIGIN_SPACE(0);
+ y[idx] = CCTK_DELTA_SPACE(1)*(j+cctk_lbnd[1]) + CCTK_ORIGIN_SPACE(1);
+ z[idx] = CCTK_DELTA_SPACE(2)*(k+cctk_lbnd[2]) + CCTK_ORIGIN_SPACE(2);
+ r[idx] = sqrt(SQR(x[idx]) + SQR(y[idx]) + SQR(z[idx]));
+ }
+ }
+ }
+
+}