![]() |
Karamelo
714599e9
Parallel Material Point Method Simulator
|
This class represents a given solid.
In this class is stored the variables related to the solid boundaries as well as all those related to the particles it contains. Moreover, it holds the evaluation of the shape functions and its derivatives. The functions member of this class are related to the Particle to Grid (P2G) and Grid to Particle (G2P) steps of the MPM algorithm. Another member is a function that updates the particles' stresses. This class holds a pointer of type Grid that points to either the local grid, as used in the Total Lagrangian MPM, or the global grid, as used in the Updated Lagrangian MPM. This class has Pointers as parent in order to be accessible from anywhere within the MPM class.
Public Member Functions | |
void | compute_deformation_gradient () |
Compute the deformation gradient directly from the grid nodes' positions. | |
void | compute_external_forces_nodes (bool) |
Compute external forces step of the Particle to Grid step of the MPM algorithm. | |
void | compute_inertia_tensor (string) |
Compute the inertia tensor necessary for the Affice PIC. | |
void | compute_internal_forces_nodes_TL () |
Compute internal forces step of the Particle to Grid step of the total Lagrangian MPM algorithm. | |
void | compute_internal_forces_nodes_UL (bool) |
Compute internal forces step of the Particle to Grid step of the updated Lagrangian MPM algorithm. | |
void | compute_mass_nodes (bool) |
Compute nodal mass step of the Particle to Grid step of the MPM algorithm. | |
void | compute_particle_acceleration () |
Update the particles' acceleration. | |
void | compute_particle_velocities_and_positions () |
Compute the particles' temporary velocities and position, part of the Grid to Particles step of the MPM algorithm. | |
void | compute_rate_deformation_gradient_TL () |
Compute the time derivative of the deformation matrix for TLMPM, when APIC is not used. | |
void | compute_rate_deformation_gradient_TL_APIC () |
Compute the time derivative of the deformation matrix for TLMPM, when APIC is used. | |
void | compute_rate_deformation_gradient_UL_APIC () |
Compute the time derivative of the deformation matrix for TLMPM, when APIC is in use. | |
void | compute_rate_deformation_gradient_UL_MUSL () |
Compute the time derivative of the deformation matrix for TLMPM, when using Modified Update Stress Last and APIC is not used. | |
void | compute_rate_deformation_gradient_UL_USL () |
Compute the time derivative of the deformation matrix for TLMPM, when using Update Stress Last and APIC is not used. | |
void | compute_velocity_nodes (bool) |
Compute nodal velocity (via momentum) step of the Particle to Grid step of the MPM algorithm. | |
void | compute_velocity_nodes_APIC (bool) |
Specific function that computes the nodal velocity (via momentum) when using Affine PIC (APIC). | |
void | copy_particle (int, int) |
void | grow (int) |
Allocate memory for the vectors used for particles or resize them. | |
void | init () |
Launch the initialization of the grid. | |
void | options (vector< string > *, vector< string >::iterator) |
Determines the material and temperature schemes used. | |
void | pack_particle (int, vector< double > &) |
Pack particles attributes into a buffer (used for generating a restart). | |
Solid (class MPM *, vector< string >) | |
The main task of the constructor is to read the input arguments and launch the creation of particles. | |
void | unpack_particle (int &, vector< int >, double[]) |
Unpack particles attributes from a buffer (used when reading a restart). | |
void | update_deformation_gradient () |
Update the deformation gradient, volume, density, and the necessary strain matrices. | |
void | update_particle_domain () |
Update the particle domain. Used with CPDI. | |
void | update_particle_velocities (double) |
void | update_stress () |
Calculate the stress, damage and temperature at each particle, and determine the maximum allowed time step. | |
Public Attributes | |
vector< Eigen::Vector3d > | a |
Particles' acceleration. | |
int | comm_n |
Number of double to pack for particle exchange between CPU. | |
vector< Eigen::Matrix3d > | D |
Symmetric part of L. | |
vector< double > | damage |
Particles' damage variable. | |
vector< double > | damage_init |
Particles' damage initiation variable. | |
vector< Eigen::Matrix3d > | Di |
Inertia tensor. | |
vector< double > | eff_plastic_strain |
Particles' effective plastic strain. | |
vector< double > | eff_plastic_strain_rate |
Particles' effective plastic strain rate. | |
vector< Eigen::Vector3d > | f |
Particles' internal forces. | |
vector< Eigen::Matrix3d > | F |
Deformation gradient matrix. | |
vector< Eigen::Matrix3d > | Fdot |
Rate of deformation gradient matrix. | |
vector< Eigen::Matrix3d > | Finv |
Inverse of the deformation gradient matrix. | |
class Grid * | grid |
Pointer to the background grid. | |
string | id |
Solid id. | |
vector< double > | ienergy |
Particles' internal energy. | |
vector< double > | J |
Determinant of the deformation matrix. | |
vector< Eigen::Matrix3d > | L |
Velocity gradient matrix. | |
vector< int > | mask |
Particles' group mask. | |
vector< double > | mass |
Particles' current mass. | |
struct Mat * | mat |
Pointer to the material. | |
double | max_p_wave_speed |
Maximum of the particle wave speed. | |
vector< Eigen::Vector3d > | mbp |
Particles' external forces times mass. | |
string | method_type |
Either tlmpm, tlcpdi, tlcpdi2, ulmpm, ulcpdi, or ulcpdi2 (all kinds of MPM supported) | |
double | mtot |
Total mass. | |
int | nc |
Number of corners per particles: \(2^{dimension}\). | |
vector< vector< int > > | neigh_np |
List of the particles neighbouring a given node. | |
vector< vector< int > > | neigh_pn |
List of the nodes neighbouring a given particle. | |
bigint | np |
Total number of particles in the domain. | |
int | np_local |
Number of local particles (in this CPU) | |
vector< int > | numneigh_np |
Number of nodes neighbouring a given node. | |
vector< int > | numneigh_pn |
Number of nodes neighbouring a given particle. | |
vector< double > | pH_regu |
Particles' regularized hydrostatic pressure. | |
vector< tagint > | ptag |
Unique identifier for particles in the system. | |
vector< Eigen::Matrix3d > | R |
Rotation matrix. | |
vector< double > | rho |
Particles' current density. | |
vector< double > | rho0 |
Particles' reference density. | |
vector< Eigen::Vector3d > | rp |
Current domain vector (CPDI1) | |
vector< Eigen::Vector3d > | rp0 |
Reference domain vector (CPDI1) | |
vector< Eigen::Matrix3d > | sigma |
Stress matrix. | |
double | solidhi [3] |
Solid global bounds. | |
double | solidsubhi [3] |
Solid local bounds. | |
vector< Eigen::Matrix3d > | strain_el |
Elastic strain matrix. | |
vector< double > | T |
Particles' temperature. | |
vector< Eigen::Vector3d > | v |
Particles' current velocity. | |
vector< Eigen::Vector3d > | v_update |
Particles' velocity at time t+dt. | |
vector< double > | vol |
Particles' current volume. | |
vector< double > | vol0 |
Particles' reference volume. | |
vector< Eigen::Matrix3d > | vol0PK1 |
Transpose of the 1st Piola-Kirchhoff matrix times vol0. | |
double | vtot |
Total volume. | |
vector< vector< double > > | wf_np |
Array of arrays (matrix) of the weight functions \(\Phi_{Ip}\) effectively the transpose of wf_pn. | |
vector< vector< double > > | wf_pn |
Array of arrays (matrix) of the weight functions \(\Phi_{pI}\). | |
vector< vector< double > > | wf_pn_corners |
Array of arrays (matrix) of the weight functions \(\Phi_{Ic}\) evaluated at the corners of the particle's domain (used in CPDI) | |
vector< vector< Eigen::Vector3d > > | wfd_np |
Array of arrays (matrix) of the derivative of the weight functions \(\partial \Phi_{Ip}/ \partial x\) effectively the transpose of wfd_pn. | |
vector< vector< Eigen::Vector3d > > | wfd_pn |
Array of arrays (matrix) of the derivative of the weight functions \(\partial \Phi_{pI}/\partial x\). | |
vector< Eigen::Vector3d > | x |
Particles' current position. | |
vector< Eigen::Vector3d > | x0 |
Particles' reference position. | |
vector< Eigen::Vector3d > | xpc |
Current position of the corners of the particles' domain (CPDI2o) | |
vector< Eigen::Vector3d > | xpc0 |
Reference position of the corners of the particles' domain (CPDI2) | |
void Solid::copy_particle | ( | int | i, |
int | j | ||
) |
Copy particle i attribute and copy it to particle j. This function is used to re-order the memory arrangment of particles. Usually this is done when particle j is deleted.
void Solid::update_particle_velocities | ( | double | alpha | ) |
Update the particles' velocities based on either PIC and/or FLIP. The argument is the ratio \(\alpha\) used between PIC and FLIP. \(\alpha = 0\) for pure PIC, \(\alpha = 1\) for pure FLIP.