Commit df254615 authored by Matthieu Schaller's avatar Matthieu Schaller
Browse files

Removed the energy statistics from the drift

parent 09aeffb8
......@@ -44,7 +44,7 @@ include_HEADERS = space.h runner.h queue.h task.h lock.h cell.h part.h const.h \
common_io.h single_io.h multipole.h map.h tools.h partition.h clocks.h parser.h \
physical_constants.h physical_constants_cgs.h potential.h version.h \
hydro_properties.h riemann.h threadpool.h cooling.h cooling_struct.h sourceterms.h \
sourceterms_struct.h
sourceterms_struct.h statistics.h
# Common source files
......@@ -53,7 +53,8 @@ AM_SOURCES = space.c runner.c queue.c task.c cell.c engine.c \
units.c common_io.c single_io.c multipole.c version.c map.c \
kernel_hydro.c tools.c part.c partition.c clocks.c parser.c \
physical_constants.c potential.c hydro_properties.c \
runner_doiact_fft.c threadpool.c cooling.c sourceterms.c
runner_doiact_fft.c threadpool.c cooling.c sourceterms.c \
statistics.c
# Include files for distribution, not installation.
nobase_noinst_HEADERS = align.h approx_math.h atomic.h cycle.h error.h inline.h kernel_hydro.h kernel_gravity.h \
......
......@@ -179,12 +179,6 @@ struct cell {
/* ID of the previous owner, e.g. runner. */
int owner;
/* Momentum of particles in cell. */
double mom[3], ang_mom[3];
/* Mass, potential, internal and kinetic energy of particles in this cell. */
double mass, e_pot, e_int, e_kin, e_rad, entropy;
/* Number of particles updated in this cell. */
int updated, g_updated;
......
......@@ -63,6 +63,7 @@
#include "runner.h"
#include "serial_io.h"
#include "single_io.h"
#include "statistics.h"
#include "timers.h"
#include "tools.h"
#include "units.h"
......@@ -1249,9 +1250,8 @@ void engine_make_external_gravity_tasks(struct engine *e) {
if (ci->nodeID != nodeID) continue;
/* If the cells is local build a self-interaction */
ci->grav_external = scheduler_addtask(sched, task_type_self,
task_subtype_external_grav, 0, 0,
ci, NULL, 0);
ci->grav_external = scheduler_addtask(
sched, task_type_self, task_subtype_external_grav, 0, 0, ci, NULL, 0);
}
}
......@@ -2462,74 +2462,60 @@ void engine_print_stats(struct engine *e) {
const ticks tic = getticks();
const struct space *s = e->s;
double e_kin = 0.0, e_int = 0.0, e_pot = 0.0, e_rad = 0.0;
double entropy = 0.0, mass = 0.0;
double mom[3] = {0.0, 0.0, 0.0}, ang_mom[3] = {0.0, 0.0, 0.0};
struct statistics stats;
bzero(&stats, sizeof(struct statistics));
/* Collect the cell data. */
for (int k = 0; k < s->nr_cells; k++)
if (s->cells_top[k].nodeID == e->nodeID) {
struct cell *c = &s->cells_top[k];
mass += c->mass;
e_kin += c->e_kin;
e_int += c->e_int;
e_pot += c->e_pot;
e_rad += c->e_rad;
entropy += c->entropy;
mom[0] += c->mom[0];
mom[1] += c->mom[1];
mom[2] += c->mom[2];
ang_mom[0] += c->ang_mom[0];
ang_mom[1] += c->ang_mom[1];
ang_mom[2] += c->ang_mom[2];
}
/* Collect the stats on this node */
stats_collect(s, &stats);
/* Aggregate the data from the different nodes. */
#ifdef WITH_MPI
{
double in[12] = {0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.};
double out[12];
out[0] = e_kin;
out[1] = e_int;
out[2] = e_pot;
out[3] = e_rad;
out[4] = mom[0];
out[5] = mom[1];
out[6] = mom[2];
out[7] = ang_mom[0];
out[8] = ang_mom[1];
out[9] = ang_mom[2];
out[10] = mass;
out[11] = entropy;
if (MPI_Reduce(out, in, 11, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD) !=
MPI_SUCCESS)
error("Failed to aggregate stats.");
e_kin = out[0];
e_int = out[1];
e_pot = out[2];
e_rad = out[3];
mom[0] = out[4];
mom[1] = out[5];
mom[2] = out[6];
ang_mom[0] = out[7];
ang_mom[1] = out[8];
ang_mom[2] = out[9];
mass = out[10];
entropy = out[11];
}
/* { */
/* double in[12] = {0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.}; */
/* double out[12]; */
/* out[0] = e_kin; */
/* out[1] = e_int; */
/* out[2] = e_pot; */
/* out[3] = e_rad; */
/* out[4] = mom[0]; */
/* out[5] = mom[1]; */
/* out[6] = mom[2]; */
/* out[7] = ang_mom[0]; */
/* out[8] = ang_mom[1]; */
/* out[9] = ang_mom[2]; */
/* out[10] = mass; */
/* out[11] = entropy; */
/* if (MPI_Reduce(out, in, 11, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD) != */
/* MPI_SUCCESS) */
/* error("Failed to aggregate stats."); */
/* e_kin = out[0]; */
/* e_int = out[1]; */
/* e_pot = out[2]; */
/* e_rad = out[3]; */
/* mom[0] = out[4]; */
/* mom[1] = out[5]; */
/* mom[2] = out[6]; */
/* ang_mom[0] = out[7]; */
/* ang_mom[1] = out[8]; */
/* ang_mom[2] = out[9]; */
/* mass = out[10]; */
/* entropy = out[11]; */
/* } */
#endif
const double e_tot = e_kin + e_int + e_pot;
// const double e_tot = e_kin + e_int + e_pot;
/* Print info */
if (e->nodeID == 0) {
fprintf(e->file_stats,
" %14e %14e %14e %14e %14e %14e %14e %14e %14e %14e %14e %14e %14e "
"%14e\n",
e->time, mass, e_tot, e_kin, e_int, e_pot, e_rad, entropy, mom[0],
mom[1], mom[2], ang_mom[0], ang_mom[1], ang_mom[2]);
fflush(e->file_stats);
}
/* if (e->nodeID == 0) { */
/* fprintf(e->file_stats, */
/* " %14e %14e %14e %14e %14e %14e %14e %14e %14e %14e %14e %14e
* %14e " */
/* "%14e\n", */
/* e->time, mass, e_tot, e_kin, e_int, e_pot, e_rad, entropy,
* mom[0], */
/* mom[1], mom[2], ang_mom[0], ang_mom[1], ang_mom[2]); */
/* fflush(e->file_stats); */
/* } */
if (e->verbose)
message("took %.3f %s.", clocks_from_ticks(getticks() - tic),
......@@ -2869,8 +2855,8 @@ void engine_drift(struct engine *e) {
e->s->nr_cells, sizeof(struct cell), 1, e);
if (e->verbose)
message("took %.3f %s (including task unskipping).", clocks_from_ticks(getticks() - tic),
clocks_getunit());
message("took %.3f %s (including task unskipping).",
clocks_from_ticks(getticks() - tic), clocks_getunit());
}
/**
......
......@@ -790,12 +790,7 @@ static void runner_do_drift(struct cell *c, struct engine *e, int drift) {
/* Drift from the last time the cell was drifted to the current time */
const double dt = (ti_current - ti_old) * timeBase;
float dx_max = 0.f, dx2_max = 0.f, h_max = 0.f;
double e_kin = 0.0, e_int = 0.0, e_pot = 0.0, e_rad = 0.0;
double entropy = 0.0, mass = 0.0;
double mom[3] = {0.0, 0.0, 0.0};
double ang_mom[3] = {0.0, 0.0, 0.0};
/* No children? */
if (!c->split) {
......@@ -820,7 +815,7 @@ static void runner_do_drift(struct cell *c, struct engine *e, int drift) {
dx2_max = (dx2_max > dx2) ? dx2_max : dx2;
}
/* Loop over all the particles in the cell (more work for these !) */
/* Loop over all the particles in the cell */
const size_t nr_parts = c->count;
for (size_t k = 0; k < nr_parts; k++) {
......@@ -839,39 +834,6 @@ static void runner_do_drift(struct cell *c, struct engine *e, int drift) {
/* Maximal smoothing length */
h_max = (h_max > p->h) ? h_max : p->h;
/* Now collect quantities for statistics */
const float half_dt =
(ti_current - (p->ti_begin + p->ti_end) / 2) * timeBase;
const double x[3] = {p->x[0], p->x[1], p->x[2]};
const float v[3] = {xp->v_full[0] + p->a_hydro[0] * half_dt,
xp->v_full[1] + p->a_hydro[1] * half_dt,
xp->v_full[2] + p->a_hydro[2] * half_dt};
const float m = hydro_get_mass(p);
/* Collect mass */
mass += m;
/* Collect momentum */
mom[0] += m * v[0];
mom[1] += m * v[1];
mom[2] += m * v[2];
/* Collect angular momentum */
ang_mom[0] += m * (x[1] * v[2] - x[2] * v[1]);
ang_mom[1] += m * (x[2] * v[0] - x[0] * v[2]);
ang_mom[2] += m * (x[0] * v[1] - x[1] * v[0]);
/* Collect energies. */
e_kin += 0.5 * m * (v[0] * v[0] + v[1] * v[1] + v[2] * v[2]);
e_pot += 0.;
e_int += m * hydro_get_internal_energy(p, half_dt);
e_rad += cooling_get_radiated_energy(xp);
/* Collect entropy */
entropy += m * hydro_get_entropy(p, half_dt);
}
/* Now, get the maximal particle motion from its square */
......@@ -892,36 +854,12 @@ static void runner_do_drift(struct cell *c, struct engine *e, int drift) {
runner_do_drift(cp, e, drift);
dx_max = max(dx_max, cp->dx_max);
h_max = max(h_max, cp->h_max);
mass += cp->mass;
e_kin += cp->e_kin;
e_int += cp->e_int;
e_pot += cp->e_pot;
e_rad += cp->e_rad;
entropy += cp->entropy;
mom[0] += cp->mom[0];
mom[1] += cp->mom[1];
mom[2] += cp->mom[2];
ang_mom[0] += cp->ang_mom[0];
ang_mom[1] += cp->ang_mom[1];
ang_mom[2] += cp->ang_mom[2];
}
}
/* Store the values */
c->h_max = h_max;
c->dx_max = dx_max;
c->mass = mass;
c->e_kin = e_kin;
c->e_int = e_int;
c->e_pot = e_pot;
c->e_rad = e_rad;
c->entropy = entropy;
c->mom[0] = mom[0];
c->mom[1] = mom[1];
c->mom[2] = mom[2];
c->ang_mom[0] = ang_mom[0];
c->ang_mom[1] = ang_mom[1];
c->ang_mom[2] = ang_mom[2];
/* Update the time of the last drift */
c->ti_old = ti_current;
......
......@@ -601,8 +601,7 @@ void space_rebuild(struct space *s, double cell_max, int verbose) {
for (size_t k = 1; k < nr_parts; k++) {
if (ind[k - 1] > ind[k]) {
error("Sort failed!");
} else if (ind[k] != cell_getid(s->cdim,
s->parts[k].x[0] * s->iwidth[0],
} else if (ind[k] != cell_getid(s->cdim, s->parts[k].x[0] * s->iwidth[0],
s->parts[k].x[1] * s->iwidth[1],
s->parts[k].x[2] * s->iwidth[2])) {
error("Incorrect indices!");
......
/*******************************************************************************
* This file is part of SWIFT.
* Copyright (c) 2016 Matthieu Schaller (matthieu.schaller@durham.ac.uk)
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU Lesser General Public License as published
* by the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* This program is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU Lesser General Public License
* along with this program. If not, see <http://www.gnu.org/licenses/>.
*
******************************************************************************/
/* Config parameters. */
#include "../config.h"
/* Some standard headers. */
#include <string.h>
/* This object's header. */
#include "statistics.h"
/* Local headers. */
#include "cooling.h"
#include "engine.h"
#include "error.h"
#include "hydro.h"
#include "threadpool.h"
/**
* @brief Information required to compute the statistics
*/
struct index_data {
const struct space *s;
int chunk_size;
int num_threads;
struct statistics *stats;
};
void stats_collect_part_mapper(void *map_data, int nr_parts, void *extra_data) {
/* Unpack the data */
struct index_data *data = (struct index_data *)extra_data;
const struct space *s = data->s;
struct part *restrict parts = (struct part *)map_data;
struct xpart *restrict xparts = s->xparts + (ptrdiff_t)(parts - s->parts);
const int chunk_size = data->chunk_size;
const int num_threads = data->num_threads;
const int ti_current = s->e->ti_current;
const double timeBase = s->e->timeBase;
/* Local accumulator */
struct statistics stats;
bzero(&stats, sizeof(struct statistics));
/* Loop over particles */
for (int k = 0; k < nr_parts; k++) {
/* Get the particle */
struct part *restrict p = &parts[k];
struct xpart *restrict xp = &xparts[k];
const float dt = (ti_current - (p->ti_begin + p->ti_end) / 2) * timeBase;
const double x[3] = {p->x[0], p->x[1], p->x[2]};
const float v[3] = {xp->v_full[0] + p->a_hydro[0] * dt,
xp->v_full[1] + p->a_hydro[1] * dt,
xp->v_full[2] + p->a_hydro[2] * dt};
const float m = hydro_get_mass(p);
/* Collect mass */
stats.mass += m;
/* Collect momentum */
stats.mom[0] += m * v[0];
stats.mom[1] += m * v[1];
stats.mom[2] += m * v[2];
/* Collect angular momentum */
stats.ang_mom[0] += m * (x[1] * v[2] - x[2] * v[1]);
stats.ang_mom[1] += m * (x[2] * v[0] - x[0] * v[2]);
stats.ang_mom[2] += m * (x[0] * v[1] - x[1] * v[0]);
/* Collect energies. */
stats.E_kin += 0.5f * m * (v[0] * v[0] + v[1] * v[1] + v[2] * v[2]);
stats.E_pot += 0.;
stats.E_int += m * hydro_get_internal_energy(p, dt);
stats.E_rad += cooling_get_radiated_energy(xp);
/* Collect entropy */
stats.entropy += m * hydro_get_entropy(p, dt);
}
/* Where do we write back ? */
const int thread_id =
((ptrdiff_t)(parts - s->parts) / chunk_size) % num_threads;
/* Now write back to memory */
data->stats[thread_id].E_kin += stats.E_kin;
data->stats[thread_id].E_int += stats.E_int;
data->stats[thread_id].E_pot += stats.E_pot;
data->stats[thread_id].E_rad += stats.E_rad;
data->stats[thread_id].entropy += stats.entropy;
data->stats[thread_id].mass += stats.mass;
data->stats[thread_id].mom[0] += stats.mom[0];
data->stats[thread_id].mom[1] += stats.mom[1];
data->stats[thread_id].mom[2] += stats.mom[2];
data->stats[thread_id].ang_mom[0] += stats.ang_mom[0];
data->stats[thread_id].ang_mom[1] += stats.ang_mom[1];
data->stats[thread_id].ang_mom[2] += stats.ang_mom[2];
}
void stats_collect(const struct space *s, struct statistics *stats) {
const int num_threads = s->e->threadpool.num_threads;
const int chunk_size = s->nr_parts / num_threads;
/* Prepare the data */
struct index_data extra_data;
extra_data.s = s;
extra_data.chunk_size = chunk_size;
extra_data.num_threads = num_threads;
extra_data.stats = malloc(sizeof(struct statistics) * num_threads);
if (extra_data.stats == NULL)
error("Impossible to allocate memory for statistics");
bzero(&extra_data.stats, sizeof(struct statistics) * num_threads);
/* Run parallel collection of statistics */
threadpool_map(&s->e->threadpool, stats_collect_part_mapper, s->parts,
s->nr_parts, sizeof(struct part), chunk_size, &extra_data);
/* Collect data from all the threads in the pool */
for (int i = 0; i < num_threads; ++i) {
stats->E_kin += extra_data.stats[i].E_kin;
stats->E_int += extra_data.stats[i].E_int;
stats->E_pot += extra_data.stats[i].E_pot;
stats->E_rad += extra_data.stats[i].E_rad;
stats->entropy += extra_data.stats[i].entropy;
stats->mass += extra_data.stats[i].mass;
stats->mom[0] += extra_data.stats[i].mom[0];
stats->mom[1] += extra_data.stats[i].mom[1];
stats->mom[2] += extra_data.stats[i].mom[2];
stats->ang_mom[0] += extra_data.stats[i].ang_mom[0];
stats->ang_mom[1] += extra_data.stats[i].ang_mom[1];
stats->ang_mom[2] += extra_data.stats[i].ang_mom[2];
}
/* Free stuff */
free(extra_data.stats);
}
/*******************************************************************************
* This file is part of SWIFT.
* Copyright (c) 2016 Matthieu Schaller (matthieu.schaller@durham.ac.uk)
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU Lesser General Public License as published
* by the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* This program is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU Lesser General Public License
* along with this program. If not, see <http://www.gnu.org/licenses/>.
*
******************************************************************************/
#ifndef SWIFT_STATISTICS_H
#define SWIFT_STATISTICS_H
/* Config parameters. */
#include "../config.h"
/* Local headers. */
#include "space.h"
/**
* @brief Quantities collected for physics statistics
*/
struct statistics {
/*! Kinetic energy */
double E_kin;
/*! Internal energy */
double E_int;
/*! Potential energy */
double E_pot;
/*! Radiative energy */
double E_rad;
/*! Entropy */
double entropy;
/*! Mass */
double mass;
/*! Momentum */
double mom[3];
/*! Angular momentum */
double ang_mom[3];
};
void stats_collect(const struct space* s, struct statistics* stats);
#endif /* SWIFT_STATISTICS_H */
Supports Markdown
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment