Modeling utilities

Description

Utilities for creating BCE particles, setting parameters via a JSON file, and output data into files with specified format.

Collaboration diagram for Modeling utilities:

Classes

class  chrono::fsi::sph::GpuTimer
 Time recorder for cuda events. More...
 

Macros

#define mF2   make_float2
 
#define mF3   make_float3
 
#define mF4   make_float4
 
#define mR2   make_Real2
 
#define mR3   make_Real3
 
#define mR4   make_Real4
 
#define mI2   make_int2
 
#define mI3   make_int3
 
#define mI4   make_int4
 
#define mU3   make_uint3
 
#define F1CAST(x)   (float*)thrust::raw_pointer_cast(&x[0])
 
#define D1CAST(x)   (double*)thrust::raw_pointer_cast(&x[0])
 
#define BCAST(x)   (bool*)thrust::raw_pointer_cast(&x[0])
 
#define I1CAST(x)   (int*)thrust::raw_pointer_cast(&x[0])
 
#define mI2CAST(x)   (int2*)thrust::raw_pointer_cast(&x[0])
 
#define mI4CAST(x)   (int4*)thrust::raw_pointer_cast(&x[0])
 
#define U1CAST(x)   (uint*)thrust::raw_pointer_cast(&x[0])
 
#define U2CAST(x)   (uint2*)thrust::raw_pointer_cast(&x[0])
 
#define U3CAST(x)   (uint3*)thrust::raw_pointer_cast(&x[0])
 
#define U4CAST(x)   (uint4*)thrust::raw_pointer_cast(&x[0])
 
#define LU1CAST(x)   (unsigned long int*)thrust::raw_pointer_cast(&x[0])
 
#define R1CAST(x)   (Real*)thrust::raw_pointer_cast(&x[0])
 
#define mR2CAST(x)   (Real2*)thrust::raw_pointer_cast(&x[0])
 
#define mR3CAST(x)   (Real3*)thrust::raw_pointer_cast(&x[0])
 
#define mR4CAST(x)   (Real4*)thrust::raw_pointer_cast(&x[0])
 
#define TCAST(x)   thrust::raw_pointer_cast(x.data())
 
#define mR3BY3CAST(x)   (Real3By3*)thrust::raw_pointer_cast(&x[0])
 
#define CUDA_KERNEL_DIM(...)   << <__VA_ARGS__>>>
 
#define cudaCheckError()
 
#define cudaMallocErrorFlag(error_flag_D)   { cudaMalloc((void**)&error_flag_D, sizeof(bool)); }
 
#define cudaFreeErrorFlag(error_flag_D)   { cudaFree(error_flag_D); }
 
#define cudaResetErrorFlag(error_flag_D)
 
#define cudaCheckErrorFlag(error_flag_D, kernel_name)
 

Functions

void chrono::fsi::sph::computeGridSize (uint n, uint blockSize, uint &numBlocks, uint &numThreads)
 Compute number of blocks and threads for calculation on GPU. More...
 
void chrono::fsi::sph::SaveParticleDataCFD (const std::string &dir, OutputLevel level, const thrust::device_vector< Real4 > &posRadD, const thrust::device_vector< Real3 > &velMasD, const thrust::device_vector< Real4 > &derivVelRhoD, const thrust::device_vector< Real4 > &rhoPresMuD, const thrust::device_vector< Real4 > &srTauMuD, const thrust::host_vector< int4 > &referenceArray, const thrust::host_vector< int4 > &referenceArrayFEA)
 Save current CFD SPH data to files. More...
 
void chrono::fsi::sph::SaveParticleDataCRM (const std::string &dir, OutputLevel level, const thrust::device_vector< Real4 > &posRadD, const thrust::device_vector< Real3 > &velMasD, const thrust::device_vector< Real4 > &derivVelRhoD, const thrust::device_vector< Real4 > &rhoPresMuD, const thrust::device_vector< Real3 > &tauXxYyZzD, const thrust::device_vector< Real3 > &tauXyXzYzD, const thrust::host_vector< int4 > &referenceArray, const thrust::host_vector< int4 > &referenceArrayFEA)
 Save current CRM SPH data to files. More...
 
void chrono::fsi::sph::SaveSolidData (const std::string &dir, double time, const thrust::device_vector< Real3 > &posRigidD, const thrust::device_vector< Real4 > &rotRigidD, const thrust::device_vector< Real3 > &velRigidD, const thrust::device_vector< Real3 > &forceRigidD, const thrust::device_vector< Real3 > &torqueRigidD, const thrust::device_vector< Real3 > &pos1DNodeD, const thrust::device_vector< Real3 > &vel1DNodeD, const thrust::device_vector< Real3 > &force1DNodeD, const thrust::device_vector< Real3 > &pos2DNodeD, const thrust::device_vector< Real3 > &vel2DNodeD, const thrust::device_vector< Real3 > &force2DNodeD)
 Save current FSI solid data. More...
 
void chrono::fsi::sph::WriteParticleFileCSV (const std::string &filename, thrust::device_vector< Real4 > &posRadD, thrust::device_vector< Real3 > &velMasD, thrust::device_vector< Real4 > &rhoPresMuD, thrust::host_vector< int4 > &referenceArray)
 Save current particle data to a CSV file. More...
 
void chrono::fsi::sph::WriteParticleFileCHPF (const std::string &filename, thrust::device_vector< Real4 > &posRadD, thrust::host_vector< int4 > &referenceArray)
 Save current particle data to a binary ChPF file. More...
 
void chrono::fsi::sph::printStruct (struct Real2 &s)
 Print a Real2 struct.
 
void chrono::fsi::sph::printStruct (struct int2 &s)
 Print a Int2 struct.
 
void chrono::fsi::sph::printStruct (struct Real3 &s)
 Print a Real3 struct.
 
void chrono::fsi::sph::printStruct (struct int3 &s)
 Print a Int3 struct.
 
void chrono::fsi::sph::printStruct (struct Real4 &s)
 Print a Real4 struct.
 
void chrono::fsi::sph::printStruct (struct int4 &s)
 Print a Int4 struct.
 
ChVector3d chrono::fsi::sph::ToChVector (const Real3 &p3)
 Convert a Real3 data structure to a ChVector3d data structure.
 
ChVector3d chrono::fsi::sph::ToChVector (const Real2 &p2)
 Convert a Real2 data structure to a ChVector3d data structure.
 
ChVector3d chrono::fsi::sph::ToChVector (const Real4 &p4)
 Convert the first 3 arguments of a Real4 data structure to a ChVector3d data structure.
 
ChQuaternion chrono::fsi::sph::ToChQuaternion (const Real4 &q4)
 Convert a Real4 data structure to a ChQuaternion data structure.
 
Real2 chrono::fsi::sph::ToReal2 (const ChVector2<> &v2)
 
Real3 chrono::fsi::sph::ToReal3 (const ChVector3<> &v3)
 Convert a ChVector data structure to a Real3 data structure.
 
Real4 chrono::fsi::sph::ToReal4 (const ChVector3d &v3, Real m)
 Convert a ChVector3d and a scalar to a Real4 data structure.
 
Real4 chrono::fsi::sph::ToReal4 (const ChQuaternion<> &q4)
 Convert a ChQuaternion data structure to a Real4 data structure.
 

Macro Definition Documentation

◆ cudaCheckError

#define cudaCheckError ( )
Value:
{ \
cudaError_t e = cudaGetLastError(); \
if (e != cudaSuccess) { \
printf("Cuda failure %s:%d: '%s'\n", __FILE__, __LINE__, cudaGetErrorString(e)); \
exit(0); \
} \
}

◆ cudaCheckErrorFlag

#define cudaCheckErrorFlag (   error_flag_D,
  kernel_name 
)
Value:
{ \
bool error_flag_H; \
cudaDeviceSynchronize(); \
cudaMemcpy(&error_flag_H, error_flag_D, sizeof(bool), cudaMemcpyDeviceToHost); \
if (error_flag_H) { \
char buffer[256]; \
sprintf(buffer, "Error flag intercepted in %s:%d from %s", __FILE__, __LINE__, kernel_name); \
printf("%s\n", buffer); \
throw std::runtime_error(buffer); \
} \
cudaError_t e = cudaGetLastError(); \
if (e != cudaSuccess) { \
char buffer[256]; \
sprintf(buffer, "CUDA failure in %s:%d from %s", __FILE__, __LINE__, cudaGetErrorString(e)); \
printf("%s\n", buffer); \
throw std::runtime_error(buffer); \
} \
}

