Skip to content
Snippets Groups Projects
Commit 7bd5c53a authored by Matthieu Schaller's avatar Matthieu Schaller
Browse files

Added files for minimal implementation of SPH

parent 97c7c027
No related branches found
No related tags found
2 merge requests!136Master,!90Improved multi-timestep SPH
......@@ -47,6 +47,8 @@ AM_SOURCES = space.c runner.c queue.c task.c cell.c engine.c \
noinst_HEADERS = approx_math.h atomic.h cycle.h error.h inline.h kernel.h vector.h \
runner_doiact.h runner_doiact_grav.h units.h intrinsics.h \
hydro.h hydro_io.h gravity.h \
hydro/Minimal/hydro.h hydro/Minimal/hydro_iact.h hydro/Minimal/hydro_io.h \
hydro/Minimal/hydro_debug.h hydro/Minimal/hydro_part.h
hydro/Default/hydro.h hydro/Default/hydro_iact.h hydro/Default/hydro_io.h \
hydro/Default/hydro_debug.h hydro/Default/hydro_part.h \
hydro/Gadget2/hydro.h hydro/Gadget2/hydro_iact.h hydro/Gadget2/hydro_io.h \
......
......@@ -66,7 +66,9 @@
#define const_iepsilon6 (const_iepsilon3* const_iepsilon3)
/* SPH variant to use */
#define GADGET2_SPH
#define MINIMAL_SPH
//#define GADGET2_SPH
//#define DEFAULT_SPH
/* System of units */
#define const_unit_length_in_cgs 1 /* 3.08567810e16 /\* 1Mpc *\/ */
......
/*******************************************************************************
* This file is part of SWIFT.
* Coypright (c) 2015 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/>.
*
******************************************************************************/
/**
* @brief Computes the hydro time-step of a given particle
*
* @param p Pointer to the particle data
* @param xp Pointer to the extended particle data
*
*/
__attribute__((always_inline)) INLINE static float hydro_compute_timestep(
struct part* p, struct xpart* xp) {
/* CFL condition */
float dt_cfl = 2.f * const_cfl * kernel_gamma * p->h / p->force.v_sig;
/* Limit change in u */
float dt_u_change = (p->force.u_dt != 0.0f)
? fabsf(const_max_u_change * p->u / p->force.u_dt)
: FLT_MAX;
return fminf(dt_cfl, fminf(dt_h_change, dt_u_change));
}
/**
* @brief Prepares a particle for the density calculation.
*
* Zeroes all the relevant arrays in preparation for the sums taking place in
* the variaous density tasks
*
* @param p The particle to act upon
*/
__attribute__((always_inline))
INLINE static void hydro_init_part(struct part* p) {
p->density.wcount = 0.f;
p->density.wcount_dh = 0.f;
p->rho = 0.f;
p->rho_dh = 0.f;
p->density.div_v = 0.f;
p->density.curl_v[0] = 0.f;
p->density.curl_v[1] = 0.f;
p->density.curl_v[2] = 0.f;
}
/**
* @brief Finishes the density calculation.
*
* Multiplies the density and number of neighbours by the appropiate constants
* and add the self-contribution term.
*
* @param p The particle to act upon
* @param time The current time
*/
__attribute__((always_inline))
INLINE static void hydro_end_density(struct part* p, float time) {
/* Some smoothing length multiples. */
const float h = p->h;
const float ih = 1.0f / h;
const float ih2 = ih * ih;
const float ih4 = ih2 * ih2;
/* Final operation on the density. */
p->rho = ih * ih2 * (p->rho + p->mass * kernel_root);
p->rho_dh = (p->rho_dh - 3.0f * p->mass * kernel_root) * ih4;
p->density.wcount =
(p->density.wcount + kernel_root) * (4.0f / 3.0 * M_PI * kernel_gamma3);
p->density.wcount_dh =
p->density.wcount_dh * ih * (4.0f / 3.0 * M_PI * kernel_gamma3);
}
/**
* @brief Prepare a particle for the force calculation.
*
* Computes viscosity term, conduction term and smoothing length gradient terms.
*
* @param p The particle to act upon
* @param xp The extended particle data to act upon
* @param time The current time
*/
__attribute__((always_inline))
INLINE static void hydro_prepare_force(struct part* p, struct xpart* xp, float time) {
/* Some smoothing length multiples. */
const float h = p->h;
const float ih = 1.0f / h;
const float ih2 = ih * ih;
const float ih4 = ih2 * ih2;
/* Pre-compute some stuff for the balsara switch. */
const float normDiv_v = fabs(p->density.div_v / p->rho * ih4);
const float normCurl_v = sqrtf(p->density.curl_v[0] * p->density.curl_v[0] +
p->density.curl_v[1] * p->density.curl_v[1] +
p->density.curl_v[2] * p->density.curl_v[2]) /
p->rho * ih4;
/* Compute this particle's sound speed. */
const float u = p->u;
const float fc = p->force.c =
sqrtf(const_hydro_gamma * (const_hydro_gamma - 1.0f) * u);
/* Compute the P/Omega/rho2. */
xp->omega = 1.0f + 0.3333333333f * h * p->rho_dh / p->rho;
p->force.POrho2 = u * (const_hydro_gamma - 1.0f) / (p->rho * xp->omega);
/* Balsara switch */
p->force.balsara = normDiv_v / (normDiv_v + normCurl_v + 0.0001f * fc * ih);
/* Viscosity parameter decay time */
const float tau = h / (2.f * const_viscosity_length * p->force.c);
/* Viscosity source term */
const float S = fmaxf(-normDiv_v, 0.f);
/* Compute the particle's viscosity parameter time derivative */
const float alpha_dot = (const_viscosity_alpha_min - p->alpha) / tau +
(const_viscosity_alpha_max - p->alpha) * S;
/* Update particle's viscosity paramter */
p->alpha += alpha_dot * (p->t_end - p->t_begin);
}
/**
* @brief Reset acceleration fields of a particle
*
* Resets all hydro acceleration and time derivative fields in preparation
* for the sums taking place in the variaous force tasks
*
* @param p The particle to act upon
*/
__attribute__((always_inline))
INLINE static void hydro_reset_acceleration(struct part* p) {
/* Reset the acceleration. */
p->a_hydro[0] = 0.0f;
p->a_hydro[1] = 0.0f;
p->a_hydro[2] = 0.0f;
/* Reset the time derivatives. */
p->force.u_dt = 0.0f;
p->h_dt = 0.0f;
p->force.v_sig = 0.0f;
}
/**
* @brief Predict additional particle fields forward in time when drifting
*
* @param p The particle
* @param xp The extended data of the particle
* @param t0 The time at the start of the drift
* @param t1 The time at the end of the drift
*/
__attribute__((always_inline)) INLINE static void hydro_predict_extra(
struct part* p, struct xpart* xp, float t0, float t1) {
float u, w;
const float dt = t1 - t0;
/* Predict internal energy */
w = p->force.u_dt / p->u * dt;
if (fabsf(w) < 0.01f) /* 1st order expansion of exp(w) */
u = p->u *=
1.0f + w * (1.0f + w * (0.5f + w * (1.0f / 6.0f + 1.0f / 24.0f * w)));
else
u = p->u *= expf(w);
/* Predict gradient term */
p->force.POrho2 = u * (const_hydro_gamma - 1.0f) / (p->rho * xp->omega);
}
/**
* @brief Finishes the force calculation.
*
* Multiplies the forces and accelerationsby the appropiate constants
*
* @param p The particle to act upon
*/
__attribute__((always_inline))
INLINE static void hydro_end_force(struct part* p) {}
/**
* @brief Kick the additional variables
*
* @param p The particle to act upon
* @param dt The time-step for this kick
*/
__attribute__((always_inline))
INLINE static void hydro_kick_extra(struct part* p, float dt) {}
/**
* @brief Converts hydro quantity of a particle
*
* Requires the density to be known
*
* @param p The particle to act upon
*/
__attribute__((always_inline))
INLINE static void hydro_convert_quantities(struct part* p) {}
/*******************************************************************************
* This file is part of SWIFT.
* Coypright (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/>.
*
******************************************************************************/
__attribute__((always_inline))
INLINE static void hydro_debug_particle(struct part* p, struct xpart* xp) {
printf(
"x=[%.3e,%.3e,%.3e], "
"v=[%.3e,%.3e,%.3e],v_full=[%.3e,%.3e,%.3e] \n a=[%.3e,%.3e,%.3e],\n "
"h=%.3e, "
"wcount=%d, m=%.3e, dh_drho=%.3e, rho=%.3e, t_begin=%.3e, t_end=%.3e\n",
p->x[0], p->x[1], p->x[2], p->v[0], p->v[1], p->v[2], xp->v_full[0],
xp->v_full[1], xp->v_full[2], p->a_hydro[0], p->a_hydro[1], p->a_hydro[2],
p->h, (int)p->density.wcount, p->mass, p->rho_dh, p->rho, p->t_begin,
p->t_end);
}
This diff is collapsed.
/*******************************************************************************
* This file is part of SWIFT.
* Coypright (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/>.
*
******************************************************************************/
/**
* @brief Reads the different particles to the HDF5 file
*
* @param h_grp The HDF5 group in which to read the arrays.
* @param N The number of particles on that MPI rank.
* @param N_total The total number of particles (only used in MPI mode)
* @param offset The offset of the particles for this MPI rank (only used in MPI
*mode)
* @param parts The particle array
*
*/
__attribute__((always_inline)) INLINE static void hydro_read_particles(
hid_t h_grp, int N, long long N_total, long long offset,
struct part* parts) {
/* Read arrays */
readArray(h_grp, "Coordinates", DOUBLE, N, 3, parts, N_total, offset, x,
COMPULSORY);
readArray(h_grp, "Velocities", FLOAT, N, 3, parts, N_total, offset, v,
COMPULSORY);
readArray(h_grp, "Masses", FLOAT, N, 1, parts, N_total, offset, mass,
COMPULSORY);
readArray(h_grp, "SmoothingLength", FLOAT, N, 1, parts, N_total, offset, h,
COMPULSORY);
readArray(h_grp, "InternalEnergy", FLOAT, N, 1, parts, N_total, offset, u,
COMPULSORY);
readArray(h_grp, "ParticleIDs", ULONGLONG, N, 1, parts, N_total, offset, id,
COMPULSORY);
readArray(h_grp, "Acceleration", FLOAT, N, 3, parts, N_total, offset, a_hydro,
OPTIONAL);
readArray(h_grp, "Density", FLOAT, N, 1, parts, N_total, offset, rho,
OPTIONAL);
}
/**
* @brief Writes the different particles to the HDF5 file
*
* @param h_grp The HDF5 group in which to write the arrays.
* @param fileName The name of the file (unsued in MPI mode).
* @param xmfFile The XMF file to write to (unused in MPI mode).
* @param N The number of particles on that MPI rank.
* @param N_total The total number of particles (only used in MPI mode)
* @param mpi_rank The MPI rank of this node (only used in MPI mode)
* @param offset The offset of the particles for this MPI rank (only used in MPI
*mode)
* @param parts The particle array
* @param us The unit system to use
*
*/
__attribute__((always_inline)) INLINE static void hydro_write_particles(
hid_t h_grp, char* fileName, FILE* xmfFile, int N, long long N_total,
int mpi_rank, long long offset, struct part* parts, struct UnitSystem* us) {
/* Write arrays */
writeArray(h_grp, fileName, xmfFile, "Coordinates", DOUBLE, N, 3, parts,
N_total, mpi_rank, offset, x, us, UNIT_CONV_LENGTH);
writeArray(h_grp, fileName, xmfFile, "Velocities", FLOAT, N, 3, parts,
N_total, mpi_rank, offset, v, us, UNIT_CONV_SPEED);
writeArray(h_grp, fileName, xmfFile, "Masses", FLOAT, N, 1, parts, N_total,
mpi_rank, offset, mass, us, UNIT_CONV_MASS);
writeArray(h_grp, fileName, xmfFile, "SmoothingLength", FLOAT, N, 1, parts,
N_total, mpi_rank, offset, h, us, UNIT_CONV_LENGTH);
writeArray(h_grp, fileName, xmfFile, "InternalEnergy", FLOAT, N, 1, parts,
N_total, mpi_rank, offset, u, us, UNIT_CONV_ENERGY_PER_UNIT_MASS);
writeArray(h_grp, fileName, xmfFile, "ParticleIDs", ULONGLONG, N, 1, parts,
N_total, mpi_rank, offset, id, us, UNIT_CONV_NO_UNITS);
writeArray(h_grp, fileName, xmfFile, "Acceleration", FLOAT, N, 3, parts,
N_total, mpi_rank, offset, a_hydro, us, UNIT_CONV_ACCELERATION);
writeArray(h_grp, fileName, xmfFile, "Density", FLOAT, N, 1, parts, N_total,
mpi_rank, offset, rho, us, UNIT_CONV_DENSITY);
}
/*******************************************************************************
* This file is part of SWIFT.
* Copyright (c) 2012 Pedro Gonnet (pedro.gonnet@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/>.
*
******************************************************************************/
/* Extra particle data not needed during the SPH loops over neighbours. */
struct xpart {
/* Old position, at last tree rebuild. */
double x_old[3];
/* Velocity at the last full step. */
float v_full[3];
/* Old density. */
float omega;
} __attribute__((aligned(xpart_align)));
/* Data of a single particle. */
struct part {
/* Particle position. */
double x[3];
/* Particle predicted velocity. */
float v[3];
/* Particle acceleration. */
float a_hydro[3];
/* Particle cutoff radius. */
float h;
/* Change in smoothing length over time. */
float h_dt;
/* Particle time of beginning of time-step. */
float t_begin;
/* Particle time of end of time-step. */
float t_end;
/* Particle internal energy. */
float u;
/* Particle density. */
float rho;
/* Derivative of the density with respect to this particle's smoothing length.
*/
float rho_dh;
/* Particle viscosity parameter */
float alpha;
/* Store density/force specific stuff. */
// union {
struct {
/* Particle velocity divergence. */
float div_v;
/* Derivative of particle number density. */
float wcount_dh;
/* Particle velocity curl. */
float curl_v[3];
/* Particle number density. */
float wcount;
} density;
struct {
/* Balsara switch */
float balsara;
/* Aggregate quantities. */
float POrho2;
/* Change in particle energy over time. */
float u_dt;
/* Signal velocity */
float v_sig;
/* Sound speed */
float c;
} force;
//};
/* Particle mass. */
float mass;
/* Particle ID. */
unsigned long long id;
/* Pointer to corresponding gravity part. */
struct gpart* gpart;
} __attribute__((aligned(part_align)));
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment