Shamrock 2025.10.0
Astrophysical Code
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Public Attributes | Static Public Attributes | List of all members
shammodels::gsph::Model< Tvec, SPHKernel > Class Template Reference

The GSPH Model class. More...

#include <shammodels/gsph/include/shammodels/gsph/Model.hpp>

+ Collaboration diagram for shammodels::gsph::Model< Tvec, SPHKernel >:

Public Types

using Tscal = shambase::VecComponent< Tvec >
 
using Kernel = SPHKernel< Tscal >
 
using Solver = Solver< Tvec, SPHKernel >
 
using SolverConfig = typename Solver::Config
 

Public Member Functions

 Model (ShamrockCtx &ctx)
 
void init ()
 Initialise the model and all the related data structures (patch scheduler in particular)
 
void init_scheduler (u32 crit_split, u32 crit_merge)
 
template<std::enable_if_t< dim==3, int > = 0>
Tvec get_box_dim_fcc_3d (Tscal dr, u32 xcnt, u32 ycnt, u32 zcnt)
 
void set_cfl_cour (Tscal cfl_cour)
 
void set_cfl_force (Tscal cfl_force)
 
void set_particle_mass (Tscal gpart_mass)
 
Tscal get_particle_mass ()
 
void resize_simulation_box (std::pair< Tvec, Tvec > box)
 
void do_vtk_dump (std::string filename, bool add_patch_world_id)
 
u64 get_total_part_count ()
 
f64 total_mass_to_part_mass (f64 totmass)
 
std::pair< Tvec, Tvec > get_ideal_fcc_box (Tscal dr, std::pair< Tvec, Tvec > box)
 
std::pair< Tvec, Tvec > get_ideal_hcp_box (Tscal dr, std::pair< Tvec, Tvec > box)
 
Tscal get_hfact ()
 
Tscal rho_h (Tscal h)
 
void add_cube_fcc_3d (Tscal dr, std::pair< Tvec, Tvec > _box)
 
void add_cube_hcp_3d (Tscal dr, std::pair< Tvec, Tvec > _box)
 
template<class T >
void apply_field_from_position (std::string field_name, const std::function< T(Tvec)> pos_to_val)
 Apply a position-dependent function to initialize a field.
 
template<class T >
void set_field_in_box (std::string field_name, T val, std::pair< Tvec, Tvec > box, u32 ivar=0)
 Set field value for particles within a box region.
 
template<class T >
void set_field_in_sphere (std::string field_name, T val, Tvec center, Tscal radius)
 Set field value for particles within a spherical region.
 
template<class T >
get_sum (std::string name)
 
SolverConfig gen_default_config ()
 
void set_solver_config (SolverConfig cfg)
 
f64 solver_logs_last_rate ()
 
u64 solver_logs_last_obj_count ()
 
shamsys::SystemMetrics solver_logs_last_system_metrics ()
 
void load_from_dump (std::string fname)
 
void dump (std::string fname)
 
TimestepLog timestep ()
 
void evolve_once ()
 
bool evolve_until (Tscal target_time, i32 niter_max=-1)
 

Public Attributes

ShamrockCtxctx
 
Solver solver
 

Static Public Attributes

static constexpr u32 dim = shambase::VectorProperties<Tvec>::dimension
 

Detailed Description

template<class Tvec, template< class > class SPHKernel>
class shammodels::gsph::Model< Tvec, SPHKernel >

The GSPH Model class.

Provides a high-level interface for setting up and running GSPH simulations. The GSPH method uses Riemann solvers at particle interfaces instead of artificial viscosity, giving sharper shock resolution.

Template Parameters
TvecVector type (e.g., f64_3)
SPHKernelKernel type (e.g., M4, M6, C2, C4, C6 for Wendland)

Definition at line 63 of file Model.hpp.

Member Typedef Documentation

◆ Kernel

template<class Tvec , template< class > class SPHKernel>
using shammodels::gsph::Model< Tvec, SPHKernel >::Kernel = SPHKernel<Tscal>

Definition at line 67 of file Model.hpp.

◆ Solver

template<class Tvec , template< class > class SPHKernel>
using shammodels::gsph::Model< Tvec, SPHKernel >::Solver = Solver<Tvec, SPHKernel>

Definition at line 69 of file Model.hpp.

◆ SolverConfig

template<class Tvec , template< class > class SPHKernel>
using shammodels::gsph::Model< Tvec, SPHKernel >::SolverConfig = typename Solver::Config

Definition at line 70 of file Model.hpp.

◆ Tscal

template<class Tvec , template< class > class SPHKernel>
using shammodels::gsph::Model< Tvec, SPHKernel >::Tscal = shambase::VecComponent<Tvec>

Definition at line 65 of file Model.hpp.

Constructor & Destructor Documentation

◆ Model()

template<class Tvec , template< class > class SPHKernel>
shammodels::gsph::Model< Tvec, SPHKernel >::Model ( ShamrockCtx ctx)
inline

Definition at line 75 of file Model.hpp.

Member Function Documentation

◆ add_cube_fcc_3d()

template<class Tvec , template< class > class SPHKernel>
void shammodels::gsph::Model< Tvec, SPHKernel >::add_cube_fcc_3d ( Tscal  dr,
std::pair< Tvec, Tvec >  _box 
)

Definition at line 99 of file Model.cpp.

◆ add_cube_hcp_3d()

template<class Tvec , template< class > class SPHKernel>
void shammodels::gsph::Model< Tvec, SPHKernel >::add_cube_hcp_3d ( Tscal  dr,
std::pair< Tvec, Tvec >  _box 
)

Definition at line 208 of file Model.cpp.

◆ apply_field_from_position()

template<class Tvec , template< class > class SPHKernel>
template<class T >
void shammodels::gsph::Model< Tvec, SPHKernel >::apply_field_from_position ( std::string  field_name,
const std::function< T(Tvec)>  pos_to_val 
)
inline

Apply a position-dependent function to initialize a field.

Sets field values by evaluating a function at each particle position. Useful for setting up spatially-varying initial conditions.

Template Parameters
TField type (e.g., Tscal for density, Tvec for velocity)
Parameters
field_nameName of the field to modify (e.g., "uint", "vxyz")
pos_to_valFunction mapping position to field value

Example:

// Set velocity as a function of position
model.apply_field_from_position<Tvec>("vxyz", [](Tvec pos) {
return Tvec{pos[0], 0.0, 0.0}; // Linear velocity profile
});

Definition at line 158 of file Model.hpp.

+ Here is the call graph for this function:

◆ do_vtk_dump()

template<class Tvec , template< class > class SPHKernel>
void shammodels::gsph::Model< Tvec, SPHKernel >::do_vtk_dump ( std::string  filename,
bool  add_patch_world_id 
)
inline

Definition at line 115 of file Model.hpp.

◆ dump()

template<class Tvec , template< class > class SPHKernel>
void shammodels::gsph::Model< Tvec, SPHKernel >::dump ( std::string  fname)
inline

Definition at line 377 of file Model.hpp.

◆ evolve_once()

template<class Tvec , template< class > class SPHKernel>
void shammodels::gsph::Model< Tvec, SPHKernel >::evolve_once ( )
inline

Definition at line 397 of file Model.hpp.

◆ evolve_until()

template<class Tvec , template< class > class SPHKernel>
bool shammodels::gsph::Model< Tvec, SPHKernel >::evolve_until ( Tscal  target_time,
i32  niter_max = -1 
)
inline

Definition at line 402 of file Model.hpp.

◆ gen_default_config()

template<class Tvec , template< class > class SPHKernel>
SolverConfig shammodels::gsph::Model< Tvec, SPHKernel >::gen_default_config ( )
inline

Definition at line 326 of file Model.hpp.

◆ get_box_dim_fcc_3d()

template<class Tvec , template< class > class SPHKernel>
template<std::enable_if_t< dim==3, int > = 0>
Tvec shammodels::gsph::Model< Tvec, SPHKernel >::get_box_dim_fcc_3d ( Tscal  dr,
u32  xcnt,
u32  ycnt,
u32  zcnt 
)
inline

Definition at line 93 of file Model.hpp.

◆ get_hfact()

template<class Tvec , template< class > class SPHKernel>
Tscal shammodels::gsph::Model< Tvec, SPHKernel >::get_hfact ( )
inline

Definition at line 126 of file Model.hpp.

◆ get_ideal_fcc_box()

template<class Tvec , template< class > class SPHKernel>
auto shammodels::gsph::Model< Tvec, SPHKernel >::get_ideal_fcc_box ( Tscal  dr,
std::pair< Tvec, Tvec >  box 
)

Definition at line 81 of file Model.cpp.

◆ get_ideal_hcp_box()

template<class Tvec , template< class > class SPHKernel>
auto shammodels::gsph::Model< Tvec, SPHKernel >::get_ideal_hcp_box ( Tscal  dr,
std::pair< Tvec, Tvec >  box 
)

Definition at line 90 of file Model.cpp.

◆ get_particle_mass()

template<class Tvec , template< class > class SPHKernel>
Tscal shammodels::gsph::Model< Tvec, SPHKernel >::get_particle_mass ( )
inline

Definition at line 109 of file Model.hpp.

◆ get_sum()

template<class Tvec , template< class > class SPHKernel>
template<class T >
T shammodels::gsph::Model< Tvec, SPHKernel >::get_sum ( std::string  name)
inline

Definition at line 306 of file Model.hpp.

◆ get_total_part_count()

template<class Tvec , template< class > class SPHKernel>
u64 shammodels::gsph::Model< Tvec, SPHKernel >::get_total_part_count ( )

Definition at line 70 of file Model.cpp.

◆ init()

template<class Tvec , template< class > class SPHKernel>
void shammodels::gsph::Model< Tvec, SPHKernel >::init ( )

Initialise the model and all the related data structures (patch scheduler in particular)

Definition at line 40 of file Model.cpp.

+ Here is the call graph for this function:

◆ init_scheduler()

template<class Tvec , template< class > class SPHKernel>
void shammodels::gsph::Model< Tvec, SPHKernel >::init_scheduler ( u32  crit_split,
u32  crit_merge 
)
inline

Old way of doing it, for backward compatibility it just overrides the values in the config before calling init()

Definition at line 86 of file Model.hpp.

+ Here is the call graph for this function:

◆ load_from_dump()

template<class Tvec , template< class > class SPHKernel>
void shammodels::gsph::Model< Tvec, SPHKernel >::load_from_dump ( std::string  fname)
inline

Definition at line 354 of file Model.hpp.

◆ resize_simulation_box()

template<class Tvec , template< class > class SPHKernel>
void shammodels::gsph::Model< Tvec, SPHKernel >::resize_simulation_box ( std::pair< Tvec, Tvec >  box)
inline

Definition at line 111 of file Model.hpp.

◆ rho_h()

template<class Tvec , template< class > class SPHKernel>
Tscal shammodels::gsph::Model< Tvec, SPHKernel >::rho_h ( Tscal  h)
inline

Definition at line 128 of file Model.hpp.

◆ set_cfl_cour()

template<class Tvec , template< class > class SPHKernel>
void shammodels::gsph::Model< Tvec, SPHKernel >::set_cfl_cour ( Tscal  cfl_cour)
inline

Definition at line 97 of file Model.hpp.

◆ set_cfl_force()

template<class Tvec , template< class > class SPHKernel>
void shammodels::gsph::Model< Tvec, SPHKernel >::set_cfl_force ( Tscal  cfl_force)
inline

Definition at line 101 of file Model.hpp.

◆ set_field_in_box()

template<class Tvec , template< class > class SPHKernel>
template<class T >
void shammodels::gsph::Model< Tvec, SPHKernel >::set_field_in_box ( std::string  field_name,
val,
std::pair< Tvec, Tvec >  box,
u32  ivar = 0 
)
inline

Set field value for particles within a box region.

Sets the specified field to a constant value for all particles whose positions fall within the given axis-aligned box. Useful for setting up discontinuous initial conditions (e.g., Sod shock tube).

