chrono::fsi::sph::ChFsiFluidSystemSPH Class Reference

Description

Physical system for an FSI-aware SPH fluid solver.

#include <ChFsiFluidSystemSPH.h>

Inheritance diagram for chrono::fsi::sph::ChFsiFluidSystemSPH:
Collaboration diagram for chrono::fsi::sph::ChFsiFluidSystemSPH:

Classes

struct  ElasticMaterialProperties
 Structure with elastic material properties. More...
 
struct  FluidProperties
 Structure with fluid properties. More...
 
struct  LinSolverParameters
 Structure with linear solver parameters (used only for implicit SPH). More...
 
struct  SPHParameters
 Structure with SPH method parameters. More...
 
struct  SplashsurfParameters
 Structure with surface reconstruction parameters. More...
 

Public Member Functions

void ReadParametersFromFile (const std::string &json_file)
 Read Chrono::FSI parameters from the specified JSON file.
 
void EnableCudaErrorCheck (bool val)
 Enable/disable CUDA error checks (default: enabled).
 
void SetInitialSpacing (double spacing)
 Set initial spacing.
 
void SetKernelMultiplier (double multiplier)
 Set multiplier for interaction length. More...
 
void SetShiftingMethod (ShiftingMethod shifting_method)
 Set the shifting method.
 
void SetContainerDim (const ChVector3d &box_dim)
 Set the fluid container dimension.
 
void SetComputationalDomain (const ChAABB &computational_AABB, BoundaryConditions bc_type)
 Set computational domain and boundary conditions on its sides. More...
 
void SetComputationalDomain (const ChAABB &computational_AABB)
 Set computational domain. More...
 
void SetActiveDomain (const ChVector3d &box_dim)
 Set dimensions of the active domain AABB. More...
 
void SetActiveDomainDelay (double duration)
 Disable use of the active domain for the given duration at the beginning of the simulation (default: 0). More...
 
void SetNumBCELayers (int num_layers)
 Set number of BCE marker layers (default: 3).
 
void SetDensity (double rho0)
 Set (initial) density.
 
void SetShiftingPPSTParameters (double push, double pull)
 Set the PPST Shifting parameters push: coefficient for the pushing term in the PPST shifting method (upon penetration with fictitious sphere) pull: coefficient for the pulling term in the PPST shifting method.
 
void SetShiftingXSPHParameters (double eps)
 Set the XSPH Shifting parameters eps: coefficient for the XSPH shifting method.
 
void SetShiftingDiffusionParameters (double A, double AFSM, double AFST)
 Set the diffusion based shifting parameters A: coefficient for the diffusion based shifting method AFSM: coefficient for the AFSM in the diffusion based shifting method AFST: coefficient for the AFST in the diffusion based shifting method.
 
void SetInitPressure (const double fzDim)
 Set prescribed initial pressure for gravity field.
 
virtual void SetGravitationalAcceleration (const ChVector3d &gravity) override
 Set gravity for the FSI syatem.
 
void SetBodyForce (const ChVector3d &force)
 Set a constant force applied to the fluid. More...
 
void SetConsistentDerivativeDiscretization (bool consistent_gradient, bool consistent_Laplacian)
 Set SPH discretization type, consistent or inconsistent.
 
void SetCohesionForce (double Fc)
 Set cohesion force of the granular material.
 
void SetSPHLinearSolver (SolverType lin_solver)
 Set the linear system solver for implicit methods.
 
void SetIntegrationScheme (IntegrationScheme scheme)
 Set the integration scheme (default: RK2). More...
 
void SetNumProximitySearchSteps (int steps)
 Set the number of steps between successive updates to neighbor lists (default: 4).
 
void SetCfdSPH (const FluidProperties &fluid_props)
 Enable solution of a CFD problem.
 
void SetElasticSPH (const ElasticMaterialProperties &mat_props)
 Enable solution of elastic SPH (for continuum representation of granular dynamics). More...
 
void CheckSPHParameters ()
 Checks the applicability of user set parameters for SPH and throws an exception if necessary.
 
void SetSPHParameters (const SPHParameters &sph_params)
 Set SPH method parameters.
 
void SetLinSolverParameters (const LinSolverParameters &linsolv_params)
 Set linear solver parameters (used only for implicit SPH).
 
void SetOutputLevel (OutputLevel output_level)
 Set simulation data output level (default: STATE_PRESSURE). More...
 
void SetBoundaryType (BoundaryMethod boundary_method)
 Set boundary treatment type (default: Adami).
 
void SetViscosityType (ViscosityMethod viscosity_method)
 Set viscosity treatment type (default: artificial unilateral).
 
void SetArtificialViscosityCoefficient (double coefficient)
 Set artificial viscosity coefficient (default: 0.02).
 
void SetKernelType (KernelType kernel_type)
 Set kernel type.
 
double GetKernelLength () const
 Return the SPH kernel length of kernel function.
 
double GetInitialSpacing () const
 Return the initial spacing of the SPH particles.
 
int GetNumBCELayers () const
 Return the number of BCE layers.
 
ChVector3d GetContainerDim () const
 Get the fluid container dimensions.
 
ChAABB GetComputationalDomain () const
 Get the computational domain.
 
double GetDensity () const
 Return density.
 
double GetViscosity () const
 Return viscosity.
 
double GetParticleMass () const
 Return SPH particle mass.
 
double GetBasePressure () const
 Return base pressure.
 
ChVector3d GetGravitationalAcceleration () const
 Return gravitational acceleration.
 
double GetSoundSpeed () const
 Return the speed of sound in the fluid phase.
 
ChVector3d GetBodyForce () const
 Return the constant force applied to the fluid (if any). More...
 
int GetNumProximitySearchSteps () const
 Get the number of steps between successive updates to neighbor lists.
 
const ChFsiParamsSPHGetParams () const
 Return the current system parameters (debugging only).
 
size_t GetNumFluidMarkers () const
 Get the current number of fluid SPH particles.
 
size_t GetNumBoundaryMarkers () const
 Get the current number of boundary BCE markers.
 
size_t GetNumRigidBodyMarkers () const
 Get the current number of rigid body BCE markers.
 
size_t GetNumFlexBodyMarkers () const
 Get the current number of flexible body BCE markers.
 
std::vector< ChVector3dGetParticlePositions () const
 Return the SPH particle positions.
 
std::vector< ChVector3dGetParticleVelocities () const
 Return the SPH particle velocities.
 
std::vector< ChVector3dGetParticleAccelerations () const
 Return the accelerations of SPH particles.
 
std::vector< ChVector3dGetParticleForces () const
 Return the forces acting on SPH particles.
 
std::vector< ChVector3dGetParticleFluidProperties () const
 Return the SPH particle fluid properties. More...
 
BoundaryMethod GetBoundaryType () const
 Return the boundary treatment type.
 
ViscosityMethod GetViscosityType () const
 Return the viscosity treatment type.
 
KernelType GetKernelType () const
 Return the kernel type.
 
void WriteParticleFile (const std::string &filename) const
 Write FSI system particle output.
 
void SaveParticleData (const std::string &dir) const
 Save current SPH particle and BCE marker data to files. More...
 
void SaveSolidData (const std::string &dir, double time) const
 Save current FSI solid data to files. More...
 
void AddSPHParticle (const ChVector3d &pos, double rho, double pres, double mu, const ChVector3d &vel=ChVector3d(0), const ChVector3d &tauXxYyZz=ChVector3d(0), const ChVector3d &tauXyXzYz=ChVector3d(0))
 Add an SPH particle with given properties to the FSI system.
 
void AddSPHParticle (const ChVector3d &pos, const ChVector3d &vel=ChVector3d(0), const ChVector3d &tauXxYyZz=ChVector3d(0), const ChVector3d &tauXyXzYz=ChVector3d(0))
 Add an SPH particle with current properties to the SPH system.
 
void AddBoxSPH (const ChVector3d &boxCenter, const ChVector3d &boxHalfDim)
 Create SPH particles in the specified box volume. More...
 
void AddBCEBoundary (const std::vector< ChVector3d > &points, const ChFramed &frame)
 Add boundary BCE markers at the specified points. More...
 
std::vector< int > FindParticlesInBox (const ChFrame<> &frame, const ChVector3d &size)
 Utility function for finding indices of SPH particles inside a given OBB. More...
 
std::vector< Real3GetPositions () const
 Extract positions of all markers (SPH and BCE).
 
std::vector< Real3GetVelocities () const
 Extract velocities of all markers (SPH and BCE).
 
std::vector< Real3GetAccelerations () const
 Extract accelerations of all markers (SPH and BCE).
 
std::vector< Real3GetForces () const
 Extract forces applied to all markers (SPH and BCE).
 
std::vector< Real3GetProperties () const
 Extract fluid properties of all markers (SPH and BCE). More...
 
std::vector< Real3GetPositions (const std::vector< int > &indices) const
 Extract positions of all markers (SPH and BCE) with indices in the provided array.
 
std::vector< Real3GetVelocities (const std::vector< int > &indices) const
 Extract velocities of all markers (SPH and BCE) with indices in the provided array.
 
std::vector< Real3GetAccelerations (const std::vector< int > &indices) const
 Extract accelerations of all markers (SPH and BCE) with indices in the provided array.
 
std::vector< Real3GetForces (const std::vector< int > &indices) const
 Extract forces applied to allmarkers (SPH and BCE) with indices in the provided array.
 
std::vector< ChVector3dCreatePointsPlate (const ChVector2d &size) const
 Create marker points on a rectangular plate of specified X-Y dimensions, assumed centered at the origin. More...
 
std::vector< ChVector3dCreatePointsBoxContainer (const ChVector3d &size, const ChVector3i &faces) const
 Create marker points for a box container of specified dimensions. More...
 
std::vector< ChVector3dCreatePointsBoxInterior (const ChVector3d &size) const
 Create interior marker points for a box of specified dimensions, assumed centered at the origin. More...
 
std::vector< ChVector3dCreatePointsBoxExterior (const ChVector3d &size) const
 Create exterior marker points for a box of specified dimensions, assumed centered at the origin. More...
 
std::vector< ChVector3dCreatePointsSphereInterior (double radius, bool polar) const
 Create interior marker points for a sphere of specified radius, assumed centered at the origin. More...
 
std::vector< ChVector3dCreatePointsSphereExterior (double radius, bool polar) const
 Create exterior marker pointss for a sphere of specified radius, assumed centered at the origin. More...
 
std::vector< ChVector3dCreatePointsCylinderInterior (double rad, double height, bool polar) const
 Create interior marker points for a cylinder of specified radius and height. More...
 
std::vector< ChVector3dCreatePointsCylinderExterior (double rad, double height, bool polar) const
 Create exterior marker points for a cylinder of specified radius and height. More...
 
std::vector< ChVector3dCreatePointsConeInterior (double rad, double height, bool polar) const
 Create interior marker points for a cone of specified radius and height. More...
 
std::vector< ChVector3dCreatePointsConeExterior (double rad, double height, bool polar) const
 Create exterior marker points for a cone of specified radius and height. More...
 
std::vector< ChVector3dCreatePointsCylinderAnnulus (double rad_inner, double rad_outer, double height, bool polar) const
 Create marker points filling a cylindrical annulus of specified radii and height. More...
 
std::vector< ChVector3dCreatePointsMesh (ChTriangleMeshConnected &mesh) const
 Create marker points filling a closed mesh. More...
 
PhysicsProblem GetPhysicsProblem () const
 
std::string GetPhysicsProblemString () const
 
std::string GetSphIntegrationSchemeString () const
 
- Public Member Functions inherited from chrono::fsi::ChFsiFluidSystem
void SetVerbose (bool verbose)
 Enable/disable verbose terminal output (default: true).
 
void SetStepSize (double step)
 Set integration step size.
 
virtual void Initialize ()
 Initialize the fluid system with no FSI support.
 
void DoStepDynamics (double step)
 Function to integrate the FSI fluid system in time.
 
double GetSimTime () const
 Get current simulation time.
 
double GetStepSize () const
 Get the integration step size.
 
double GetRtf () const
 Get current estimated RTF (real time factor).
 
double GetTimerStep () const
 Return the time in seconds for fluid dynamics over the last step.
 

Friends

class ChFsiSystemSPH
 
class ChFsiInterfaceSPH
 
class ChFsiProblemSPH
 
class ChFsiSplashsurfSPH
 

Additional Inherited Members

- Protected Member Functions inherited from chrono::fsi::ChFsiFluidSystem
void UseNodeDirections (bool val)
 
- Protected Attributes inherited from chrono::fsi::ChFsiFluidSystem
bool m_verbose
 enable/disable m_verbose terminal output
 
std::string m_outdir
 output directory
 
bool m_use_node_directions
 use FEA node directions
 
bool m_is_initialized
 set to true once the Initialize function is called
 
double m_step
 time step for fluid dynamics
 
unsigned int m_frame
 current simulation frame
 

Member Function Documentation

◆ AddBCEBoundary()

void chrono::fsi::sph::ChFsiFluidSystemSPH::AddBCEBoundary ( const std::vector< ChVector3d > &  points,
const ChFramed frame 
)

Add boundary BCE markers at the specified points.

The points are assumed to be provided relative to the specified frame. These BCE markers are not associated with a particular FSI body and, as such, cannot be used to extract fluid forces and moments. If fluid reaction forces are needed, create an FSI body with the desirted geometry or list of BCE points and add it through the contianing FSI system.

◆ AddBoxSPH()

void chrono::fsi::sph::ChFsiFluidSystemSPH::AddBoxSPH ( const ChVector3d boxCenter,
const ChVector3d boxHalfDim 
)

Create SPH particles in the specified box volume.

The SPH particles are created on a uniform grid with resolution equal to the FSI initial separation.

◆ CreatePointsBoxContainer()

std::vector< ChVector3d > chrono::fsi::sph::ChFsiFluidSystemSPH::CreatePointsBoxContainer ( const ChVector3d size,
const ChVector3i faces 
) const

Create marker points for a box container of specified dimensions.

The box volume is assumed to be centered at the origin. The 'faces' input vector specifies which faces of the container are to be created: for each direction, a value of -1 indicates the face in the negative direction, a value of +1 indicates the face in the positive direction, and a value of 2 indicates both faces. Setting a value of 0 does not create container faces in that direction. Markers are created in a number of layers corresponding to system parameters.

◆ CreatePointsBoxExterior()

std::vector< ChVector3d > chrono::fsi::sph::ChFsiFluidSystemSPH::CreatePointsBoxExterior ( const ChVector3d size) const

Create exterior marker points for a box of specified dimensions, assumed centered at the origin.

Markers are created outside the box, in a number of layers corresponding to system parameters.

◆ CreatePointsBoxInterior()

std::vector< ChVector3d > chrono::fsi::sph::ChFsiFluidSystemSPH::CreatePointsBoxInterior ( const ChVector3d size) const

Create interior marker points for a box of specified dimensions, assumed centered at the origin.

Markers are created inside the box, in a number of layers corresponding to system parameters.

◆ CreatePointsConeExterior()

std::vector< ChVector3d > chrono::fsi::sph::ChFsiFluidSystemSPH::CreatePointsConeExterior ( double  rad,
double  height,
bool  polar 
) const

Create exterior marker points for a cone of specified radius and height.

The cone is assumed centered at the origin and aligned with the Z axis. Markers are created outside the cone, in a number of layers corresponding to system parameters. Markers are created using cylinderical coordinates (polar=true), or else on a uniform Cartesian grid.

◆ CreatePointsConeInterior()

std::vector< ChVector3d > chrono::fsi::sph::ChFsiFluidSystemSPH::CreatePointsConeInterior ( double  rad,
double  height,
bool  polar 
) const

Create interior marker points for a cone of specified radius and height.

The cone is assumed centered at the origin and aligned with the Z axis. Markers are created inside the cone, in a number of layers corresponding to system parameters. Markers are created using cylinderical coordinates (polar=true), or else on a uniform Cartesian grid.

◆ CreatePointsCylinderAnnulus()

std::vector< ChVector3d > chrono::fsi::sph::ChFsiFluidSystemSPH::CreatePointsCylinderAnnulus ( double  rad_inner,
double  rad_outer,
double  height,
bool  polar 
) const

Create marker points filling a cylindrical annulus of specified radii and height.

The cylinder annulus is assumed centered at the origin and aligned with the Z axis. Markers are created using cylindrical coordinates (polar=true), or else on a uniform Cartesian grid.

◆ CreatePointsCylinderExterior()

std::vector< ChVector3d > chrono::fsi::sph::ChFsiFluidSystemSPH::CreatePointsCylinderExterior ( double  rad,
double  height,
bool  polar 
) const

Create exterior marker points for a cylinder of specified radius and height.

The cylinder is assumed centered at the origin and aligned with the Z axis. Markers are created outside the cylinder, in a number of layers corresponding to system parameters. Markers are created using cylindrical coordinates (polar=true), or else on a uniform Cartesian grid.

◆ CreatePointsCylinderInterior()

std::vector< ChVector3d > chrono::fsi::sph::ChFsiFluidSystemSPH::CreatePointsCylinderInterior ( double  rad,
double  height,
bool  polar 
) const

Create interior marker points for a cylinder of specified radius and height.

The cylinder is assumed centered at the origin and aligned with the Z axis. Markers are created inside the cylinder, in a number of layers corresponding to system parameters. Markers are created using cylindrical coordinates (polar=true), or else on a uniform Cartesian grid.

◆ CreatePointsMesh()

std::vector< ChVector3d > chrono::fsi::sph::ChFsiFluidSystemSPH::CreatePointsMesh ( ChTriangleMeshConnected mesh) const

Create marker points filling a closed mesh.

Markers are created on a Cartesian grid with a separation corresponding to system parameters.

Calculate uu parameter and test bound

ray intersection

◆ CreatePointsPlate()

std::vector< ChVector3d > chrono::fsi::sph::ChFsiFluidSystemSPH::CreatePointsPlate ( const ChVector2d size) const

Create marker points on a rectangular plate of specified X-Y dimensions, assumed centered at the origin.

Markers are created in a number of layers (in the negative Z direction) corresponding to system parameters.

◆ CreatePointsSphereExterior()

std::vector< ChVector3d > chrono::fsi::sph::ChFsiFluidSystemSPH::CreatePointsSphereExterior ( double  radius,
bool  polar 
) const

Create exterior marker pointss for a sphere of specified radius, assumed centered at the origin.

Markers are created outside the sphere, in a number of layers corresponding to system parameters. Markers are created using spherical coordinates (polar=true), or else on a uniform Cartesian grid.

◆ CreatePointsSphereInterior()

std::vector< ChVector3d > chrono::fsi::sph::ChFsiFluidSystemSPH::CreatePointsSphereInterior ( double  radius,
bool  polar 
) const

Create interior marker points for a sphere of specified radius, assumed centered at the origin.

Markers are created inside the sphere, in a number of layers corresponding to system parameters. Markers are created using spherical coordinates (polar=true), or else on a uniform Cartesian grid.

◆ FindParticlesInBox()

std::vector< int > chrono::fsi::sph::ChFsiFluidSystemSPH::FindParticlesInBox ( const ChFrame<> &  frame,
const ChVector3d size 
)

Utility function for finding indices of SPH particles inside a given OBB.

The object-oriented box, of specified size, is assumed centered at the origin of the provided frame and aligned with the axes of that frame.

◆ GetBodyForce()

ChVector3d chrono::fsi::sph::ChFsiFluidSystemSPH::GetBodyForce ( ) const

Return the constant force applied to the fluid (if any).

TODO: RENAME

◆ GetParticleFluidProperties()

std::vector< ChVector3d > chrono::fsi::sph::ChFsiFluidSystemSPH::GetParticleFluidProperties ( ) const

Return the SPH particle fluid properties.

For each SPH particle, the 3-dimensional vector contains density, pressure, and viscosity.

◆ GetProperties()

