Shamrock 2025.10.0
Astrophysical Code
Loading...
Searching...
No Matches
pyGSPHModel.cpp
Go to the documentation of this file.
1// -------------------------------------------------------//
2//
3// SHAMROCK code for hydrodynamics
4// Copyright (c) 2021-2026 Timothée David--Cléris <tim.shamrock@proton.me>
5// SPDX-License-Identifier: CeCILL Free Software License Agreement v2.1
6// Shamrock is licensed under the CeCILL 2.1 License, see LICENSE for more information
7//
8// -------------------------------------------------------//
9
27#include "shambase/memory.hpp"
35#include <pybind11/cast.h>
36#include <pybind11/numpy.h>
37#include <memory>
38
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;
42
43 using Tscal = shambase::VecComponent<Tvec>;
44
45 using T = Model<Tvec, SPHKernel>;
46 using TConfig = typename T::SolverConfig;
47
48 shamlog_debug_ln("[Py]", "registering class :", name_config, typeid(T).name());
49 shamlog_debug_ln("[Py]", "registering class :", name_model, typeid(T).name());
50
51 py::class_<TConfig> config_cls(m, name_config.c_str());
52
53 shammodels::common::add_json_defs<TConfig>(config_cls);
54
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)
58 // Riemann solver config
59 .def(
60 "set_riemann_iterative",
61 [](TConfig &self, Tscal tol, u32 max_iter) {
62 self.set_riemann_iterative(tol, max_iter);
63 },
64 py::kw_only(),
65 py::arg("tolerance") = Tscal{1e-6},
66 py::arg("max_iter") = 20,
67 R"==(
68 Set iterative Riemann solver (van Leer 1997).
69
70 This is the most accurate but slower Riemann solver.
71 Uses Newton-Raphson iteration to find the pressure in the star region.
72
73 Parameters
74 ----------
75 tolerance : float
76 Convergence tolerance for Newton-Raphson iteration (default: 1e-6)
77 max_iter : int
78 Maximum number of iterations (default: 20)
79)==")
80 .def(
81 "set_riemann_hllc",
82 [](TConfig &self) {
83 self.set_riemann_hllc();
84 },
85 R"==(
86 Set HLLC approximate Riemann solver.
87
88 Fast approximate Riemann solver that captures contact discontinuities.
89 Recommended for general use - good balance of accuracy and speed.
90)==")
91 // Reconstruction config
92 .def(
93 "set_reconstruct_piecewise_constant",
94 [](TConfig &self) {
95 self.set_reconstruct_piecewise_constant();
96 },
97 R"==(
98 Set first-order piecewise constant reconstruction.
99
100 Sets all gradients to zero. Most diffusive but most stable.
101 Good for very strong shocks or initial testing.
102)==")
103 // EOS config
104 .def(
105 "set_eos_adiabatic",
106 [](TConfig &self, Tscal gamma) {
107 self.set_eos_adiabatic(gamma);
108 },
109 py::arg("gamma"),
110 R"==(
111 Set adiabatic equation of state: P = (\gamma-1) \rho u
112
113 Parameters
114 ----------
115 gamma : float
116 Adiabatic index (e.g., 5/3 for monatomic gas, 7/5 for diatomic)
117)==")
118 .def(
119 "set_eos_isothermal",
120 [](TConfig &self, Tscal cs) {
121 self.set_eos_isothermal(cs);
122 },
123 py::arg("cs"),
124 R"==(
125 Set isothermal equation of state: P = cs^2 \rho
126
127 Parameters
128 ----------
129 cs : float
130 Sound speed
131)==")
132 // Boundary config
133 .def("set_boundary_free", &TConfig::set_boundary_free)
134 .def("set_boundary_periodic", &TConfig::set_boundary_periodic)
135 // External forces
136 .def(
137 "add_ext_force_point_mass",
138 [](TConfig &self, Tscal central_mass, Tscal Racc) {
139 self.add_ext_force_point_mass(central_mass, Racc);
140 },
141 py::kw_only(),
142 py::arg("central_mass"),
143 py::arg("Racc"))
144 // Units
145 .def("set_units", &TConfig::set_units)
146 // CFL
147 .def(
148 "set_cfl_cour",
149 [](TConfig &self, Tscal cfl_cour) {
150 self.cfl_config.cfl_cour = cfl_cour;
151 })
152 .def(
153 "set_cfl_force",
154 [](TConfig &self, Tscal cfl_force) {
155 self.cfl_config.cfl_force = cfl_force;
156 })
157 .def(
158 "set_particle_mass",
159 [](TConfig &self, Tscal gpart_mass) {
160 self.gpart_mass = gpart_mass;
161 })
162 .def(
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;
167 },
168 py::kw_only(),
169 py::arg("split_load_value"),
170 py::arg("merge_load_value"));
171
172 py::class_<T>(m, name_model.c_str())
173 .def(py::init([](ShamrockCtx &ctx) {
174 return std::make_unique<T>(ctx);
175 }))
176 .def("init", &T::init)
177 .def("init_scheduler", &T::init_scheduler)
178 .def("evolve_once", &T::evolve_once)
179 .def(
180 "evolve_until",
181 [](T &self, f64 target_time, i32 niter_max) {
182 return self.evolve_until(target_time, niter_max);
183 },
184 py::arg("target_time"),
185 py::kw_only(),
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)
194 .def(
195 "get_box_dim_fcc_3d",
196 [](T &self, f64 dr, u32 xcnt, u32 ycnt, u32 zcnt) {
197 return self.get_box_dim_fcc_3d(dr, xcnt, ycnt, zcnt);
198 })
199 .def(
200 "get_ideal_fcc_box",
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});
203 })
204 .def(
205 "get_ideal_hcp_box",
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});
208 })
209 .def(
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});
213 })
214 .def(
215 "add_cube_fcc_3d",
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});
218 })
219 .def(
220 "add_cube_hcp_3d",
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});
223 })
224 .def("get_total_part_count", &T::get_total_part_count)
225 .def("total_mass_to_part_mass", &T::total_mass_to_part_mass)
226 .def(
227 "set_field_in_box",
228 [](T &self,
229 std::string field_name,
230 std::string field_type,
231 pybind11::object value,
232 f64_3 box_min,
233 f64_3 box_max,
234 u32 ivar) {
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);
244 } else {
246 "unknown field type: " + field_type + ". Valid types: f64, f64_3, u32");
247 }
248 },
249 py::arg("field_name"),
250 py::arg("field_type"),
251 py::arg("value"),
252 py::arg("box_min"),
253 py::arg("box_max"),
254 py::kw_only(),
255 py::arg("ivar") = 0,
256 R"==(
257 Set field value for particles within a box region.
258
259 Useful for setting up discontinuous initial conditions like Sod shock tube.
260
261 Parameters
262 ----------
263 field_name : str
264 Name of the field to set (e.g., "vxyz", "uint", "hpart")
265 field_type : str
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)
269 box_min : tuple
270 Minimum corner of the box (x, y, z)
271 box_max : tuple
272 Maximum corner of the box (x, y, z)
273 ivar : int
274 Variable index for multi-component fields (default: 0)
275
276 Examples
277 --------
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))
282)==")
283 .def(
284 "set_field_in_sphere",
285 [](T &self,
286 std::string field_name,
287 std::string field_type,
288 pybind11::object value,
289 f64_3 center,
290 f64 radius) {
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);
297 } else {
299 "unknown field type");
300 }
301 },
302 py::arg("field_name"),
303 py::arg("field_type"),
304 py::arg("value"),
305 py::arg("center"),
306 py::arg("radius"),
307 R"==(
308 Set field value for particles within a spherical region.
309
310 Useful for setting up point-source initial conditions like Sedov blast.
311
312 Parameters
313 ----------
314 field_name : str
315 Name of the field to set (e.g., "uint")
316 field_type : str
317 Type of the field: "f64" or "f64_3"
318 value : float or tuple
319 Value to set (type must match field_type)
320 center : tuple
321 Center of the sphere (x, y, z)
322 radius : float
323 Radius of the sphere
324
325 Examples
326 --------
327 >>> # Sedov blast: inject energy in central sphere
328 >>> model.set_field_in_sphere("uint", "f64", u_blast, (0,0,0), r_blast)
329)==")
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>)
332 .def(
333 "get_sum",
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));
339 } else {
341 "unknown field type");
342 }
343 })
344 .def(
345 "gen_default_config",
346 [](T &self) {
347 return self.gen_default_config();
348 })
349 .def(
350 "get_current_config",
351 [](T &self) {
352 return self.solver.solver_config;
353 })
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)
358 .def(
359 "get_time",
360 [](T &self) {
361 return self.solver.solver_config.get_time();
362 })
363 .def(
364 "get_dt",
365 [](T &self) {
366 return self.solver.solver_config.get_dt();
367 })
368 .def(
369 "set_time",
370 [](T &self, Tscal t) {
371 return self.solver.solver_config.set_time(t);
372 })
373 .def(
374 "set_next_dt",
375 [](T &self, Tscal dt) {
376 return self.solver.solver_config.set_next_dt(dt);
377 })
378 .def(
379 "load_from_dump",
380 &T::load_from_dump,
381 py::arg("filename"),
382 R"==(
383 Load simulation state from a Shamrock dump file.
384
385 Uses the shared ShamrockDump mechanism (same as SPH).
386
387 Parameters
388 ----------
389 filename : str
390 Path to the dump file
391
392 Example
393 -------
394 >>> model.load_from_dump("checkpoint.shamrock")
395)==")
396 .def(
397 "dump",
398 &T::dump,
399 py::arg("filename"),
400 R"==(
401 Write simulation state to a Shamrock dump file.
402
403 Uses the shared ShamrockDump mechanism (same as SPH).
404
405 Parameters
406 ----------
407 filename : str
408 Path to the dump file
409
410 Example
411 -------
412 >>> model.dump("checkpoint.shamrock")
413)==");
414}
415
416using namespace shammodels::gsph;
417
418Register_pymod(pygsphmodel) {
419
420 py::module mgsph = m.def_submodule("model_gsph", "Shamrock GSPH (Godunov SPH) solver");
421
422 using namespace shammodels::gsph;
423
424 // Register GSPH models for different kernels
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");
431
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");
438
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>>>;
446
447 m.def(
448 "get_Model_GSPH",
449 [](ShamrockCtx &ctx, std::string vector_type, std::string kernel) -> VariantGSPHModelBind {
450 VariantGSPHModelBind ret;
451
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);
464 } else {
466 "unknown combination of representation and kernel");
467 }
468
469 return ret;
470 },
471 py::kw_only(),
472 py::arg("context"),
473 py::arg("vector_type") = "f64_3",
474 py::arg("sph_kernel") = "M4",
475 R"==(
476 Create a GSPH (Godunov SPH) model.
477
478 GSPH uses Riemann solvers at particle interfaces instead of artificial viscosity,
479 giving sharper shock resolution.
480
481 Parameters
482 ----------
483 context : ShamrockCtx
484 Shamrock context
485 vector_type : str
486 Vector type, e.g., "f64_3" for 3D double precision (default: "f64_3")
487 sph_kernel : str
488 SPH kernel type: "M4" (cubic spline, default), "M6", "M8" (quintic spline),
489 "C2", "C4", "C6" (Wendland kernels)
490
491 Returns
492 -------
493 GSPHModel
494 A GSPH model instance
495
496 Examples
497 --------
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)
504)==");
505}
MPI scheduler.
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
The GSPH Model class.
Definition Model.hpp:63
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...
sph kernels
Functions related to the MPI communicator.