Template Parameters
TField type (e.g., Tscal for scalars, Tvec for vectors)
Parameters
field_nameName of the field to modify (e.g., "uint", "vxyz")
valValue to set for particles in the region
boxBounding box as (min_corner, max_corner)
ivarVariable index for multi-variable fields (default: 0)

Example:

// Sod shock tube: set left state internal energy
model.set_field_in_box("uint", u_left, {box_min, interface_pos});
// Set right state
model.set_field_in_box("uint", u_right, {interface_pos, box_max});

Definition at line 215 of file Model.hpp.

+ Here is the call graph for this function:

◆ set_field_in_sphere()

template<class Tvec , template< class > class SPHKernel>
template<class T >
void shammodels::gsph::Model< Tvec, SPHKernel >::set_field_in_sphere ( std::string  field_name,
val,
Tvec  center,
Tscal  radius 
)
inline

Set field value for particles within a spherical region.

Sets the specified field to a constant value for all particles whose positions fall within the given sphere. Useful for setting up point-source initial conditions (e.g., Sedov blast).

Template Parameters
TField type (must be single-variable, e.g., Tscal)
Parameters
field_nameName of the field to modify (e.g., "uint")
valValue to set for particles in the region
centerCenter of the sphere
radiusRadius of the sphere

Example:

// Sedov blast: inject energy in central sphere
Tscal blast_energy_per_particle = E_blast / n_particles_in_sphere;
model.set_field_in_sphere("uint", blast_energy_per_particle, origin, r_blast);

Definition at line 274 of file Model.hpp.

+ Here is the call graph for this function:

◆ set_particle_mass()

template<class Tvec , template< class > class SPHKernel>
void shammodels::gsph::Model< Tvec, SPHKernel >::set_particle_mass ( Tscal  gpart_mass)
inline

Definition at line 105 of file Model.hpp.

◆ set_solver_config()

template<class Tvec , template< class > class SPHKernel>
void shammodels::gsph::Model< Tvec, SPHKernel >::set_solver_config ( SolverConfig  cfg)
inline

Definition at line 335 of file Model.hpp.

◆ solver_logs_last_obj_count()

template<class Tvec , template< class > class SPHKernel>
u64 shammodels::gsph::Model< Tvec, SPHKernel >::solver_logs_last_obj_count ( )
inline

Definition at line 345 of file Model.hpp.

◆ solver_logs_last_rate()

template<class Tvec , template< class > class SPHKernel>
f64 shammodels::gsph::Model< Tvec, SPHKernel >::solver_logs_last_rate ( )
inline

Definition at line 344 of file Model.hpp.

◆ solver_logs_last_system_metrics()

template<class Tvec , template< class > class SPHKernel>
shamsys::SystemMetrics shammodels::gsph::Model< Tvec, SPHKernel >::solver_logs_last_system_metrics ( )
inline

Definition at line 346 of file Model.hpp.

◆ timestep()

template<class Tvec , template< class > class SPHKernel>
TimestepLog shammodels::gsph::Model< Tvec, SPHKernel >::timestep ( )
inline

Definition at line 395 of file Model.hpp.

◆ total_mass_to_part_mass()

template<class Tvec , template< class > class SPHKernel>
f64 shammodels::gsph::Model< Tvec, SPHKernel >::total_mass_to_part_mass ( f64  totmass)

Definition at line 76 of file Model.cpp.

Member Data Documentation

◆ ctx

template<class Tvec , template< class > class SPHKernel>
ShamrockCtx& shammodels::gsph::Model< Tvec, SPHKernel >::ctx

Definition at line 72 of file Model.hpp.

◆ dim

template<class Tvec , template< class > class SPHKernel>
constexpr u32 shammodels::gsph::Model< Tvec, SPHKernel >::dim = shambase::VectorProperties<Tvec>::dimension
staticconstexpr

Definition at line 66 of file Model.hpp.

◆ solver

template<class Tvec , template< class > class SPHKernel>
Solver shammodels::gsph::Model< Tvec, SPHKernel >::solver

Definition at line 73 of file Model.hpp.


The documentation for this class was generated from the following files: