diff options
author | schnetter <> | 2003-08-10 10:52:00 +0000 |
---|---|---|
committer | schnetter <> | 2003-08-10 10:52:00 +0000 |
commit | de6ebb94b4f69684d5cf21204bded61ff3e5b094 (patch) | |
tree | 411b1b8d453234b2638349e21a84e7af2c59cc6e /Carpet | |
parent | 4f971f6a5a474559ae876affe86edbd1272cc68e (diff) |
Fix bug in three timelevel initialisation:
Fix bug in three timelevel initialisation:
Change order of time level flipping and restricting.
ggf.cc: Remove commented out code.
darcs-hash:20030810105209-07bb3-40c12c5fb838cd56a99f26c8ca71529bfe11a900.gz
Diffstat (limited to 'Carpet')
-rw-r--r-- | Carpet/Carpet/src/Comm.cc | 4 | ||||
-rw-r--r-- | Carpet/Carpet/src/Initialise.cc | 19 | ||||
-rw-r--r-- | Carpet/Carpet/src/Restrict.cc | 10 | ||||
-rw-r--r-- | Carpet/CarpetLib/src/ggf.cc | 6 |
4 files changed, 15 insertions, 24 deletions
diff --git a/Carpet/Carpet/src/Comm.cc b/Carpet/Carpet/src/Comm.cc index c445bb3e9..e00efbfa6 100644 --- a/Carpet/Carpet/src/Comm.cc +++ b/Carpet/Carpet/src/Comm.cc @@ -10,7 +10,7 @@ #include "carpet.hh" extern "C" { - static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/Carpet/src/Comm.cc,v 1.20 2003/07/20 21:03:43 schnetter Exp $"; + static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/Carpet/src/Comm.cc,v 1.21 2003/08/10 12:52:09 schnetter Exp $"; CCTK_FILEVERSION(Carpet_Carpet_Comm_cc); } @@ -75,7 +75,7 @@ namespace Carpet { // user) const CCTK_REAL time = (cgh->cctk_time - cctk_initial_time) / delta_time; #if 0 - const CCTK_REAL time1 = tt->time (tl, reflevel, mglevel); + const CCTK_REAL time1 = tt->time (tl, reflevel, mglevel) / tt->get_delta (0, 0); const CCTK_REAL time2 = (cgh->cctk_time - cctk_initial_time) / delta_time; assert (fabs((time1 - time2) / (fabs(time1) + fabs(time2) + fabs(cgh->cctk_delta_time))) < 1e-12); #endif diff --git a/Carpet/Carpet/src/Initialise.cc b/Carpet/Carpet/src/Initialise.cc index 49b1b8da1..d7baabef3 100644 --- a/Carpet/Carpet/src/Initialise.cc +++ b/Carpet/Carpet/src/Initialise.cc @@ -12,7 +12,7 @@ #include "carpet.hh" extern "C" { - static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/Carpet/src/Initialise.cc,v 1.30 2003/06/18 18:24:27 schnetter Exp $"; + static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/Carpet/src/Initialise.cc,v 1.31 2003/08/10 12:52:09 schnetter Exp $"; CCTK_FILEVERSION(Carpet_Carpet_Initialise_cc); } @@ -229,7 +229,6 @@ namespace Carpet { BEGIN_REFLEVEL_LOOP(cgh) { BEGIN_MGLEVEL_LOOP(cgh) { - // Cycle time levels (ignore arrays) cout << "3TL rl=" << reflevel << " cycling" << endl; CycleTimeLevels (cgh); @@ -304,15 +303,6 @@ namespace Carpet { BEGIN_REVERSE_REFLEVEL_LOOP(cgh) { BEGIN_MGLEVEL_LOOP(cgh) { - cgh->cctk_time = cctk_initial_time + cgh->cctk_delta_time / cgh->cctk_timefac; - - // Restrict - cout << "3TL rl=" << reflevel << " restricting" << endl; - Restrict (cgh); - - Waypoint ("%*sScheduling POSTRESTRICT", 2*reflevel, ""); - CCTK_ScheduleTraverse ("POSTRESTRICT", cgh, CallFunction); - // Flip time levels cout << "3TL rl=" << reflevel << " flipping" << endl; FlipTimeLevels (cgh); @@ -330,6 +320,13 @@ namespace Carpet { << " time=" << tt->get_time (reflevel, mglevel) << " time=" << cgh->cctk_time / cgh->cctk_delta_time << endl; + // Restrict + cout << "3TL rl=" << reflevel << " restricting" << endl; + Restrict (cgh); + + Waypoint ("%*sScheduling POSTRESTRICT", 2*reflevel, ""); + CCTK_ScheduleTraverse ("POSTRESTRICT", cgh, CallFunction); + // Cycle time levels cout << "3TL rl=" << reflevel << " cycling" << endl; CycleTimeLevels (cgh); diff --git a/Carpet/Carpet/src/Restrict.cc b/Carpet/Carpet/src/Restrict.cc index e392af3ca..1467cdd58 100644 --- a/Carpet/Carpet/src/Restrict.cc +++ b/Carpet/Carpet/src/Restrict.cc @@ -10,7 +10,7 @@ #include "carpet.hh" extern "C" { - static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/Carpet/src/Restrict.cc,v 1.19 2003/07/20 21:03:43 schnetter Exp $"; + static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/Carpet/src/Restrict.cc,v 1.20 2003/08/10 12:52:09 schnetter Exp $"; CCTK_FILEVERSION(Carpet_Carpet_Restrict_cc); } @@ -40,11 +40,9 @@ namespace Carpet { // use background time here (which may not be modified by // the user) const CCTK_REAL time = tt->time (tl, reflevel, mglevel); - if (tl==0) { - const CCTK_REAL time1 = tt->time (tl, reflevel, mglevel); - const CCTK_REAL time2 = cgh->cctk_time / cgh->cctk_delta_time; - assert (fabs((time1 - time2) / (fabs(time1) + fabs(time2) + fabs(cgh->cctk_delta_time))) < 1e-12); - } + const CCTK_REAL time1 = tt->time (0, reflevel, mglevel) / tt->get_delta (0, 0); + const CCTK_REAL time2 = cgh->cctk_time / cgh->cctk_delta_time; + assert (fabs((time1 - time2) / (fabs(time1) + fabs(time2) + fabs(cgh->cctk_delta_time))) < 1e-12); if (reflevel < hh->reflevels()-1) { diff --git a/Carpet/CarpetLib/src/ggf.cc b/Carpet/CarpetLib/src/ggf.cc index ba33f8a48..d8b39cec6 100644 --- a/Carpet/CarpetLib/src/ggf.cc +++ b/Carpet/CarpetLib/src/ggf.cc @@ -1,4 +1,4 @@ -// $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/ggf.cc,v 1.25 2003/06/18 18:24:28 schnetter Exp $ +// $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetLib/src/ggf.cc,v 1.26 2003/08/10 12:52:09 schnetter Exp $ #include <assert.h> #include <stdlib.h> @@ -478,10 +478,6 @@ void ggf<D>::ref_restrict (int tl, int rl, int c, int ml, CCTK_REAL time) { // Require same times - // SHH: removed assert and added warning -// if (t.get_time(rl,ml) != t.get_time(rl+1,ml)) { -// printf("WARNING: Time on rl %d is %d, time on rl %d is %d\n",rl,t.get_time(rl,ml),rl+1,t.get_time(rl+1,ml)); -// } assert (t.get_time(rl,ml) == t.get_time(rl+1,ml)); const vector<int> tl2s(1,tl); intercat (tl ,rl ,c,ml, &dh<D>::dboxes::recv_ref_fine, |