Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Feature/multiphase material #675

Merged
merged 14 commits into from
Jul 25, 2020
Merged
Show file tree
Hide file tree
Changes from 3 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
25 changes: 13 additions & 12 deletions include/mesh.tcc
Original file line number Diff line number Diff line change
Expand Up @@ -476,8 +476,9 @@ bool mpm::Mesh<Tdim>::generate_material_points(
unsigned before_generation = this->nparticles();
bool checks = false;
// Get material
std::vector<std::shared_ptr<mpm::Material<Tdim>>> material;
for (auto m_id : material_ids) material.emplace_back(materials_.at(m_id));
std::vector<std::shared_ptr<mpm::Material<Tdim>>> materials;
for (auto m_id : material_ids)
materials.emplace_back(materials_.at(m_id));

// If set id is -1, use all cells
auto cset = (cset_id == -1) ? this->cells_ : cell_sets_.at(cset_id);
Expand All @@ -501,8 +502,8 @@ bool mpm::Mesh<Tdim>::generate_material_points(
status = this->add_particle(particle, checks);
if (status) {
map_particles_[pid]->assign_cell(*citr);
for (unsigned phase = 0; phase < material.size(); phase++)
map_particles_[pid]->assign_material(material[phase], phase);
for (unsigned phase = 0; phase < materials.size(); phase++)
map_particles_[pid]->assign_material(materials[phase], phase);
pids.emplace_back(pid);
} else
throw std::runtime_error("Generate particles in mesh failed");
Expand Down Expand Up @@ -543,8 +544,8 @@ bool mpm::Mesh<Tdim>::create_particles(
// Particle ids
std::vector<mpm::Index> pids;
// Get material
std::vector<std::shared_ptr<mpm::Material<Tdim>>> material;
for (auto m_id : material_ids) material.emplace_back(materials_.at(m_id));
std::vector<std::shared_ptr<mpm::Material<Tdim>>> materials;
for (auto m_id : material_ids) materials.emplace_back(materials_.at(m_id));
// Check if particle coordinates is empty
if (coordinates.empty())
throw std::runtime_error("List of coordinates is empty");
Expand All @@ -563,8 +564,8 @@ bool mpm::Mesh<Tdim>::create_particles(

// If insertion is successful
if (insert_status) {
for (unsigned phase = 0; phase < material.size(); phase++)
map_particles_[pid]->assign_material(material[phase], phase);
for (unsigned phase = 0; phase < materials.size(); phase++)
map_particles_[pid]->assign_material(materials[phase], phase);
pids.emplace_back(pid);
} else
throw std::runtime_error("Addition of particle to mesh failed!");
Expand Down Expand Up @@ -1724,9 +1725,9 @@ void mpm::Mesh<Tdim>::inject_particles(double current_time) {
unsigned pid = this->nparticles();
bool checks = false;
// Get material
std::vector<std::shared_ptr<mpm::Material<Tdim>>> material;
std::vector<std::shared_ptr<mpm::Material<Tdim>>> materials;
for (auto m_id : injection.material_ids)
material.emplace_back(materials_.at(m_id));
materials.emplace_back(materials_.at(m_id));

// Check if duration is within the current time
if (injection.start_time <= current_time &&
Expand Down Expand Up @@ -1760,8 +1761,8 @@ void mpm::Mesh<Tdim>::inject_particles(double current_time) {
unsigned status = this->add_particle(particle, checks);
if (status) {
map_particles_[pid]->assign_cell(*citr);
for (unsigned phase = 0; phase < material.size(); phase++)
map_particles_[pid]->assign_material(material[phase], phase);
for (unsigned phase = 0; phase < materials.size(); phase++)
map_particles_[pid]->assign_material(materials[phase], phase);
++pid;
injected_particles.emplace_back(particle);
}
Expand Down
12 changes: 12 additions & 0 deletions include/particles/particle.h
Original file line number Diff line number Diff line change
Expand Up @@ -289,8 +289,20 @@ class Particle : public ParticleBase<Tdim> {
//! Return neighbour ids
std::vector<mpm::Index> neighbours() const override { return neighbours_; };

protected:
//! Initialise particle material container
//! \details This function allocate memory and initialise the material related
//! containers according to the particle phase, i.e. solid or fluid particle
//! has phase_size = 1, whereas two-phase (solid-fluid) or three-phase
//! (solid-water-air) particle have phase_size = 2 and 3, respectively.
//! \param[in] phase_size The material phase size
void initialise_material(unsigned phase_size = 1);

private:
//! Compute strain rate
//! \param[in] dn_dx The spatial gradient of shape function
//! \param[in] phase Index to indicate phase
//! \retval strain rate at particle inside a cell
inline Eigen::Matrix<double, 6, 1> compute_strain_rate(
const Eigen::MatrixXd& dn_dx, unsigned phase) noexcept;

Expand Down
30 changes: 20 additions & 10 deletions include/particles/particle.tcc
Original file line number Diff line number Diff line change
Expand Up @@ -7,10 +7,8 @@ mpm::Particle<Tdim>::Particle(Index id, const VectorDim& coord)
cell_ = nullptr;
// Nodes
nodes_.clear();
// Set material pointer to null
material_.emplace_back(nullptr);
material_id_.emplace_back(std::numeric_limits<unsigned>::max());
state_variables_.emplace_back(mpm::dense_map());
// Set material containers
this->initialise_material(1);
// Logger
std::string logger =
"particle" + std::to_string(Tdim) + "d::" + std::to_string(id);
Expand All @@ -24,9 +22,8 @@ mpm::Particle<Tdim>::Particle(Index id, const VectorDim& coord, bool status)
this->initialise();
cell_ = nullptr;
nodes_.clear();
material_.emplace_back(nullptr);
material_id_.emplace_back(std::numeric_limits<unsigned>::max());
state_variables_.emplace_back(mpm::dense_map());
// Set material containers
this->initialise_material(1);
//! Logger
std::string logger =
"particle" + std::to_string(Tdim) + "d::" + std::to_string(id);
Expand Down Expand Up @@ -253,6 +250,18 @@ void mpm::Particle<Tdim>::initialise() {
this->properties_["displacements"] = [&]() { return displacement(); };
}

//! Initialise particle material container
template <unsigned Tdim>
void mpm::Particle<Tdim>::initialise_material(unsigned phase_size) {
material_.resize(phase_size);
material_id_.resize(phase_size);
state_variables_.resize(phase_size);
std::fill(material_.begin(), material_.end(), nullptr);
std::fill(material_id_.begin(), material_id_.end(),
std::numeric_limits<unsigned>::max());
std::fill(state_variables_.begin(), state_variables_.end(), mpm::dense_map());
}

//! Assign material state variables from neighbour particle
template <unsigned Tdim>
bool mpm::Particle<Tdim>::assign_material_state_vars(
Expand Down Expand Up @@ -386,9 +395,10 @@ bool mpm::Particle<Tdim>::assign_material(
try {
// Check if material is valid and properties are set
if (material != nullptr) {
material_[phase] = material;
material_id_[phase] = material_[phase]->id();
state_variables_[phase] = material_[phase]->initialise_state_variables();
material_.at(phase) = material;
material_id_.at(phase) = material_[phase]->id();
state_variables_.at(phase) =
material_[phase]->initialise_state_variables();
status = true;
} else {
throw std::runtime_error("Material is undefined!");
Expand Down
4 changes: 2 additions & 2 deletions tests/mesh_test_2d.cc
Original file line number Diff line number Diff line change
Expand Up @@ -578,7 +578,7 @@ TEST_CASE("Mesh is checked for 2D case", "[mesh][2D]") {
// Gauss point generation
Json jgen;
jgen["type"] = "gauss";
jgen["material_id"] = mid;
jgen["material_id"] = {0};
kks32 marked this conversation as resolved.
Show resolved Hide resolved
jgen["cset_id"] = 1;
jgen["particle_type"] = "P2D";
jgen["check_duplicates"] = false;
Expand All @@ -595,7 +595,7 @@ TEST_CASE("Mesh is checked for 2D case", "[mesh][2D]") {
// Gauss point generation
Json jgen;
jgen["type"] = "inject";
jgen["material_id"] = mid;
jgen["material_id"] = {0};
jgen["cset_id"] = 1;
jgen["particle_type"] = "P2D";
jgen["check_duplicates"] = false;
Expand Down