summaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorAnton Khirnov <anton@khirnov.net>2018-12-28 13:09:13 +0100
committerAnton Khirnov <anton@khirnov.net>2018-12-28 13:09:13 +0100
commit30816fb1931f8ea3a70215daf77636eb425cc91a (patch)
tree4d7653b87a91606249e1138b6423ba1b044a3df8
parent133c5ff7ebce3bb4301a841bf9de9ba85fca23a2 (diff)
Make more solver parameters runtime-configurable.
-rw-r--r--param.ccl30
-rw-r--r--src/maximal_slicing_axi_mg.c28
2 files changed, 52 insertions, 6 deletions
diff --git a/param.ccl b/param.ccl
index 0795ac4..204a942 100644
--- a/param.ccl
+++ b/param.ccl
@@ -14,6 +14,36 @@ CCTK_INT fd_stencil "finite differencing stencil"
} 1
RESTRICTED:
+CCTK_INT maxiter "maximum number of solver iterations"
+{
+ 1: :: ""
+} 64
+
+RESTRICTED:
+CCTK_INT nb_cycles "number of coarse-correction cycles per level"
+{
+ 1: :: ""
+} 1
+
+RESTRICTED:
+CCTK_INT nb_relax_pre "number of relaxation steps before coarse-grid correction"
+{
+ 0: :: ""
+} 2
+
+RESTRICTED:
+CCTK_INT nb_relax_post "number of relaxation steps after coarse-grid correction"
+{
+ 0: :: ""
+} 2
+
+RESTRICTED:
+CCTK_REAL tol_residual "maximum absolute value of the residual"
+{
+ 0: :: ""
+} 1e-12
+
+RESTRICTED:
CCTK_INT stats_every "print elliptic solver stats every <count> coarsest-level steps"
{
0: :: ""
diff --git a/src/maximal_slicing_axi_mg.c b/src/maximal_slicing_axi_mg.c
index 0fe9a54..93130d7 100644
--- a/src/maximal_slicing_axi_mg.c
+++ b/src/maximal_slicing_axi_mg.c
@@ -78,6 +78,11 @@ typedef struct MSMGContext {
cGH *gh;
int fd_stencil;
+ int maxiter;
+ int nb_cycles;
+ int nb_relax_post;
+ int nb_relax_pre;
+ double tol_residual;
CoordPatch *patches;
int nb_patches;
@@ -236,9 +241,11 @@ static CoordPatch *get_coord_patch(MSMGContext *ms, int level)
cp->solver->boundaries[MG2D_BOUNDARY_0U]->type = MG2D_BC_TYPE_FIXVAL;
cp->solver->boundaries[MG2D_BOUNDARY_1U]->type = MG2D_BC_TYPE_FIXVAL;
- cp->solver->maxiter = 32;
- cp->solver->tol = 1e-10;
- cp->solver->nb_cycles = 2;
+ cp->solver->maxiter = ms->maxiter;
+ cp->solver->tol = ms->tol_residual;
+ cp->solver->nb_cycles = ms->nb_cycles;
+ cp->solver->nb_relax_post = ms->nb_relax_post;
+ cp->solver->nb_relax_pre = ms->nb_relax_pre;
cp->solver->opaque = ms;
cp->solver->log_callback = log_callback;
@@ -318,7 +325,9 @@ static void print_stats(MSMGContext *ms)
ms->log_level = orig_log_level;
}
-static int context_init(cGH *gh, int fd_stencil, const char *loglevel_str,
+static int context_init(cGH *gh, int fd_stencil, int maxiter, int nb_cycles,
+ int nb_relax_pre, int nb_relax_post, double tol_residual,
+ const char *loglevel_str,
MSMGContext **ctx)
{
MSMGContext *ms;
@@ -329,7 +338,12 @@ static int context_init(cGH *gh, int fd_stencil, const char *loglevel_str,
return -ENOMEM;
ms->gh = gh;
- ms->fd_stencil = fd_stencil;
+ ms->fd_stencil = fd_stencil;
+ ms->maxiter = maxiter;
+ ms->nb_cycles = nb_cycles;
+ ms->nb_relax_pre = nb_relax_pre;
+ ms->nb_relax_post = nb_relax_post;
+ ms->tol_residual = tol_residual;
for (int i = 0; i < ARRAY_ELEMS(log_levels); i++) {
if (!strcmp(loglevel_str, log_levels[i].str)) {
@@ -869,7 +883,9 @@ void msa_mg_init(CCTK_ARGUMENTS)
DECLARE_CCTK_PARAMETERS;
int ret;
- ret = context_init(cctkGH, fd_stencil, loglevel, &ms);
+ ret = context_init(cctkGH, fd_stencil, maxiter, nb_cycles,
+ nb_relax_pre, nb_relax_post, tol_residual,
+ loglevel, &ms);
if (ret < 0)
CCTK_WARN(0, "Error initializing the solver context");
}