summaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorAnton Khirnov <anton@khirnov.net>2019-07-25 10:56:04 +0200
committerAnton Khirnov <anton@khirnov.net>2019-07-25 10:56:04 +0200
commit43c18a92f8430a26f56b790baea2fdc1607dd02b (patch)
tree04b1084407e6d9ffef66b40c4964e506fd1bc779
parent68e2c0b3166527e95b66c024fcc06395d059c6bf (diff)
Fix syncing the solution values.
-rw-r--r--interface.ccl5
-rw-r--r--schedule.ccl4
-rw-r--r--src/qms.c12
3 files changed, 15 insertions, 6 deletions
diff --git a/interface.ccl b/interface.ccl
index 37e039a..f742897 100644
--- a/interface.ccl
+++ b/interface.ccl
@@ -15,8 +15,9 @@ CCTK_INT FUNCTION MoLNumIntegratorSubsteps ()
USES FUNCTION MoLNumIntegratorSubsteps
public:
-REAL W_val0 TYPE=GF TIMELEVELS=1
-REAL W_val1 TYPE=GF TIMELEVELS=1
+REAL W_mg TYPE=GF TIMELEVELS=2
+REAL W_val0 TYPE=GF TIMELEVELS=1 tags='Prolongation="None"'
+REAL W_val1 TYPE=GF TIMELEVELS=1 tags='Prolongation="None"'
REAL W_val0_time TYPE=ARRAY DIM=1 SIZE=32 DISTRIB=constant
REAL W_val1_time TYPE=ARRAY DIM=1 SIZE=32 DISTRIB=constant
diff --git a/schedule.ccl b/schedule.ccl
index 70e90c0..52075bf 100644
--- a/schedule.ccl
+++ b/schedule.ccl
@@ -22,10 +22,11 @@ if (CCTK_EQUALS(lapse_source, "QMS_MG")) {
SYNC: W_pred1
SYNC: W_val0
SYNC: W_val1
+ SYNC: W_mg
} "Quasimaximal slicing solve W"
SCHEDULE qms_mg_sync IN CCTK_POSTSTEP BEFORE qms_mg_solve {
- SYNC: W_val1
+ SYNC: W_mg
LANG: C
} ""
@@ -41,6 +42,7 @@ if (CCTK_EQUALS(lapse_source, "QMS_MG")) {
# LANG: C
#} "Quasimaximal slicing"
+ STORAGE: W_mg[2]
STORAGE: W_pred0
STORAGE: W_pred1
STORAGE: W_pred0_time
diff --git a/src/qms.c b/src/qms.c
index de2bb12..162bd47 100644
--- a/src/qms.c
+++ b/src/qms.c
@@ -899,7 +899,7 @@ void qms_mg_solve(CCTK_ARGUMENTS)
double *dst = bnd->val + j * bnd->val_stride;
for (ptrdiff_t i = -j; i < (ptrdiff_t)cp->solver->local_size[0] + j; i++) {
const ptrdiff_t idx = CCTK_GFINDEX3D(cctkGH, ABS(i) + cp->offset_left[0], cp->y_idx, cp->offset_left[1] + cp->solver->local_size[1] - 1 + j);
- dst[i] = W_val1[idx];
+ dst[i] = W_mg[idx];
}
}
}
@@ -910,7 +910,7 @@ void qms_mg_solve(CCTK_ARGUMENTS)
double *dst = bnd->val + j * bnd->val_stride;
for (ptrdiff_t i = -j; i < (ptrdiff_t)cp->solver->local_size[1] + j; i++) {
const ptrdiff_t idx = CCTK_GFINDEX3D(cctkGH, cp->offset_left[1] + cp->solver->local_size[0] - 1 + j, cp->y_idx, cp->offset_left[1] + ABS(i));
- dst[i] = W_val1[idx];
+ dst[i] = W_mg[idx];
}
}
}
@@ -921,12 +921,18 @@ void qms_mg_solve(CCTK_ARGUMENTS)
if (ret < 0)
CCTK_WARN(0, "Error solving the quasi-maximal slicing equation");
+ {
+ double *W_mg_1 = CCTK_VarDataPtr(cctkGH, 1, "QuasiMaximalSlicingMG::W_mg");
+ solution_to_grid(cp, cp->solver, W_mg);
+ memcpy(W_mg_1, W_mg, grid_size * sizeof(*W_mg_1));
+ }
+
/* add the solution to the list of past solutions */
{
memcpy(W_val0, W_val1, sizeof(*W_val0) * grid_size);
W_val0_time[reflevel] = W_val1_time[reflevel];
- solution_to_grid(cp, cp->solver, W_val1);
+ memcpy(W_val1, W_mg, sizeof(*W_val1) * grid_size);
W_val1_time[reflevel] = ms->gh->cctk_time;
}