Description
Physical system for and SPH FSI-aware fluid solver.
#include <ChFluidSystemSPH.h>
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... | |
Public Member Functions | |
void | ReadParametersFromFile (const std::string &json_file) |
Read Chrono::FSI parameters from the specified JSON file. | |
void | SetInitialSpacing (double spacing) |
Set initial spacing. | |
void | SetKernelMultiplier (double multiplier) |
Set multiplier for interaction length. More... | |
void | SetContainerDim (const ChVector3d &boxDim) |
Set the fluid container dimension. | |
void | SetBoundaries (const ChVector3d &cMin, const ChVector3d &cMax) |
Set periodic boundary condition for fluid. | |
void | SetActiveDomain (const ChVector3d &boxHalfDim) |
Set half-dimensions of the active domain. 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 | 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 | SetSPHMethod (SPHMethod SPH_method) |
Set the SPH method and, optionally, the linear solver type. | |
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 (BoundaryType boundary_type) |
Set boundary treatment type (default: Adami). | |
void | SetViscosityType (ViscosityType viscosity_type) |
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. | |
virtual void | Initialize (unsigned int num_fsi_bodies, unsigned int num_fsi_nodes1D, unsigned int num_fsi_elements1D, unsigned int num_fsi_nodes2D, unsigned int num_fsi_elements2D, const std::vector< FsiBodyState > &body_states, const std::vector< FsiMeshState > &mesh1D_states, const std::vector< FsiMeshState > &mesh2D_states) override |
Initialize the SPH fluid system with FSI support. | |
virtual void | Initialize () override |
Initialize the SPH fluid system with no FSI support. | |
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 |
Set the fluid container dimension. | |
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 sph::SimParams & | GetParams () 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< ChVector3d > | GetParticlePositions () const |
Return the SPH particle positions. | |
std::vector< ChVector3d > | GetParticleVelocities () const |
Return the SPH particle velocities. | |
std::vector< ChVector3d > | GetParticleAccelerations () const |
Return the accelerations of SPH particles. | |
std::vector< ChVector3d > | GetParticleForces () const |
Return the forces acting on SPH particles. | |
std::vector< ChVector3d > | GetParticleFluidProperties () const |
Return the SPH particle fluid properties. More... | |
BoundaryType | GetBoundaryType () const |
Return the boundary treatment type. | |
ViscosityType | GetViscosityType () const |
Return the viscosity treatment type. | |
KernelType | GetKernelType () const |
Return the kernel type. | |
void | WriteParticleFile (const std::string &filename, OutputMode mode) 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 | AddPlateBCE (std::shared_ptr< ChBody > body, const ChFrame<> &frame, const ChVector2d &size) |
Add BCE markers for a rectangular plate of specified X-Y dimensions and associate them with the given body. More... | |
void | AddBoxContainerBCE (std::shared_ptr< ChBody > body, const ChFrame<> &frame, const ChVector3d &size, const ChVector3i faces) |
Add BCE markers for a box container of specified dimensions and associate them with the given body. More... | |
size_t | AddBoxBCE (std::shared_ptr< ChBody > body, const ChFrame<> &frame, const ChVector3d &size, bool solid) |
Add BCE markers for a box of specified dimensions and associate them with the given body. More... | |
size_t | AddSphereBCE (std::shared_ptr< ChBody > body, const ChFrame<> &frame, double radius, bool solid, bool polar=true) |
Add BCE markers for a sphere of specified radius and associate them with the given body. More... | |
size_t | AddCylinderBCE (std::shared_ptr< ChBody > body, const ChFrame<> &frame, double radius, double height, bool solid, bool polar=true) |
Add BCE markers for a cylinder of specified radius and height and associate them with the given body. More... | |
size_t | AddConeBCE (std::shared_ptr< ChBody > body, const ChFrame<> &frame, double radius, double height, bool solid, bool polar=true) |
Add BCE markers for a cone of specified radius and height and associate them with the given body. More... | |
size_t | AddCylinderAnnulusBCE (std::shared_ptr< ChBody > body, const ChFrame<> &frame, double radius_inner, double radius_outer, double height, bool polar=true) |
Add BCE markers for a cylindrical annulus of specified radii and height and associate them with the given body. More... | |
size_t | AddPointsBCE (std::shared_ptr< ChBody > body, const std::vector< ChVector3d > &points, const ChFrame<> &rel_frame, bool solid) |
Add BCE markers from a set of points and associate them with the given body. More... | |
void | SetBcePattern1D (BcePatternMesh1D pattern, bool remove_center=false) |
Set the BCE marker pattern for 1D flexible solids for subsequent calls to AddFsiMesh. More... | |
void | SetBcePattern2D (BcePatternMesh2D pattern, bool remove_center=false) |
Set the BCE marker pattern for 2D flexible solids for subsequent calls to AddFsiMesh. More... | |
thrust::device_vector< int > | FindParticlesInBox (const ChFrame<> &frame, const ChVector3d &size) |
Utility function for finding indices of SPH particles inside a given OBB. More... | |
thrust::device_vector< sph::Real3 > | GetPositions () |
Extract positions of all markers (SPH and BCE). More... | |
thrust::device_vector< sph::Real3 > | GetVelocities () |
Extract velocities of all markers (SPH and BCE). More... | |
thrust::device_vector< sph::Real3 > | GetAccelerations () |
Extract accelerations of all markers (SPH and BCE). More... | |
thrust::device_vector< sph::Real3 > | GetForces () |
Extract forces applied to all markers (SPH and BCE). More... | |
thrust::device_vector< sph::Real3 > | GetProperties () |
Extract fluid properties of all markers (SPH and BCE). More... | |
thrust::device_vector< sph::Real3 > | GetPositions (const thrust::device_vector< int > &indices) |
Extract positions of all markers (SPH and BCE) with indices in the provided array. More... | |
thrust::device_vector< sph::Real3 > | GetVelocities (const thrust::device_vector< int > &indices) |
Extract velocities of all markers (SPH and BCE) with indices in the provided array. More... | |
thrust::device_vector< sph::Real3 > | GetAccelerations (const thrust::device_vector< int > &indices) |
Extract accelerations of all markers (SPH and BCE) with indices in the provided array. More... | |
thrust::device_vector< sph::Real3 > | GetForces (const thrust::device_vector< int > &indices) |
Extract forces applied to allmarkers (SPH and BCE) with indices in the provided array. More... | |
void | CreateBCE_Plate (const ChVector2d &size, std::vector< ChVector3d > &bce) |
Create BCE markers on a rectangular plate of specified X-Y dimensions, assumed centered at the origin. More... | |
void | CreateBCE_BoxInterior (const ChVector3d &size, std::vector< ChVector3d > &bce) |
Create BCE interior markers for a box of specified dimensions, assumed centered at the origin. More... | |
void | CreateBCE_BoxExterior (const ChVector3d &size, std::vector< ChVector3d > &bce) |
Create exterior BCE markers for a box of specified dimensions, assumed centered at the origin. More... | |
void | CreateBCE_SphereInterior (double radius, bool polar, std::vector< ChVector3d > &bce) |
Create interior BCE markers for a sphere of specified radius, assumed centered at the origin. More... | |
void | CreateBCE_SphereExterior (double radius, bool polar, std::vector< ChVector3d > &bce) |
Create exterior BCE markers for a sphere of specified radius, assumed centered at the origin. More... | |
void | CreateBCE_CylinderInterior (double rad, double height, bool polar, std::vector< ChVector3d > &bce) |
Create interior BCE markers for a cylinder of specified radius and height. More... | |
void | CreateBCE_CylinderExterior (double rad, double height, bool polar, std::vector< ChVector3d > &bce) |
Create exterior BCE markers for a cylinder of specified radius and height. More... | |
void | CreateBCE_ConeInterior (double rad, double height, bool polar, std::vector< ChVector3d > &bce) |
Create interior BCE particles for a cone of specified radius and height. More... | |
void | CreateBCE_ConeExterior (double rad, double height, bool polar, std::vector< ChVector3d > &bce) |
Create exterior BCE particles for a cone of specified radius and height. More... | |
PhysicsProblem | GetPhysicsProblem () const |
std::string | GetPhysicsProblemString () const |
std::string | GetSphMethodTypeString () const |
Public Member Functions inherited from chrono::fsi::ChFluidSystem | |
void | SetVerbose (bool verbose) |
Enable/disable verbose terminal output (default: true). | |
void | SetStepSize (double step) |
Set integration step size. | |
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. | |
Static Public Member Functions | |
static void | CreatePoints_CylinderAnnulus (double rad_inner, double rad_outer, double height, bool polar, double delta, std::vector< ChVector3d > &points) |
Create points at suggested separation 'delta' filling a cylindrical annulus of specified radii and height. More... | |
static void | CreatePoints_Mesh (ChTriangleMeshConnected &mesh, double delta, std::vector< ChVector3d > &points) |
Create points at suggested separation 'delta' filling a closed mesh. More... | |
Friends | |
class | ChFsiSystemSPH |
class | ChFsiInterfaceSPH |
class | ChFsiVisualizationGL |
class | ChFsiVisualizationVSG |
Additional Inherited Members | |
Protected Attributes inherited from chrono::fsi::ChFluidSystem | |
bool | m_verbose |
enable/disable m_verbose terminal output | |
std::string | m_outdir |
output directory | |
bool | m_is_initialized |
set to true once the Initialize function is called | |
double | m_step |
time step for fluid dynamics | |
double | m_time |
current fluid dynamics simulation time | |
Member Function Documentation
◆ AddBoxBCE()
size_t chrono::fsi::ChFluidSystemSPH::AddBoxBCE | ( | std::shared_ptr< ChBody > | body, |
const ChFrame<> & | frame, | ||
const ChVector3d & | size, | ||
bool | solid | ||
) |
Add BCE markers for a box of specified dimensions and associate them with the given body.
The box is assumed to be centered at the origin of the provided frame and aligned with its axes. BCE markers are created inside the box if solid=true, and outside the box otherwise. BCE markers are created in a number of layers corresponding to system parameters.
◆ AddBoxContainerBCE()
void chrono::fsi::ChFluidSystemSPH::AddBoxContainerBCE | ( | std::shared_ptr< ChBody > | body, |
const ChFrame<> & | frame, | ||
const ChVector3d & | size, | ||
const ChVector3i | faces | ||
) |
Add BCE markers for a box container of specified dimensions and associate them with the given body.
The center of the box volume is at the origin of the given frame and the the container is aligned with the frame axes. Such a container is assumed to be used as a boundary. 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. BCE markers are created in a number of layers corresponding to system parameters. Such a container is assumed to be used as a fixed boundary and the associated body is not tracked in FSI.
◆ AddBoxSPH()
void chrono::fsi::ChFluidSystemSPH::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.
◆ AddConeBCE()
size_t chrono::fsi::ChFluidSystemSPH::AddConeBCE | ( | std::shared_ptr< ChBody > | body, |
const ChFrame<> & | frame, | ||
double | radius, | ||
double | height, | ||
bool | solid, | ||
bool | polar = true |
||
) |
Add BCE markers for a cone of specified radius and height and associate them with the given body.
The cone is assumed centered at the origin of the provided frame and aligned with its Z axis. BCE markers are created inside the cone if solid=true, and outside the cone otherwise. BCE markers are created in a number of layers corresponding to system parameters. BCE markers are created using cylinderical coordinates (default), or else on a uniform Cartesian grid.
◆ AddCylinderAnnulusBCE()
size_t chrono::fsi::ChFluidSystemSPH::AddCylinderAnnulusBCE | ( | std::shared_ptr< ChBody > | body, |
const ChFrame<> & | frame, | ||
double | radius_inner, | ||
double | radius_outer, | ||
double | height, | ||
bool | polar = true |
||
) |
Add BCE markers for a cylindrical annulus of specified radii and height and associate them with the given body.
The cylindrical annulus is assumed centered at the origin of the provided frame and aligned with its Z axis. BCE markers are created in a number of layers corresponding to system parameters. BCE markers are created using cylinderical coordinates (default), or else on a uniform Cartesian grid. Such a cylindrical annulus is assumed to be used on a solid body.
◆ AddCylinderBCE()
size_t chrono::fsi::ChFluidSystemSPH::AddCylinderBCE | ( | std::shared_ptr< ChBody > | body, |
const ChFrame<> & | frame, | ||
double | radius, | ||
double | height, | ||
bool | solid, | ||
bool | polar = true |
||
) |
Add BCE markers for a cylinder of specified radius and height and associate them with the given body.
The cylinder is assumed centered at the origin of the provided frame and aligned with its Z axis. BCE markers are created inside the cylinder if solid=true, and outside the cylinder otherwise. BCE markers are created in a number of layers corresponding to system parameters. BCE markers are created using cylinderical coordinates (default), or else on a uniform Cartesian grid.
◆ AddPlateBCE()
void chrono::fsi::ChFluidSystemSPH::AddPlateBCE | ( | std::shared_ptr< ChBody > | body, |
const ChFrame<> & | frame, | ||
const ChVector2d & | size | ||
) |
Add BCE markers for a rectangular plate of specified X-Y dimensions and associate them with the given body.
BCE markers are created in a number of layers corresponding to system parameters. X-Y BCE layers are created in the negative Z direction of the plate orientation frame. Such a plate is assumed to be used as boundary.
◆ AddPointsBCE()
size_t chrono::fsi::ChFluidSystemSPH::AddPointsBCE | ( | std::shared_ptr< ChBody > | body, |
const std::vector< ChVector3d > & | points, | ||
const ChFrame<> & | rel_frame, | ||
bool | solid | ||
) |
Add BCE markers from a set of points and associate them with the given body.
The points are assumed to be provided relative to the specified frame. The BCE markers are created in the absolute coordinate frame.
◆ AddSphereBCE()
size_t chrono::fsi::ChFluidSystemSPH::AddSphereBCE | ( | std::shared_ptr< ChBody > | body, |
const ChFrame<> & | frame, | ||
double | radius, | ||
bool | solid, | ||
bool | polar = true |
||
) |
Add BCE markers for a sphere of specified radius and associate them with the given body.
The sphere is assumed to be centered at the origin of the provided frame. BCE markers are created inside the sphere if solid=true, and outside the sphere otherwise. BCE markers are created in a number of layers corresponding to system parameters. BCE markers are created using spherical coordinates (default), or else on a uniform Cartesian grid.
◆ CreateBCE_BoxExterior()
void chrono::fsi::ChFluidSystemSPH::CreateBCE_BoxExterior | ( | const ChVector3d & | size, |
std::vector< ChVector3d > & | bce | ||
) |
Create exterior BCE markers for a box of specified dimensions, assumed centered at the origin.
BCE markers are created outside the box, in a number of layers corresponding to system parameters.
◆ CreateBCE_BoxInterior()
void chrono::fsi::ChFluidSystemSPH::CreateBCE_BoxInterior | ( | const ChVector3d & | size, |
std::vector< ChVector3d > & | bce | ||
) |
Create BCE interior markers for a box of specified dimensions, assumed centered at the origin.
BCE markers are created inside the box, in a number of layers corresponding to system parameters.
◆ CreateBCE_ConeExterior()
void chrono::fsi::ChFluidSystemSPH::CreateBCE_ConeExterior | ( | double | rad, |
double | height, | ||
bool | polar, | ||
std::vector< ChVector3d > & | bce | ||
) |
Create exterior BCE particles for a cone of specified radius and height.
The cone is assumed centered at the origin and aligned with the Z axis. BCE markers are created outside the cone, in a number of layers corresponding to system parameters. BCE markers are created using cylinderical coordinates (polar=true), or else on a uniform Cartesian grid.
◆ CreateBCE_ConeInterior()
void chrono::fsi::ChFluidSystemSPH::CreateBCE_ConeInterior | ( | double | rad, |
double | height, | ||
bool | polar, | ||
std::vector< ChVector3d > & | bce | ||
) |
Create interior BCE particles for a cone of specified radius and height.
The cone is assumed centered at the origin and aligned with the Z axis. BCE markers are created inside the cone, in a number of layers corresponding to system parameters. BCE markers are created using cylinderical coordinates (polar=true), or else on a uniform Cartesian grid.
◆ CreateBCE_CylinderExterior()
void chrono::fsi::ChFluidSystemSPH::CreateBCE_CylinderExterior | ( | double | rad, |
double | height, | ||
bool | polar, | ||
std::vector< ChVector3d > & | bce | ||
) |
Create exterior BCE markers for a cylinder of specified radius and height.
The cylinder is assumed centered at the origin and aligned with the Z axis. BCE markers are created outside the cylinder, in a number of layers corresponding to system parameters. BCE markers are created using cylindrical coordinates (polar=true), or else on a uniform Cartesian grid.
◆ CreateBCE_CylinderInterior()
void chrono::fsi::ChFluidSystemSPH::CreateBCE_CylinderInterior | ( | double | rad, |
double | height, | ||
bool | polar, | ||
std::vector< ChVector3d > & | bce | ||
) |
Create interior BCE markers for a cylinder of specified radius and height.
The cylinder is assumed centered at the origin and aligned with the Z axis. BCE markers are created inside the cylinder, in a number of layers corresponding to system parameters. BCE markers are created using cylindrical coordinates (polar=true), or else on a uniform Cartesian grid.
◆ CreateBCE_Plate()
void chrono::fsi::ChFluidSystemSPH::CreateBCE_Plate | ( | const ChVector2d & | size, |
std::vector< ChVector3d > & | bce | ||
) |
Create BCE markers on a rectangular plate of specified X-Y dimensions, assumed centered at the origin.
BCE markers are created in a number of layers corresponding to system parameters. BCE layers are created in the negative Z direction.
◆ CreateBCE_SphereExterior()
void chrono::fsi::ChFluidSystemSPH::CreateBCE_SphereExterior | ( | double | radius, |
bool | polar, | ||
std::vector< ChVector3d > & | bce | ||
) |
Create exterior BCE markers for a sphere of specified radius, assumed centered at the origin.
BCE markers are created outside the sphere, in a number of layers corresponding to system parameters. BCE markers are created using spherical coordinates (polar=true), or else on a uniform Cartesian grid.
◆ CreateBCE_SphereInterior()
void chrono::fsi::ChFluidSystemSPH::CreateBCE_SphereInterior | ( | double | radius, |
bool | polar, | ||
std::vector< ChVector3d > & | bce | ||
) |
Create interior BCE markers for a sphere of specified radius, assumed centered at the origin.
BCE markers are created inside the sphere, in a number of layers corresponding to system parameters. BCE markers are created using spherical coordinates (polar=true), or else on a uniform Cartesian grid.
◆ CreatePoints_CylinderAnnulus()
|
static |
Create points at suggested separation 'delta' filling a cylindrical annulus of specified radii and height.
The cylinder annulus is assumed centered at the origin and aligned with the Z axis. BCE markers are created using cylindrical coordinates (polar=true), or else on a uniform Cartesian grid.
◆ CreatePoints_Mesh()
|
static |
Create points at suggested separation 'delta' filling a closed mesh.
Calculate uu parameter and test bound
ray intersection
◆ FindParticlesInBox()
thrust::device_vector< int > chrono::fsi::ChFluidSystemSPH::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. The return value is a device thrust vector.
◆ GetAccelerations() [1/2]
thrust::device_vector< Real3 > chrono::fsi::ChFluidSystemSPH::GetAccelerations | ( | ) |
Extract accelerations of all markers (SPH and BCE).
The return value is a device thrust vector.
◆ GetAccelerations() [2/2]
thrust::device_vector< Real3 > chrono::fsi::ChFluidSystemSPH::GetAccelerations | ( | const thrust::device_vector< int > & | indices | ) |
Extract accelerations of all markers (SPH and BCE) with indices in the provided array.
The return value is a device thrust vector.
◆ GetBodyForce()
ChVector3d chrono::fsi::ChFluidSystemSPH::GetBodyForce | ( | ) | const |
Return the constant force applied to the fluid (if any).
TODO: RENAME
◆ GetForces() [1/2]
thrust::device_vector< Real3 > chrono::fsi::ChFluidSystemSPH::GetForces | ( | ) |
Extract forces applied to all markers (SPH and BCE).
The return value is a device thrust vector.
◆ GetForces() [2/2]
thrust::device_vector< Real3 > chrono::fsi::ChFluidSystemSPH::GetForces | ( | const thrust::device_vector< int > & | indices | ) |
Extract forces applied to allmarkers (SPH and BCE) with indices in the provided array.
The return value is a device thrust vector.
◆ GetParticleFluidProperties()
std::vector< ChVector3d > chrono::fsi::ChFluidSystemSPH::GetParticleFluidProperties | ( | ) | const |
Return the SPH particle fluid properties.
For each SPH particle, the 3-dimensional vector contains density, pressure, and viscosity.
◆ GetPositions() [1/2]
thrust::device_vector< Real3 > chrono::fsi::ChFluidSystemSPH::GetPositions | ( | ) |
Extract positions of all markers (SPH and BCE).
The return value is a device thrust vector.
◆ GetPositions() [2/2]
thrust::device_vector< Real3 > chrono::fsi::ChFluidSystemSPH::GetPositions | ( | const thrust::device_vector< int > & | indices | ) |
Extract positions of all markers (SPH and BCE) with indices in the provided array.
The return value is a device thrust vector.
◆ GetProperties()
thrust::device_vector< Real3 > chrono::fsi::ChFluidSystemSPH::GetProperties | ( | ) |
Extract fluid properties of all markers (SPH and BCE).
For each SPH particle, the 3-dimensional vector contains density, pressure, and viscosity.
◆ GetVelocities() [1/2]
thrust::device_vector< Real3 > chrono::fsi::ChFluidSystemSPH::GetVelocities | ( | ) |
Extract velocities of all markers (SPH and BCE).
The return value is a device thrust vector.
◆ GetVelocities() [2/2]
thrust::device_vector< Real3 > chrono::fsi::ChFluidSystemSPH::GetVelocities | ( | const thrust::device_vector< int > & | indices | ) |
Extract velocities of all markers (SPH and BCE) with indices in the provided array.
The return value is a device thrust vector.
◆ SaveParticleData()
void chrono::fsi::ChFluidSystemSPH::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::ChFluidSystemSPH::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::ChFluidSystemSPH::SetActiveDomain | ( | const ChVector3d & | boxHalfDim | ) |
Set half-dimensions of the active domain.
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 actual (CFD) simulations, but rather oinly when Chrono::FSI is used for continuum representation of granular dynamics (in terramechanics).
◆ SetActiveDomainDelay()
void chrono::fsi::ChFluidSystemSPH::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.
◆ SetBcePattern1D()
void chrono::fsi::ChFluidSystemSPH::SetBcePattern1D | ( | BcePatternMesh1D | pattern, |
bool | remove_center = false |
||
) |
Set the BCE marker pattern for 1D flexible solids for subsequent calls to AddFsiMesh.
By default, a full set of BCE markers is used across each section, including a central marker.
- Parameters
-
pattern marker pattern in cross-section remove_center eliminate markers on center line
◆ SetBcePattern2D()
void chrono::fsi::ChFluidSystemSPH::SetBcePattern2D | ( | BcePatternMesh2D | pattern, |
bool | remove_center = false |
||
) |
Set the BCE marker pattern for 2D flexible solids for subsequent calls to AddFsiMesh.
By default, BCE markers are created centered on the mesh surface, with a layer of BCEs on the surface.
- Parameters
-
pattern pattern of marker locations along normal remove_center eliminate markers on surface
◆ SetBodyForce()
void chrono::fsi::ChFluidSystemSPH::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.
◆ SetElasticSPH()
void chrono::fsi::ChFluidSystemSPH::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.
◆ SetKernelMultiplier()
void chrono::fsi::ChFluidSystemSPH::SetKernelMultiplier | ( | double | multiplier | ) |
Set multiplier for interaction length.
h = multiplier * initial_spacing.
◆ SetOutputLevel()
void chrono::fsi::ChFluidSystemSPH::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/ChFluidSystemSPH.h
- /builds/uwsbel/chrono/src/chrono_fsi/sph/ChFluidSystemSPH.cpp