aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorErik Schnetter <schnetter@cct.lsu.edu>2007-03-11 20:40:00 +0000
committerErik Schnetter <schnetter@cct.lsu.edu>2007-03-11 20:40:00 +0000
commitebee4f847dfd6405b0e711424bac4e58af5b2ec5 (patch)
treeda6357687e1eb94d91ee9b5e6a163d36d8150957
parent1fb9122f23b55ad4bacd7dee3e2752ca50547a50 (diff)
CarpetRegrid2: Regrid only if the tracked objects have moved sufficiently
Regrid only if the tracked objects have moved a certain minimum distance. darcs-hash:20070311204027-dae7b-7c39e3c51d6b617188b477aed73facd5fe7e6fe8.gz
-rw-r--r--Carpet/CarpetRegrid2/interface.ccl7
-rw-r--r--Carpet/CarpetRegrid2/param.ccl70
-rw-r--r--Carpet/CarpetRegrid2/schedule.ccl2
-rw-r--r--Carpet/CarpetRegrid2/src/initialise.cc6
-rw-r--r--Carpet/CarpetRegrid2/src/regrid.cc35
5 files changed, 119 insertions, 1 deletions
diff --git a/Carpet/CarpetRegrid2/interface.ccl b/Carpet/CarpetRegrid2/interface.ccl
index 186ea206e..7c3a4e90f 100644
--- a/Carpet/CarpetRegrid2/interface.ccl
+++ b/Carpet/CarpetRegrid2/interface.ccl
@@ -92,3 +92,10 @@ CCTK_REAL radii[10] TYPE=array DIM=1 SIZE=30 DISTRIB=constant
{
radius
} "Radii of refined regions for each level"
+
+PRIVATE:
+
+CCTK_REAL old_positions[10] TYPE=scalar
+{
+ old_position_x old_position_y old_position_z
+} "Old positions of refined regions"
diff --git a/Carpet/CarpetRegrid2/param.ccl b/Carpet/CarpetRegrid2/param.ccl
index b4b088405..97a42ce57 100644
--- a/Carpet/CarpetRegrid2/param.ccl
+++ b/Carpet/CarpetRegrid2/param.ccl
@@ -83,6 +83,13 @@ CCTK_REAL radius_1[30] "Radius of refined region for this level of this centre"
+CCTK_REAL movement_threshold_1 "Minimum movement to trigger a regridding" STEERABLE=always
+{
+ 0:* :: ""
+} 0.0
+
+
+
################################################################################
CCTK_INT num_levels_2 "Number of refinement levels for this centre" STEERABLE=always
@@ -116,6 +123,13 @@ CCTK_REAL radius_2[30] "Radius of refined region for this level of this centre"
+CCTK_REAL movement_threshold_2 "Minimum movement to trigger a regridding" STEERABLE=always
+{
+ 0:* :: ""
+} 0.0
+
+
+
################################################################################
CCTK_INT num_levels_3 "Number of refinement levels for this centre" STEERABLE=always
@@ -149,6 +163,13 @@ CCTK_REAL radius_3[30] "Radius of refined region for this level of this centre"
+CCTK_REAL movement_threshold_3 "Minimum movement to trigger a regridding" STEERABLE=always
+{
+ 0:* :: ""
+} 0.0
+
+
+
################################################################################
CCTK_INT num_levels_4 "Number of refinement levels for this centre" STEERABLE=always
@@ -182,6 +203,13 @@ CCTK_REAL radius_4[30] "Radius of refined region for this level of this centre"
+CCTK_REAL movement_threshold_4 "Minimum movement to trigger a regridding" STEERABLE=always
+{
+ 0:* :: ""
+} 0.0
+
+
+
################################################################################
CCTK_INT num_levels_5 "Number of refinement levels for this centre" STEERABLE=always
@@ -215,6 +243,13 @@ CCTK_REAL radius_5[30] "Radius of refined region for this level of this centre"
+CCTK_REAL movement_threshold_5 "Minimum movement to trigger a regridding" STEERABLE=always
+{
+ 0:* :: ""
+} 0.0
+
+
+
################################################################################
CCTK_INT num_levels_6 "Number of refinement levels for this centre" STEERABLE=always
@@ -248,6 +283,13 @@ CCTK_REAL radius_6[30] "Radius of refined region for this level of this centre"
+CCTK_REAL movement_threshold_6 "Minimum movement to trigger a regridding" STEERABLE=always
+{
+ 0:* :: ""
+} 0.0
+
+
+
################################################################################
CCTK_INT num_levels_7 "Number of refinement levels for this centre" STEERABLE=always
@@ -281,6 +323,13 @@ CCTK_REAL radius_7[30] "Radius of refined region for this level of this centre"
+CCTK_REAL movement_threshold_7 "Minimum movement to trigger a regridding" STEERABLE=always
+{
+ 0:* :: ""
+} 0.0
+
+
+
################################################################################
CCTK_INT num_levels_8 "Number of refinement levels for this centre" STEERABLE=always
@@ -314,6 +363,13 @@ CCTK_REAL radius_8[30] "Radius of refined region for this level of this centre"
+CCTK_REAL movement_threshold_8 "Minimum movement to trigger a regridding" STEERABLE=always
+{
+ 0:* :: ""
+} 0.0
+
+
+
################################################################################
CCTK_INT num_levels_9 "Number of refinement levels for this centre" STEERABLE=always
@@ -347,6 +403,13 @@ CCTK_REAL radius_9[30] "Radius of refined region for this level of this centre"
+CCTK_REAL movement_threshold_9 "Minimum movement to trigger a regridding" STEERABLE=always
+{
+ 0:* :: ""
+} 0.0
+
+
+
################################################################################
CCTK_INT num_levels_10 "Number of refinement levels for this centre" STEERABLE=always
@@ -377,3 +440,10 @@ CCTK_REAL radius_10[30] "Radius of refined region for this level of this centre"
{
0:* :: ""
} 1.0
+
+
+
+CCTK_REAL movement_threshold_10 "Minimum movement to trigger a regridding" STEERABLE=always
+{
+ 0:* :: ""
+} 0.0
diff --git a/Carpet/CarpetRegrid2/schedule.ccl b/Carpet/CarpetRegrid2/schedule.ccl
index 6299ba209..431e7d7d9 100644
--- a/Carpet/CarpetRegrid2/schedule.ccl
+++ b/Carpet/CarpetRegrid2/schedule.ccl
@@ -1,7 +1,7 @@
# Schedule definitions for thorn CarpetRegrid2
STORAGE: last_iteration last_map
-STORAGE: num_levels positions radii
+STORAGE: num_levels positions radii old_positions
diff --git a/Carpet/CarpetRegrid2/src/initialise.cc b/Carpet/CarpetRegrid2/src/initialise.cc
index 3904a1d71..b2b7104db 100644
--- a/Carpet/CarpetRegrid2/src/initialise.cc
+++ b/Carpet/CarpetRegrid2/src/initialise.cc
@@ -134,6 +134,12 @@ namespace CarpetRegrid2 {
radius[index2 (lsh, rl, 9)] = radius_10[rl];
}
}
+
+ for (int n = 0; n < num_centres; ++ n) {
+ old_position_x[n] = position_x[n];
+ old_position_y[n] = position_y[n];
+ old_position_z[n] = position_z[n];
+ }
}
} // namespace CarpetRegrid2
diff --git a/Carpet/CarpetRegrid2/src/regrid.cc b/Carpet/CarpetRegrid2/src/regrid.cc
index 7ac61c7b4..9ee669564 100644
--- a/Carpet/CarpetRegrid2/src/regrid.cc
+++ b/Carpet/CarpetRegrid2/src/regrid.cc
@@ -561,6 +561,34 @@ namespace CarpetRegrid2 {
Carpet::map > * last_map))));
}
}
+
+ if (do_recompose and cctk_iteration != 0) {
+ // Regrid only if the positions have changed sufficiently
+ do_recompose = false;
+ for (int n = 0; n < num_centres; ++ n) {
+ CCTK_REAL const dist2 =
+ pow (position_x[n] - old_position_x[n], 2) +
+ pow (position_y[n] - old_position_y[n], 2) +
+ pow (position_z[n] - old_position_z[n], 2);
+ CCTK_REAL mindist;
+ switch (n) {
+ case 0: mindist = movement_threshold_1; break;
+ case 1: mindist = movement_threshold_2; break;
+ case 2: mindist = movement_threshold_3; break;
+ case 3: mindist = movement_threshold_4; break;
+ case 4: mindist = movement_threshold_5; break;
+ case 5: mindist = movement_threshold_6; break;
+ case 6: mindist = movement_threshold_7; break;
+ case 7: mindist = movement_threshold_8; break;
+ case 8: mindist = movement_threshold_9; break;
+ case 9: mindist = movement_threshold_10; break;
+ default: assert (0);
+ }
+ do_recompose = dist2 >= pow (mindist, 2);
+ if (do_recompose) break;
+ } // for n
+ }
+
if (do_recompose) {
* last_iteration = cctk_iteration;
* last_map = Carpet::map;
@@ -583,6 +611,13 @@ namespace CarpetRegrid2 {
// Make multigrid aware
Carpet::MakeMultigridBoxes (cctkGH, regss, regsss);
+ // Remember current positions
+ for (int n = 0; n < num_centres; ++ n) {
+ old_position_x[n] = position_x[n];
+ old_position_y[n] = position_y[n];
+ old_position_z[n] = position_z[n];
+ }
+
} // if do_recompose
return do_recompose;