std::vector< Real3 > chrono::fsi::sph::ChFsiFluidSystemSPH::GetProperties ( ) const

Extract fluid properties of all markers (SPH and BCE).

For each SPH particle, the 3-dimensional vector contains density, pressure, and viscosity.

◆ SaveParticleData()

void chrono::fsi::sph::ChFsiFluidSystemSPH::SaveParticleData ( const std::string &  dir) const

Save current SPH particle and BCE marker data to files.

This function creates three CSV files for SPH particles, boundary BCE markers, and solid BCE markers data.

◆ SaveSolidData()

void chrono::fsi::sph::ChFsiFluidSystemSPH::SaveSolidData ( const std::string &  dir,
double  time 
) const

Save current FSI solid data to files.

This function creates CSV files for force and torque on rigid bodies and flexible nodes.

◆ SetActiveDomain()

void chrono::fsi::sph::ChFsiFluidSystemSPH::SetActiveDomain ( const ChVector3d box_dim)

Set dimensions of the active domain AABB.

This value activates only those SPH particles that are within an AABB of the specified size from an object interacting with the "fluid" phase. Note that this setting should not be used for CFD simulations, but rather only when solving problems using the CRM (continuum representation of granular dynamics) for terramechanics simulations.

◆ SetActiveDomainDelay()

void chrono::fsi::sph::ChFsiFluidSystemSPH::SetActiveDomainDelay ( double  duration)

Disable use of the active domain for the given duration at the beginning of the simulation (default: 0).

This parameter is used for settling operations where all particles must be active through the settling process.

◆ SetBodyForce()

void chrono::fsi::sph::ChFsiFluidSystemSPH::SetBodyForce ( const ChVector3d force)

Set a constant force applied to the fluid.

Solid bodies are not explicitly affected by this force, but they are affected indirectly through the fluid.

◆ SetComputationalDomain() [1/2]

void chrono::fsi::sph::ChFsiFluidSystemSPH::SetComputationalDomain ( const ChAABB computational_AABB)

Set computational domain.

Note that this version leaves the setting for BC type unchanged.

◆ SetComputationalDomain() [2/2]

void chrono::fsi::sph::ChFsiFluidSystemSPH::SetComputationalDomain ( const ChAABB computational_AABB,
BoundaryConditions  bc_type 
)

Set computational domain and boundary conditions on its sides.

bc_type indicates the types of BCs imposed in the three directions of the computational domain. By default, no special boundary conditions are imposed in any direction (BCType::NONE).

◆ SetElasticSPH()

void chrono::fsi::sph::ChFsiFluidSystemSPH::SetElasticSPH ( const ElasticMaterialProperties mat_props)

Enable solution of elastic SPH (for continuum representation of granular dynamics).

By default, a ChSystemFSI solves an SPH fluid dynamics problem.

◆ SetIntegrationScheme()

void chrono::fsi::sph::ChFsiFluidSystemSPH::SetIntegrationScheme ( IntegrationScheme  scheme)

Set the integration scheme (default: RK2).

All explicit integration schemes use a so-called Weakly-Compressible SPH formulation, based on an equation of state that relates pressure to density.

◆ SetKernelMultiplier()

void chrono::fsi::sph::ChFsiFluidSystemSPH::SetKernelMultiplier ( double  multiplier)

Set multiplier for interaction length.

h = multiplier * initial_spacing.

◆ SetOutputLevel()

void chrono::fsi::sph::ChFsiFluidSystemSPH::SetOutputLevel ( OutputLevel  output_level)

Set simulation data output level (default: STATE_PRESSURE).

Options:

  • STATE marker state, velocity, and acceleration
  • STATE_PRESSURE STATE plus density and pressure
  • CFD_FULL STATE_PRESSURE plus various CFD parameters
  • CRM_FULL STATE_PRESSURE plus normal and shear stress

The documentation for this class was generated from the following files:
  • /builds/uwsbel/chrono/src/chrono_fsi/sph/ChFsiFluidSystemSPH.h
  • /builds/uwsbel/chrono/src/chrono_fsi/sph/ChFsiFluidSystemSPH.cpp