35#include <pybind11/cast.h>
36#include <pybind11/numpy.h>
39template<
class Tvec,
template<
class>
class SPHKernel>
40void add_gsph_instance(py::module &m, std::string name_config, std::string name_model) {
41 using namespace shammodels::gsph;
43 using Tscal = shambase::VecComponent<Tvec>;
46 using TConfig =
typename T::SolverConfig;
48 shamlog_debug_ln(
"[Py]",
"registering class :", name_config,
typeid(T).name());
49 shamlog_debug_ln(
"[Py]",
"registering class :", name_model,
typeid(T).name());
51 py::class_<TConfig> config_cls(m, name_config.c_str());
53 shammodels::common::add_json_defs<TConfig>(config_cls);
55 config_cls.def(
"print_status", &TConfig::print_status)
56 .def(
"set_tree_reduction_level", &TConfig::set_tree_reduction_level)
57 .def(
"set_two_stage_search", &TConfig::set_two_stage_search)
60 "set_riemann_iterative",
61 [](TConfig &self, Tscal tol,
u32 max_iter) {
62 self.set_riemann_iterative(tol, max_iter);
65 py::arg(
"tolerance") = Tscal{1e-6},
66 py::arg(
"max_iter") = 20,
68 Set iterative Riemann solver (van Leer 1997).
70 This is the most accurate but slower Riemann solver.
71 Uses Newton-Raphson iteration to find the pressure in the star region.
76 Convergence tolerance for Newton-Raphson iteration (default: 1e-6)
78 Maximum number of iterations (default: 20)
83 self.set_riemann_hllc();
86 Set HLLC approximate Riemann solver.
88 Fast approximate Riemann solver that captures contact discontinuities.
89 Recommended for general use - good balance of accuracy and speed.
93 "set_reconstruct_piecewise_constant",
95 self.set_reconstruct_piecewise_constant();
98 Set first-order piecewise constant reconstruction.
100 Sets all gradients to zero. Most diffusive but most stable.
101 Good for very strong shocks or initial testing.
106 [](TConfig &self, Tscal gamma) {
107 self.set_eos_adiabatic(gamma);
111 Set adiabatic equation of state: P = (\gamma-1) \rho u
116 Adiabatic index (e.g., 5/3 for monatomic gas, 7/5 for diatomic)
119 "set_eos_isothermal",
120 [](TConfig &self, Tscal cs) {
121 self.set_eos_isothermal(cs);
125 Set isothermal equation of state: P = cs^2 \rho
133 .def(
"set_boundary_free", &TConfig::set_boundary_free)
134 .def(
"set_boundary_periodic", &TConfig::set_boundary_periodic)
137 "add_ext_force_point_mass",
138 [](TConfig &self, Tscal central_mass, Tscal Racc) {
139 self.add_ext_force_point_mass(central_mass, Racc);
142 py::arg(
"central_mass"),
145 .def(
"set_units", &TConfig::set_units)
149 [](TConfig &self, Tscal cfl_cour) {
150 self.cfl_config.cfl_cour = cfl_cour;
154 [](TConfig &self, Tscal cfl_force) {
155 self.cfl_config.cfl_force = cfl_force;
159 [](TConfig &self, Tscal gpart_mass) {
160 self.gpart_mass = gpart_mass;
163 "set_scheduler_config",
164 [](TConfig &self,
u64 split_crit,
u64 merge_crit) {
165 self.scheduler_conf.split_load_value = split_crit;
166 self.scheduler_conf.merge_load_value = merge_crit;
169 py::arg(
"split_load_value"),
170 py::arg(
"merge_load_value"));
172 py::class_<T>(m, name_model.c_str())
174 return std::make_unique<T>(ctx);
176 .def(
"init", &T::init)
177 .def(
"init_scheduler", &T::init_scheduler)
178 .def(
"evolve_once", &T::evolve_once)
181 [](T &self,
f64 target_time,
i32 niter_max) {
182 return self.evolve_until(target_time, niter_max);
184 py::arg(
"target_time"),
186 py::arg(
"niter_max") = -1)
187 .def(
"timestep", &T::timestep)
188 .def(
"set_cfl_cour", &T::set_cfl_cour, py::arg(
"cfl_cour"))
189 .def(
"set_cfl_force", &T::set_cfl_force, py::arg(
"cfl_force"))
190 .def(
"set_particle_mass", &T::set_particle_mass, py::arg(
"gpart_mass"))
191 .def(
"get_particle_mass", &T::get_particle_mass)
192 .def(
"rho_h", &T::rho_h)
193 .def(
"get_hfact", &T::get_hfact)
195 "get_box_dim_fcc_3d",
197 return self.get_box_dim_fcc_3d(dr, xcnt, ycnt, zcnt);
201 [](T &self,
f64 dr, f64_3 box_min, f64_3 box_max) {
202 return self.get_ideal_fcc_box(dr, {box_min, box_max});
206 [](T &self,
f64 dr, f64_3 box_min, f64_3 box_max) {
207 return self.get_ideal_hcp_box(dr, {box_min, box_max});
210 "resize_simulation_box",
211 [](T &self, f64_3 box_min, f64_3 box_max) {
212 return self.resize_simulation_box({box_min, box_max});
216 [](T &self,
f64 dr, f64_3 box_min, f64_3 box_max) {
217 return self.add_cube_fcc_3d(dr, {box_min, box_max});
221 [](T &self,
f64 dr, f64_3 box_min, f64_3 box_max) {
222 return self.add_cube_hcp_3d(dr, {box_min, box_max});
224 .def(
"get_total_part_count", &T::get_total_part_count)
225 .def(
"total_mass_to_part_mass", &T::total_mass_to_part_mass)
229 std::string field_name,
230 std::string field_type,
231 pybind11::object value,
235 if (field_type ==
"f64") {
236 f64 val = value.cast<
f64>();
237 self.set_field_in_box(field_name, val, {box_min, box_max}, ivar);
238 }
else if (field_type ==
"f64_3") {
239 f64_3 val = value.cast<f64_3>();
240 self.set_field_in_box(field_name, val, {box_min, box_max}, ivar);
241 }
else if (field_type ==
"u32") {
242 u32 val = value.cast<
u32>();
243 self.set_field_in_box(field_name, val, {box_min, box_max}, ivar);
246 "unknown field type: " + field_type +
". Valid types: f64, f64_3, u32");
249 py::arg(
"field_name"),
250 py::arg(
"field_type"),
257 Set field value for particles within a box region.
259 Useful for setting up discontinuous initial conditions like Sod shock tube.
264 Name of the field to set (e.g., "vxyz", "uint", "hpart")
266 Type of the field: "f64", "f64_3", or "u32"
267 value : float, tuple, or int
268 Value to set (type must match field_type)
270 Minimum corner of the box (x, y, z)
272 Maximum corner of the box (x, y, z)
274 Variable index for multi-component fields (default: 0)
278 >>> # Sod shock tube: set left state internal energy
279 >>> model.set_field_in_box("uint", "f64", u_left, (-1,-1,-1), (0,1,1))
280 >>> # Set right state
281 >>> model.set_field_in_box("uint", "f64", u_right, (0,-1,-1), (1,1,1))
284 "set_field_in_sphere",
286 std::string field_name,
287 std::string field_type,
288 pybind11::object value,
291 if (field_type ==
"f64") {
292 f64 val = value.cast<
f64>();
293 self.set_field_in_sphere(field_name, val, center, radius);
294 }
else if (field_type ==
"f64_3") {
295 f64_3 val = value.cast<f64_3>();
296 self.set_field_in_sphere(field_name, val, center, radius);
299 "unknown field type");
302 py::arg(
"field_name"),
303 py::arg(
"field_type"),
308 Set field value for particles within a spherical region.
310 Useful for setting up point-source initial conditions like Sedov blast.
315 Name of the field to set (e.g., "uint")
317 Type of the field: "f64" or "f64_3"
318 value : float or tuple
319 Value to set (type must match field_type)
321 Center of the sphere (x, y, z)
327 >>> # Sedov blast: inject energy in central sphere
328 >>> model.set_field_in_sphere("uint", "f64", u_blast, (0,0,0), r_blast)
330 .def("apply_field_from_position_f64_3", &T::template apply_field_from_position<f64_3>)
331 .def(
"apply_field_from_position_f64", &T::template apply_field_from_position<f64>)
334 [](T &self, std::string field_name, std::string field_type) {
335 if (field_type ==
"f64") {
336 return py::cast(self.template get_sum<f64>(field_name));
337 }
else if (field_type ==
"f64_3") {
338 return py::cast(self.template get_sum<f64_3>(field_name));
341 "unknown field type");
345 "gen_default_config",
347 return self.gen_default_config();
350 "get_current_config",
352 return self.solver.solver_config;
354 .def(
"set_solver_config", &T::set_solver_config)
355 .def(
"do_vtk_dump", &T::do_vtk_dump)
356 .def(
"solver_logs_last_rate", &T::solver_logs_last_rate)
357 .def(
"solver_logs_last_obj_count", &T::solver_logs_last_obj_count)
361 return self.solver.solver_config.get_time();
366 return self.solver.solver_config.get_dt();
370 [](T &self, Tscal t) {
371 return self.solver.solver_config.set_time(t);
375 [](T &self, Tscal dt) {
376 return self.solver.solver_config.set_next_dt(dt);
383 Load simulation state from a Shamrock dump file.
385 Uses the shared ShamrockDump mechanism (same as SPH).
390 Path to the dump file
394 >>> model.load_from_dump("checkpoint.shamrock")
401 Write simulation state to a Shamrock dump file.
403 Uses the shared ShamrockDump mechanism (same as SPH).
408 Path to the dump file
412 >>> model.dump("checkpoint.shamrock")
416using namespace shammodels::gsph;
420 py::module mgsph = m.def_submodule(
"model_gsph",
"Shamrock GSPH (Godunov SPH) solver");
422 using namespace shammodels::gsph;
425 add_gsph_instance<f64_3, shammath::M4>(
426 mgsph,
"GSPHModel_f64_3_M4_SolverConfig",
"GSPHModel_f64_3_M4");
427 add_gsph_instance<f64_3, shammath::M6>(
428 mgsph,
"GSPHModel_f64_3_M6_SolverConfig",
"GSPHModel_f64_3_M6");
429 add_gsph_instance<f64_3, shammath::M8>(
430 mgsph,
"GSPHModel_f64_3_M8_SolverConfig",
"GSPHModel_f64_3_M8");
432 add_gsph_instance<f64_3, shammath::C2>(
433 mgsph,
"GSPHModel_f64_3_C2_SolverConfig",
"GSPHModel_f64_3_C2");
434 add_gsph_instance<f64_3, shammath::C4>(
435 mgsph,
"GSPHModel_f64_3_C4_SolverConfig",
"GSPHModel_f64_3_C4");
436 add_gsph_instance<f64_3, shammath::C6>(
437 mgsph,
"GSPHModel_f64_3_C6_SolverConfig",
"GSPHModel_f64_3_C6");
439 using VariantGSPHModelBind = std::variant<
440 std::unique_ptr<Model<f64_3, shammath::M4>>,
441 std::unique_ptr<Model<f64_3, shammath::M6>>,
442 std::unique_ptr<Model<f64_3, shammath::M8>>,
443 std::unique_ptr<Model<f64_3, shammath::C2>>,
444 std::unique_ptr<Model<f64_3, shammath::C4>>,
445 std::unique_ptr<Model<f64_3, shammath::C6>>>;
449 [](
ShamrockCtx &ctx, std::string vector_type, std::string kernel) -> VariantGSPHModelBind {
450 VariantGSPHModelBind ret;
452 if (vector_type ==
"f64_3" && kernel ==
"M4") {
453 ret = std::make_unique<Model<f64_3, shammath::M4>>(ctx);
454 }
else if (vector_type ==
"f64_3" && kernel ==
"M6") {
455 ret = std::make_unique<Model<f64_3, shammath::M6>>(ctx);
456 }
else if (vector_type ==
"f64_3" && kernel ==
"M8") {
457 ret = std::make_unique<Model<f64_3, shammath::M8>>(ctx);
458 }
else if (vector_type ==
"f64_3" && kernel ==
"C2") {
459 ret = std::make_unique<Model<f64_3, shammath::C2>>(ctx);
460 }
else if (vector_type ==
"f64_3" && kernel ==
"C4") {
461 ret = std::make_unique<Model<f64_3, shammath::C4>>(ctx);
462 }
else if (vector_type ==
"f64_3" && kernel ==
"C6") {
463 ret = std::make_unique<Model<f64_3, shammath::C6>>(ctx);
466 "unknown combination of representation and kernel");
473 py::arg(
"vector_type") =
"f64_3",
474 py::arg(
"sph_kernel") =
"M4",
476 Create a GSPH (Godunov SPH) model.
478 GSPH uses Riemann solvers at particle interfaces instead of artificial viscosity,
479 giving sharper shock resolution.
483 context : ShamrockCtx
486 Vector type, e.g., "f64_3" for 3D double precision (default: "f64_3")
488 SPH kernel type: "M4" (cubic spline, default), "M6", "M8" (quintic spline),
489 "C2", "C4", "C6" (Wendland kernels)
494 A GSPH model instance
498 >>> ctx = shamrock.ShamrockCtx()
499 >>> model = shamrock.get_Model_GSPH(context=ctx) # Uses M4 kernel by default
500 >>> config = model.gen_default_config()
501 >>> config.set_riemann_hllc()
502 >>> config.set_eos_adiabatic(1.4)
503 >>> model.set_solver_config(config)
double f64
Alias for double.
std::uint32_t u32
32 bit unsigned integer
std::uint64_t u64
64 bit unsigned integer
std::int32_t i32
32 bit integer
This header file contains utility functions related to exception handling in the code.
GSPH Model class - high-level interface for GSPH simulations.
void throw_with_loc(std::string message, SourceLocation loc=SourceLocation{})
Throw an exception and append the source location to it.
Pybind11 include and definitions.
#define Register_pymod(placeholdername)
Register a python module init function using static initialisation.
Utilities to convert JSON objects to Python objects and vice versa. TODO: try to convert directly wit...
Functions related to the MPI communicator.