◆ cudaResetErrorFlag

#define cudaResetErrorFlag (   error_flag_D)
Value:
{ \
bool error_flag_H = false; \
cudaMemcpy(error_flag_D, &error_flag_H, sizeof(bool), cudaMemcpyHostToDevice); \
}

Function Documentation

◆ computeGridSize()

void chrono::fsi::sph::computeGridSize ( uint  n,
uint  blockSize,
uint numBlocks,
uint numThreads 
)

Compute number of blocks and threads for calculation on GPU.

This function calculates the number of blocks and threads for a given number of elements based on the blockSize.

Parameters
ntotal number of elements
blockSizeblock size (threads per block)
numBlocksnumber of blocks [output]
numThreadsnumber of threads [output]

◆ SaveParticleDataCFD()

void chrono::fsi::sph::SaveParticleDataCFD ( const std::string &  dir,
OutputLevel  level,
const thrust::device_vector< Real4 > &  posRadD,
const thrust::device_vector< Real3 > &  velMasD,
const thrust::device_vector< Real4 > &  derivVelRhoD,
const thrust::device_vector< Real4 > &  rhoPresMuD,
const thrust::device_vector< Real4 > &  srTauMuD,
const thrust::host_vector< int4 > &  referenceArray,
const thrust::host_vector< int4 > &  referenceArrayFEA 
)

Save current CFD SPH data to files.

Create separate files to write fluid, boundary BCE, rigid BCE, and flex BCE marker information. The amount of data saved for each marker is controlled by the specified OutputLevel (e.g., for level CFD_FULL, the output includes shear rate).

◆ SaveParticleDataCRM()

void chrono::fsi::sph::SaveParticleDataCRM ( const std::string &  dir,
OutputLevel  level,
const thrust::device_vector< Real4 > &  posRadD,
const thrust::device_vector< Real3 > &  velMasD,
const thrust::device_vector< Real4 > &  derivVelRhoD,
const thrust::device_vector< Real4 > &  rhoPresMuD,
const thrust::device_vector< Real3 > &  tauXxYyZzD,
const thrust::device_vector< Real3 > &  tauXyXzYzD,
const thrust::host_vector< int4 > &  referenceArray,
const thrust::host_vector< int4 > &  referenceArrayFEA 
)

Save current CRM SPH data to files.

Create separate files to write fluid, boundary BCE, rigid BCE, and flex BCE marker information. The amount of data saved for each marker is controlled by the specified OutputLevel (e.g., for level CRM_FULL, the output includes stress).

◆ SaveSolidData()

void chrono::fsi::sph::SaveSolidData ( const std::string &  dir,
double  time,
const thrust::device_vector< Real3 > &  posRigidD,
const thrust::device_vector< Real4 > &  rotRigidD,
const thrust::device_vector< Real3 > &  velRigidD,
const thrust::device_vector< Real3 > &  forceRigidD,
const thrust::device_vector< Real3 > &  torqueRigidD,
const thrust::device_vector< Real3 > &  pos1DNodeD,
const thrust::device_vector< Real3 > &  vel1DNodeD,
const thrust::device_vector< Real3 > &  force1DNodeD,
const thrust::device_vector< Real3 > &  pos2DNodeD,
const thrust::device_vector< Real3 > &  vel2DNodeD,
const thrust::device_vector< Real3 > &  force2DNodeD 
)

Save current FSI solid data.

Append states and fluid forces at current time for all solids in the FSI problem.

◆ WriteParticleFileCHPF()

void chrono::fsi::sph::WriteParticleFileCHPF ( const std::string &  filename,
thrust::device_vector< Real4 > &  posRadD,
thrust::host_vector< int4 > &  referenceArray 
)

Save current particle data to a binary ChPF file.

Write particle positions only.

◆ WriteParticleFileCSV()

void chrono::fsi::sph::WriteParticleFileCSV ( const std::string &  filename,
thrust::device_vector< Real4 > &  posRadD,
thrust::device_vector< Real3 > &  velMasD,
thrust::device_vector< Real4 > &  rhoPresMuD,
thrust::host_vector< int4 > &  referenceArray 
)

Save current particle data to a CSV file.

Write particle positions, velocities, rho, pressure, and mu.