summaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorAnton Khirnov <anton@khirnov.net>2019-08-13 17:52:25 +0200
committerAnton Khirnov <anton@khirnov.net>2019-08-13 17:52:25 +0200
commit203d9c252b454c5af5ab69aaa65debba438a3e70 (patch)
tree4640555fd351432ec2681652b1abdeec3dfb5d11
parenta3a65f085dbb604b19263587f69659fe6c5627dd (diff)
Print timing stats.
-rw-r--r--param.ccl2
-rw-r--r--schedule.ccl4
-rw-r--r--src/qms.c120
3 files changed, 125 insertions, 1 deletions
diff --git a/param.ccl b/param.ccl
index bab5dfe..b7d9bcf 100644
--- a/param.ccl
+++ b/param.ccl
@@ -80,7 +80,7 @@ RESTRICTED:
CCTK_INT stats_every "print elliptic solver stats every <count> coarsest-level steps" STEERABLE=always
{
0: :: ""
-} 0
+} 32
RESTRICTED:
STRING loglevel "logging level" STEERABLE=always
diff --git a/schedule.ccl b/schedule.ccl
index 21df6c5..e8241c1 100644
--- a/schedule.ccl
+++ b/schedule.ccl
@@ -36,6 +36,10 @@ if (CCTK_EQUALS(lapse_source, "QMS_MG")) {
LANG: C
} ""
+ SCHEDULE qms_mg_terminate_print_stats IN CCTK_TERMINATE {
+ LANG: C
+ } ""
+
#SCHEDULE quasimaximal_slicing_axi IN MoL_PseudoEvolution {
# LANG: C
#} "Quasimaximal slicing"
diff --git a/src/qms.c b/src/qms.c
index 8ba8b1e..7fbde0b 100644
--- a/src/qms.c
+++ b/src/qms.c
@@ -62,6 +62,25 @@ typedef struct QMSMGContext {
CoordPatch *patches;
int nb_patches;
+
+ /* timings */
+ int64_t time_solve;
+ int64_t count_solve;
+
+ int64_t time_export;
+ int64_t count_export;
+
+ int64_t time_mg2d;
+ int64_t count_mg2d;
+
+ int64_t time_fill;
+ int64_t count_fill;
+
+ int64_t time_bnd;
+ int64_t count_bnd;
+
+ int64_t time_eval;
+ int64_t count_eval;
} QMSMGContext;
static const struct {
@@ -76,6 +95,14 @@ static const struct {
{ "debug", MG2D_LOG_DEBUG },
};
+#include <sys/time.h>
+static inline int64_t gettime(void)
+{
+ struct timeval tv;
+ gettimeofday(&tv, NULL);
+ return (int64_t)tv.tv_sec * 1000000 + tv.tv_usec;
+}
+
static QMSMGContext *qms_context;
static int ctz(int a)
@@ -243,6 +270,46 @@ static CoordPatch *get_coord_patch(QMSMGContext *ms, int level)
return cp;
}
+static void coord_patch_free(CoordPatch *cp)
+{
+ mg2d_solver_free(&cp->solver);
+
+ free(cp->filter);
+}
+
+static void print_stats(QMSMGContext *ms)
+{
+ int orig_log_level;
+ int64_t total = ms->time_solve + ms->time_eval;
+
+ fprintf(stderr, "%2.2f%% eval: %ld runs; %g s total; %g ms avg per run\n",
+ (double)ms->time_eval * 100.0 / total, ms->count_eval,
+ ms->time_eval / 1e6, ms->time_eval / 1e3 / ms->count_eval);
+
+ fprintf(stderr, "%2.2f%% solve: %ld runs; %g s total; %g ms avg per run\n",
+ (double)ms->time_solve * 100.0 / total, ms->count_solve,
+ ms->time_solve / 1e6, ms->time_solve / 1e3 / ms->count_solve);
+ fprintf(stderr, " %2.2f%% fill coeffs %2.2f%% boundaries %2.2f%% mg2d %2.2f%% export\n",
+ (double)ms->time_fill * 100.0 / ms->time_solve,
+ (double)ms->time_bnd * 100.0 / ms->time_solve,
+ (double)ms->time_mg2d * 100.0 / ms->time_solve,
+ (double)ms->time_export * 100.0 / ms->time_solve);
+
+ orig_log_level = ms->log_level;
+ ms->log_level = MG2D_LOG_VERBOSE;
+
+ for (int i = 0; i < ms->nb_patches; i++) {
+ CoordPatch *cp = get_coord_patch(ms, i);
+ char indent_buf[64];
+
+ snprintf(indent_buf, sizeof(indent_buf), " [%d] ", cp->level);
+
+ if (cp->solver)
+ mg2d_print_stats(cp->solver, indent_buf);
+ }
+ ms->log_level = orig_log_level;
+}
+
#define FD4(arr, stride, h_inv) \
((-1.0 * arr[2 * stride] + 8.0 * arr[1 * stride] - 8.0 * arr[-1 * stride] + 1.0 * arr[-2 * stride]) * h_inv / 12.0)
#define FD24(arr, stride, h_inv) \
@@ -914,6 +981,7 @@ void qms_mg_solve(CCTK_ARGUMENTS)
const double time = ms->gh->cctk_time;
int64_t grid_size;
+ int64_t total_start, start;
int ret;
@@ -926,12 +994,20 @@ void qms_mg_solve(CCTK_ARGUMENTS)
fabs(time - W_val1_time[ms->solve_level]) > 1e-13)
return;
+ total_start = gettime();
+
cp = get_coord_patch(ms, reflevel);
grid_size = cctkGH->cctk_lsh[0] * cctkGH->cctk_lsh[1] * cctkGH->cctk_lsh[2];
+ start = gettime();
+
fill_eq_coeffs(ms, cp, cp->solver);
+ ms->time_fill += gettime() - start;
+ ms->count_fill++;
+
+ start = gettime();
if (reflevel > ms->solve_level) {
/* fill the solver boundary conditions */
if (cp->solver->local_start[1] + cp->solver->local_size[1] == cp->solver->domain_size) {
@@ -967,12 +1043,18 @@ void qms_mg_solve(CCTK_ARGUMENTS)
CCTK_WARN(0, "Error setting the initial guess");
}
}
+ ms->time_bnd += gettime() - start;
+ ms->count_bnd++;
fprintf(stderr, "%d qms solve: time %g\n", reflevel, time);
+ start = gettime();
ret = mg2d_solve(cp->solver);
+ ms->time_mg2d += gettime() - start;
+ ms->count_mg2d++;
if (ret < 0)
CCTK_WARN(0, "Error solving the quasi-maximal slicing equation");
+ start = gettime();
{
double *W_mg_1 = CCTK_VarDataPtr(cctkGH, 1, "QuasiMaximalSlicingMG::W_mg");
solution_to_grid(cp, cp->solver, W_mg);
@@ -1017,6 +1099,15 @@ void qms_mg_solve(CCTK_ARGUMENTS)
W_pred1_time[reflevel] = time;
}
+ ms->time_export += gettime() - start;
+ ms->count_export++;
+
+ ms->time_solve += gettime() - total_start;
+ ms->count_solve++;
+
+ if (stats_every > 0 && reflevel == ms->solve_level && ms->count_solve > 0 && ms->count_eval > 0 &&
+ !(ms->count_solve % stats_every))
+ print_stats(ms);
}
void qms_mg_sync(CCTK_ARGUMENTS)
@@ -1034,6 +1125,8 @@ void qms_mg_eval(CCTK_ARGUMENTS)
const double time = cctkGH->cctk_time;
const int64_t grid_size = cctkGH->cctk_lsh[0] * cctkGH->cctk_lsh[1] * cctkGH->cctk_lsh[2];
+ int64_t start;
+
int i, ret;
ms = qms_context;
@@ -1041,6 +1134,7 @@ void qms_mg_eval(CCTK_ARGUMENTS)
if (reflevel < ms->solve_level || time >= switchoff_time + 1.0)
return;
+ start = gettime();
{
double *u0 = W_pred0;
double *u1 = W_pred1;
@@ -1070,6 +1164,8 @@ void qms_mg_eval(CCTK_ARGUMENTS)
W[i] = (u1[i] * fact1 + u0[i] * fact0) * fact;
}
+ ms->time_eval += gettime() - start;
+ ms->count_eval++;
}
static int context_init(cGH *gh, int solve_level, int fd_stencil, int maxiter, int exact_size, int nb_cycles,
@@ -1115,6 +1211,21 @@ static int context_init(cGH *gh, int solve_level, int fd_stencil, int maxiter, i
return 0;
}
+static void context_free(QMSMGContext **pms)
+{
+ QMSMGContext *ms = *pms;
+
+ if (!ms)
+ return;
+
+ for (int i = 0; i < ms->nb_patches; i++)
+ coord_patch_free(&ms->patches[i]);
+ free(ms->patches);
+
+ free(ms);
+ *pms = NULL;
+}
+
void qms_mg_init(CCTK_ARGUMENTS)
{
DECLARE_CCTK_ARGUMENTS;
@@ -1169,3 +1280,12 @@ void qms_mg_inithist(CCTK_ARGUMENTS)
memset(W_val0, 0, sizeof(double) * grid_size);
memset(W_val1, 0, sizeof(double) * grid_size);
}
+
+void qms_mg_terminate_print_stats(CCTK_ARGUMENTS)
+{
+ const int reflevel = ctz(cctkGH->cctk_levfac[0]);
+ if (reflevel == 0 && qms_context) {
+ print_stats(qms_context);
+ context_free(&qms_context);
+ }
+}