![]() |
Shamrock 2025.10.0
Astrophysical Code
|
The GSPH Solver class. More...
#include <shammodels/gsph/include/shammodels/gsph/Solver.hpp>
Collaboration diagram for shammodels::gsph::Solver< Tvec, SPHKernel >:Public Types | |
| using | Tscal = shambase::VecComponent< Tvec > |
| using | Kernel = SPHKernel< Tscal > |
| using | Config = SolverConfig< Tvec, SPHKernel > |
| using | u_morton = u32 |
| using | GhostHandle = GSPHGhostHandler< Tvec > |
| using | GhostHandleCache = typename GhostHandle::CacheMap |
| using | RTree = typename Config::RTree |
Public Member Functions | |
| PatchScheduler & | scheduler () |
| void | init_required_fields () |
| void | gen_serial_patch_tree () |
| void | reset_serial_patch_tree () |
| void | gen_ghost_handler (Tscal time_val) |
| void | reset_ghost_handler () |
| void | build_ghost_cache () |
| void | clear_ghost_cache () |
| void | merge_position_ghost () |
| void | build_merged_pos_trees () |
| void | clear_merged_pos_trees () |
| void | compute_presteps_rint () |
| void | reset_presteps_rint () |
| void | start_neighbors_cache () |
| void | reset_neighbors_cache () |
| void | gsph_prestep (Tscal time_val, Tscal dt) |
| void | apply_position_boundary (Tscal time_val) |
| void | do_predictor_leapfrog (Tscal dt) |
| void | init_ghost_layout () |
| void | communicate_merge_ghosts_fields () |
| void | reset_merge_ghosts_fields () |
| void | compute_omega () |
| void | compute_eos_fields () |
| void | reset_eos_fields () |
| void | copy_eos_to_patchdata () |
| Copy EOS fields from solvergraph to patchdata for persistence. | |
| void | compute_gradients () |
| Compute gradients for MUSCL reconstruction. | |
| void | prepare_corrector () |
| void | update_derivs () |
| Update derivatives using GSPH Riemann solver. | |
| Tscal | compute_dt_cfl () |
| Compute CFL timestep constraint. | |
| bool | apply_corrector (Tscal dt, u64 Npart_all) |
| void | update_sync_load_values () |
| Solver (ShamrockCtx &context) | |
| void | init_solver_graph () |
| void | vtk_do_dump (std::string filename, bool add_patch_world_id) |
| void | print_timestep_logs () |
| TimestepLog | evolve_once () |
| Tscal | evolve_once_time_expl (Tscal t_current, Tscal dt_input) |
| bool | evolve_until (Tscal target_time, i32 niter_max=-1) |
Public Attributes | |
| ShamrockCtx & | context |
| SolverStorage< Tvec, u_morton > | storage {} |
| Config | solver_config |
| sph::SolverLog | solve_logs |
Static Public Attributes | |
| static constexpr u32 | dim = shambase::VectorProperties<Tvec>::dimension |
| static constexpr Tscal | Rkern = Kernel::Rkern |
The GSPH Solver class.
Implements the Godunov SPH method using Riemann solvers at particle interfaces instead of artificial viscosity.
| Tvec | Vector type (e.g., f64_3) |
| SPHKernel | Kernel type (e.g., M4, M6, C2, C4, C6) |
Definition at line 68 of file Solver.hpp.
| using shammodels::gsph::Solver< Tvec, SPHKernel >::Config = SolverConfig<Tvec, SPHKernel> |
Definition at line 74 of file Solver.hpp.
| using shammodels::gsph::Solver< Tvec, SPHKernel >::GhostHandle = GSPHGhostHandler<Tvec> |
Definition at line 95 of file Solver.hpp.
| using shammodels::gsph::Solver< Tvec, SPHKernel >::GhostHandleCache = typename GhostHandle::CacheMap |
Definition at line 96 of file Solver.hpp.
| using shammodels::gsph::Solver< Tvec, SPHKernel >::Kernel = SPHKernel<Tscal> |
Definition at line 72 of file Solver.hpp.
| using shammodels::gsph::Solver< Tvec, SPHKernel >::RTree = typename Config::RTree |
Definition at line 109 of file Solver.hpp.
| using shammodels::gsph::Solver< Tvec, SPHKernel >::Tscal = shambase::VecComponent<Tvec> |
Definition at line 70 of file Solver.hpp.
| using shammodels::gsph::Solver< Tvec, SPHKernel >::u_morton = u32 |
Definition at line 76 of file Solver.hpp.
|
inline |
Definition at line 180 of file Solver.hpp.
| bool shammodels::gsph::Solver< Tvec, Kern >::apply_corrector | ( | Tscal | dt, |
| u64 | Npart_all | ||
| ) |
Definition at line 1645 of file Solver.cpp.
| void shammodels::gsph::Solver< Tvec, Kern >::apply_position_boundary | ( | Tscal | time_val | ) |
Definition at line 516 of file Solver.cpp.
| void shammodels::gsph::Solver< Tvec, Kern >::build_ghost_cache | ( | ) |
Definition at line 185 of file Solver.cpp.
| void shammodels::gsph::Solver< Tvec, Kern >::build_merged_pos_trees | ( | ) |
Definition at line 249 of file Solver.cpp.
| void shammodels::gsph::Solver< Tvec, Kern >::clear_ghost_cache | ( | ) |
Definition at line 198 of file Solver.cpp.
| void shammodels::gsph::Solver< Tvec, Kern >::clear_merged_pos_trees | ( | ) |
Definition at line 288 of file Solver.cpp.
| void shammodels::gsph::Solver< Tvec, Kern >::communicate_merge_ghosts_fields | ( | ) |
Definition at line 643 of file Solver.cpp.
| shammodels::gsph::Solver< Tvec, Kern >::Tscal shammodels::gsph::Solver< Tvec, Kern >::compute_dt_cfl | ( | ) |
Compute CFL timestep constraint.
Computes timestep from:
Definition at line 1541 of file Solver.cpp.
Here is the call graph for this function:| void shammodels::gsph::Solver< Tvec, Kern >::compute_eos_fields | ( | ) |
Definition at line 1116 of file Solver.cpp.
| void shammodels::gsph::Solver< Tvec, Kern >::compute_gradients | ( | ) |
Compute gradients for MUSCL reconstruction.
Computes density, pressure, and velocity gradients for each particle using SPH kernel gradient summation. Only called when MUSCL reconstruction is enabled (reconstruct_config.is_muscl() == true).
Reference: Cha & Whitworth (2003)
Definition at line 1286 of file Solver.cpp.
Here is the call graph for this function:| void shammodels::gsph::Solver< Tvec, Kern >::compute_omega | ( | ) |
Definition at line 827 of file Solver.cpp.
| void shammodels::gsph::Solver< Tvec, Kern >::compute_presteps_rint | ( | ) |
Definition at line 294 of file Solver.cpp.
| void shammodels::gsph::Solver< Tvec, Kern >::copy_eos_to_patchdata | ( | ) |
Copy EOS fields from solvergraph to patchdata for persistence.
Copies density, pressure, and soundspeed from the temporary solvergraph fields to the persistent patchdata layout. This ensures thermodynamic state is preserved across simulation restarts and available for VTK output.
Definition at line 1224 of file Solver.cpp.
Here is the call graph for this function:| void shammodels::gsph::Solver< Tvec, Kern >::do_predictor_leapfrog | ( | Tscal | dt | ) |
Definition at line 563 of file Solver.cpp.
| shammodels::gsph::TimestepLog shammodels::gsph::Solver< Tvec, Kern >::evolve_once | ( | ) |
patch_rank_owner is automatically updated since it is just a lambda
Definition at line 1716 of file Solver.cpp.
Here is the call graph for this function:
|
inline |
Definition at line 197 of file Solver.hpp.
|
inline |
Definition at line 204 of file Solver.hpp.
| void shammodels::gsph::Solver< Tvec, Kern >::gen_ghost_handler | ( | Tscal | time_val | ) |
Definition at line 135 of file Solver.cpp.
| void shammodels::gsph::Solver< Tvec, Kern >::gen_serial_patch_tree | ( | ) |
Definition at line 126 of file Solver.cpp.
| void shammodels::gsph::Solver< Tvec, Kern >::gsph_prestep | ( | Tscal | time_val, |
| Tscal | dt | ||
| ) |
Definition at line 509 of file Solver.cpp.
| void shammodels::gsph::Solver< Tvec, Kern >::init_ghost_layout | ( | ) |
Definition at line 625 of file Solver.cpp.
|
inline |
Definition at line 88 of file Solver.hpp.
| void shammodels::gsph::Solver< Tvec, Kern >::init_solver_graph | ( | ) |
Definition at line 70 of file Solver.cpp.
| void shammodels::gsph::Solver< Tvec, Kern >::merge_position_ghost | ( | ) |
Definition at line 204 of file Solver.cpp.
| void shammodels::gsph::Solver< Tvec, Kern >::prepare_corrector | ( | ) |
Definition at line 1469 of file Solver.cpp.
|
inline |
Definition at line 186 of file Solver.hpp.
| void shammodels::gsph::Solver< Tvec, Kern >::reset_eos_fields | ( | ) |
Definition at line 1219 of file Solver.cpp.
|
inline |
Definition at line 99 of file Solver.hpp.
| void shammodels::gsph::Solver< Tvec, Kern >::reset_merge_ghosts_fields | ( | ) |
Definition at line 822 of file Solver.cpp.
| void shammodels::gsph::Solver< Tvec, Kern >::reset_neighbors_cache | ( | ) |
Definition at line 504 of file Solver.cpp.
| void shammodels::gsph::Solver< Tvec, Kern >::reset_presteps_rint | ( | ) |
Definition at line 328 of file Solver.cpp.
|
inline |
Definition at line 92 of file Solver.hpp.
|
inline |
Definition at line 81 of file Solver.hpp.
| void shammodels::gsph::Solver< Tvec, Kern >::start_neighbors_cache | ( | ) |
Definition at line 333 of file Solver.cpp.
| void shammodels::gsph::Solver< Tvec, Kern >::update_derivs | ( | ) |
Update derivatives using GSPH Riemann solver.
This is the core GSPH step: for each particle pair, solve the 1D Riemann problem and compute forces from the interface pressure p*.
Definition at line 1534 of file Solver.cpp.
Here is the call graph for this function:| void shammodels::gsph::Solver< Tvec, Kern >::update_sync_load_values | ( | ) |
Definition at line 1713 of file Solver.cpp.
| void shammodels::gsph::Solver< Tvec, Kern >::vtk_do_dump | ( | std::string | filename, |
| bool | add_patch_world_id | ||
| ) |
Definition at line 119 of file Solver.cpp.
| ShamrockCtx& shammodels::gsph::Solver< Tvec, SPHKernel >::context |
Definition at line 80 of file Solver.hpp.
|
staticconstexpr |
Definition at line 71 of file Solver.hpp.
|
staticconstexpr |
Definition at line 78 of file Solver.hpp.
| sph::SolverLog shammodels::gsph::Solver< Tvec, SPHKernel >::solve_logs |
Definition at line 86 of file Solver.hpp.
| Config shammodels::gsph::Solver< Tvec, SPHKernel >::solver_config |
Definition at line 85 of file Solver.hpp.
| SolverStorage<Tvec, u_morton> shammodels::gsph::Solver< Tvec, SPHKernel >::storage {} |
Definition at line 83 of file Solver.hpp.