PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
Loading...
Searching...
No Matches
material::PmbMaterial Class Reference

A class providing methods to compute energy density and force of peridynamic material. More...

#include <material.h>

Inheritance diagram for material::PmbMaterial:
Collaboration diagram for material::PmbMaterial:

Public Member Functions

 PmbMaterial (inp::MaterialDeck &deck, const size_t &dim, const double &horizon)
 Constructor.
 
bool isStateActive () const override
 Returns true if state-based potential is active.
 
std::pair< double, double > getBondEF (const double &r, const double &s, bool &fs, const bool &break_bonds) const override
 Returns energy and force between bond due to pairwise interaction.
 
std::pair< double, double > getBondEF (const double &r, const double &s, bool &fs, const double &mx, const double &thetax) const override
 Returns energy and force between bond due to state-based model.
 
util::Point getBondForceDirection (const util::Point &dx, const util::Point &du) const override
 Returns the unit vector along which bond-force acts.
 
double getS (const util::Point &dx, const util::Point &du) const override
 Returns the bond strain.
 
double getSc (const double &r) const override
 Returns critical bond strain.
 
double getDensity () const override
 Returns the density of the material.
 
double getInfFn (const double &r) const override
 Returns the value of influence function.
 
double getMoment (const size_t &i) const override
 Returns the moment of influence function.
 
double getHorizon () const override
 Returns horizon.
 
inp::MatData computeMaterialProperties (const size_t &dim) const override
 Computes elastic and fracture material properties and returns the data.
 
std::string printStr (int nt, int lvl) const override
 Returns the string containing printable information about the object.
 
void print (int nt, int lvl) const override
 Prints the information about the object.
 
void print () const override
 Prints the information about the object.
 
- Public Member Functions inherited from material::Material
 Material (std::string name="")
 Constructor.
 
virtual ~Material ()
 Destructor.
 
std::string name ()
 Returns name of the material.
 
size_t getDimension ()
 Returns dimension of the problem.
 
bool isPlaneStrain ()
 Returns plane-strain condition.
 

Private Member Functions

void computeParameters (inp::MaterialDeck &deck, const size_t &dim)
 Compute material model parameters.
 

Private Attributes

double d_horizon
 Horizon.
 
double d_density
 Density.
 
Material parameters
double d_c
 Parameter C.
 
double d_s0
 Parameter \( \beta \).
 

Detailed Description

A class providing methods to compute energy density and force of peridynamic material.

Definition at line 665 of file material.h.

Constructor & Destructor Documentation

◆ PmbMaterial()

material::PmbMaterial::PmbMaterial ( inp::MaterialDeck deck,
const size_t &  dim,
const double &  horizon 
)
inline

Constructor.

Parameters
deckInput deck which contains user-specified information
dimDimension
horizonHorizon

Definition at line 674 of file material.h.

675 : Material("PMBBond"), d_horizon(horizon), d_density(deck.d_density),
676 d_c(0.), d_s0(0.) {
677
678 // set global fields
679 if (dimension != dim)
680 dimension = dim;
681
682 if (is_plane_strain != deck.d_isPlaneStrain)
684
685 // create influence function
686 if (deck.d_influenceFnType == 0) {
687 if (influence_fn == nullptr)
688 influence_fn = std::make_shared<material::ConstInfluenceFn>(
689 deck.d_influenceFnParams, dim);
690 }
691 else if (deck.d_influenceFnType == 1) {
692 if (influence_fn == nullptr)
693 influence_fn = std::make_shared<material::LinearInfluenceFn>(
694 deck.d_influenceFnParams, dim);
695 }
696 else if (deck.d_influenceFnType == 2) {
697 if (influence_fn == nullptr)
698 influence_fn = std::make_shared<material::GaussianInfluenceFn>(
699 deck.d_influenceFnParams, dim);
700 }
701 else {
702 std::cerr << "Error: Influence function type = "
703 << deck.d_influenceFnType
704 << " is invalid.\n";
705 exit(1);
706 }
707
708 // check if we need to compute the material parameters
710 computeParameters(deck, dim);
711 else {
712 d_c = deck.d_bondPotentialParams[0];
713 d_s0 = deck.d_bondPotentialParams[1];
714 }
715 };
Material(std::string name="")
Constructor.
Definition material.h:83
double d_density
Density.
Definition material.h:1005
void computeParameters(inp::MaterialDeck &deck, const size_t &dim)
Compute material model parameters.
Definition material.h:919
double d_s0
Parameter .
Definition material.h:1016
double d_horizon
Horizon.
Definition material.h:1002
double d_c
Parameter C.
Definition material.h:1013
bool is_plane_strain
Is plane-stress condition active.
Definition material.h:29
size_t dimension
Dimension of the domain.
Definition material.h:26
std::shared_ptr< material::BaseInfluenceFn > influence_fn
Store pointer to influence function globally.
Definition material.h:32
bool d_isPlaneStrain
Indicates if the 2-d simulation is of plane-strain type (thick material) or plane-stress type (thin m...
double d_density
Density of material.
size_t d_influenceFnType
Type of influence function.
std::vector< double > d_bondPotentialParams
List of parameters for pairwise potential.
bool d_computeParamsFromElastic
Compute Peridynamic material properties from elastic properties.
std::vector< double > d_influenceFnParams
List of parameters for influence function.

References computeParameters(), inp::MaterialDeck::d_bondPotentialParams, d_c, inp::MaterialDeck::d_computeParamsFromElastic, inp::MaterialDeck::d_influenceFnParams, inp::MaterialDeck::d_influenceFnType, inp::MaterialDeck::d_isPlaneStrain, and d_s0.

Here is the call graph for this function:

Member Function Documentation

◆ computeMaterialProperties()

inp::MatData material::PmbMaterial::computeMaterialProperties ( const size_t &  dim) const
inlineoverridevirtual

Computes elastic and fracture material properties and returns the data.

Parameters
dimDimension of the problem
Returns
inp::MatData Material data

Implements material::Material.

Definition at line 846 of file material.h.

846 {
847
848 auto data = inp::MatData();
849
850 // set Poisson's ratio to 1/4
851 data.d_nu = 0.25;
852
853 // compute peridynamic parameters
854 if (dim == 2) {
855 data.d_lambda = d_c * (M_PI * std::pow(d_horizon, 3.0)) / 24.0;
856 data.d_E = data.toELambda(data.d_lambda, data.d_nu);
857 data.d_Gc = d_s0 * d_s0 * (9.0 * data.d_E * d_horizon) / (5.0 * M_PI);
858 } else if (dim == 3) {
859 data.d_lambda = d_c * (M_PI * std::pow(d_horizon, 3.0)) / 24.0;
860 data.d_E = data.toELambda(data.d_lambda, data.d_nu);
861 data.d_Gc = d_s0 * d_s0 * (9.0 * data.d_E * d_horizon) / (5.0 * M_PI);
862 }
863 data.d_mu = data.d_lambda;
864 data.d_G = data.d_lambda;
865 data.d_K =
866 data.toK(data.d_E, data.d_nu);
867 data.d_KIc = data.toKIc(
868 data.d_Gc, data.d_nu, data.d_E);
869
870 return data;
871 };
Structure for elastic properties and fracture properties.

References d_c, d_horizon, and d_s0.

◆ computeParameters()

void material::PmbMaterial::computeParameters ( inp::MaterialDeck deck,
const size_t &  dim 
)
inlineprivate

Compute material model parameters.

Parameters
deckMaterialDeck
dimDimension of the domain

Definition at line 919 of file material.h.

919 {
920 //
921 // Need following elastic and fracture properties
922 // 1. E or K
923 // 2. Gc or KIc
924 // For bond-based, Poisson's ratio is fixed to 1/4
925 //
926 if (util::isLess(deck.d_matData.d_E, 0.) &&
927 util::isLess(deck.d_matData.d_K, 0.)) {
928 std::cerr << "Error: Require either Young's modulus E or Bulk modulus K"
929 " to compute the RNP bond-based peridynamic parameters.\n";
930 exit(1);
931 }
932 if (util::isGreater(deck.d_matData.d_E, 0.) &&
933 util::isGreater(deck.d_matData.d_K, 0.)) {
934 std::cout << "Warning: Both Young's modulus E and Bulk modulus K are "
935 "provided.\n";
936 std::cout << "Warning: To compute the RNP bond-based peridynamic "
937 "parameters, we only require one of those.\n";
938 std::cout
939 << "Warning: Selecting Young's modulus to compute parameters.\n";
940 }
941
942 if (util::isLess(deck.d_matData.d_Gc, 0.) &&
943 util::isLess(deck.d_matData.d_KIc, 0.)) {
944 std::cerr << "Error: Require either critical energy release rate Gc or "
945 "critical stress intensity factor KIc to compute the RNP "
946 "bond-based peridynamic parameters.\n";
947 exit(1);
948 } else if (util::isGreater(deck.d_matData.d_Gc, 0.) &&
949 util::isGreater(deck.d_matData.d_KIc, 0.)) {
950 std::cout << "Warning: Both critical energy release rate Gc and critical "
951 "stress intensity factor KIc are provided.\n";
952 std::cout << "Warning: To compute the RNP bond-based peridynamic "
953 "parameters, we only require one of those.\n";
954 std::cout << "Warning: Selecting critical energy release rate Gc to "
955 "compute parameters.\n";
956 }
957
958 // set Poisson's ratio to 1/4
959 deck.d_matData.d_nu = 0.25;
960
961 // compute E if not provided or K if not provided
962 if (deck.d_matData.d_E > 0.)
963 deck.d_matData.d_K =
964 deck.d_matData.toK(deck.d_matData.d_E, deck.d_matData.d_nu);
965
966 if (deck.d_matData.d_K > 0. && deck.d_matData.d_E < 0.)
967 deck.d_matData.d_E =
968 deck.d_matData.toE(deck.d_matData.d_K, deck.d_matData.d_nu);
969
970 if (deck.d_matData.d_Gc > 0.)
971 deck.d_matData.d_KIc = deck.d_matData.toKIc(
972 deck.d_matData.d_Gc, deck.d_matData.d_nu, deck.d_matData.d_E);
973
974 if (deck.d_matData.d_KIc > 0. && deck.d_matData.d_Gc < 0.)
975 deck.d_matData.d_Gc = deck.d_matData.toGc(
976 deck.d_matData.d_KIc, deck.d_matData.d_nu, deck.d_matData.d_E);
977
978 // compute lame parameter
979 deck.d_matData.d_lambda =
981 deck.d_matData.d_G =
982 deck.d_matData.toGE(deck.d_matData.d_E, deck.d_matData.d_nu);
983 deck.d_matData.d_mu = deck.d_matData.d_G;
984
985 // compute peridynamic parameters
986 if (dim == 2) {
987 // Ha, Bobaru 2010 "Studies of dynamic crack propagation and crack branching
988 // with peridynamics"
989 d_c = 24.0 * deck.d_matData.d_E / (M_PI * std::pow(d_horizon, 3.0) *
990 (1. - deck.d_matData.d_nu));
991 d_s0 = std::sqrt(5.0 * M_PI * deck.d_matData.d_Gc /
992 (9.0 * deck.d_matData.d_E * d_horizon));
993 } else if (dim == 3) {
994 d_c = 24.0 * deck.d_matData.d_lambda / (M_PI * std::pow(d_horizon, 3.0));
995 d_s0 = std::sqrt(5.0 * M_PI * deck.d_matData.d_Gc /
996 (9.0 * deck.d_matData.d_E * d_horizon));
997 }
998 };
bool isGreater(const double &a, const double &b)
Returns true if a > b.
Definition function.cpp:15
bool isLess(const double &a, const double &b)
Returns true if a < b.
Definition function.cpp:20
double toGc(double KIc, double nu, double E)
Compute critical energy release rate Gc from critical stress-intensity factor KIc,...
double toLambdaE(double E, double nu)
Compute Lame first parameter lambda from Young's modulus E and Poisson's ratio nu.
double d_mu
Lame second parameter.
double toKIc(double Gc, double nu, double E)
Compute critical stress-intensity factor KIc from critical energy release rate Gc,...
double d_KIc
Critical stress intensity factor.
double toK(double E, double nu)
Compute Bulk modulus K from Young's modulus K and Poisson's ratio nu.
double d_K
Bulk modulus.
double d_lambda
Lame first parameter.
double toGE(double E, double nu)
Compute shear modulus from Young's modulus E and Poisson's ratio nu.
double d_nu
Poisson's ratio.
double d_G
Shear modulus or Lame second parameter.
double d_E
Young's elastic modulus.
double d_Gc
Critical energy release rate.
double toE(double K, double nu)
Compute Young's modulus E from Bulk modulus K and Poisson's ratio nu.
inp::MatData d_matData
List of elastic and fracture properties.

References d_c, inp::MatData::d_E, inp::MatData::d_G, inp::MatData::d_Gc, d_horizon, inp::MatData::d_K, inp::MatData::d_KIc, inp::MatData::d_lambda, inp::MaterialDeck::d_matData, inp::MatData::d_mu, inp::MatData::d_nu, d_s0, util::isGreater(), util::isLess(), inp::MatData::toE(), inp::MatData::toGc(), inp::MatData::toGE(), inp::MatData::toK(), inp::MatData::toKIc(), and inp::MatData::toLambdaE().

Referenced by PmbMaterial().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ getBondEF() [1/2]

std::pair< double, double > material::PmbMaterial::getBondEF ( const double &  r,
const double &  s,
bool &  fs,
const bool &  break_bonds 
) const
inlineoverridevirtual

Returns energy and force between bond due to pairwise interaction.

Parameters
rReference (initial) bond length
sBond strain
fsBond fracture state
break_bondsFlag to specify whether bonds are allowed to break or not
Returns
value Pair of energy and force

Implements material::Material.

Definition at line 732 of file material.h.

734 {
735
736 if (!break_bonds)
737 return std::make_pair(getInfFn(r) * 0.5 * d_c * s *
738 s * r,
739 getInfFn(r) * d_c * s);
740
741 // check if fracture state of the bond need to be updated
742 if (!fs && util::isGreater(std::abs(s), d_s0 + 1.0e-10))
743 fs = true;
744
745 // if bond is not fractured, return energy and force from nonlinear
746 // potential otherwise return energy of fractured bond, and zero force
747 if (!fs)
748 return std::make_pair(getInfFn(r) * 0.5 * d_c * s *
749 s * r,
750 getInfFn(r) * d_c * s);
751 else
752 return std::make_pair(
753 getInfFn(r) * 0.5 * d_c * d_s0 * d_s0 * r, 0.);
754 };
double getInfFn(const double &r) const override
Returns the value of influence function.
Definition material.h:815

References d_c, d_s0, getInfFn(), and util::isGreater().

Referenced by getBondEF().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ getBondEF() [2/2]

std::pair< double, double > material::PmbMaterial::getBondEF ( const double &  r,
const double &  s,
bool &  fs,
const double &  mx,
const double &  thetax 
) const
inlineoverridevirtual

Returns energy and force between bond due to state-based model.

Parameters
rReference (initial) bond length
sBond strain
fsBond fracture state
mxWeighted volume at node
thetaxDilation
Returns
value Pair of energy and force

Implements material::Material.

Definition at line 767 of file material.h.

768 {
769
770 return this->getBondEF(r, s, fs, true);
771 };
std::pair< double, double > getBondEF(const double &r, const double &s, bool &fs, const bool &break_bonds) const override
Returns energy and force between bond due to pairwise interaction.
Definition material.h:732

References getBondEF().

Here is the call graph for this function:

◆ getBondForceDirection()

util::Point material::PmbMaterial::getBondForceDirection ( const util::Point dx,
const util::Point du 
) const
inlineoverridevirtual

Returns the unit vector along which bond-force acts.

Parameters
dxReference bond vector
duDifference of displacement
Returns
vector Unit vector

Implements material::Material.

Definition at line 780 of file material.h.

781 {
782 return (dx + du) / (dx + du).length();
783 };

◆ getDensity()

double material::PmbMaterial::getDensity ( ) const
inlineoverridevirtual

Returns the density of the material.

Returns
density Density of the material

Implements material::Material.

Definition at line 807 of file material.h.

807{ return d_density; };

References d_density.

◆ getHorizon()

double material::PmbMaterial::getHorizon ( ) const
inlineoverridevirtual

Returns horizon.

Returns
horizon Horizon

Implements material::Material.

Definition at line 837 of file material.h.

837{ return d_horizon; };

References d_horizon.

◆ getInfFn()

double material::PmbMaterial::getInfFn ( const double &  r) const
inlineoverridevirtual

Returns the value of influence function.

Parameters
rReference (initial) bond length
Returns
value Influence function at r

Implements material::Material.

Definition at line 815 of file material.h.

815 {
816 return getGlobalInfFn(r / d_horizon);
817 };
double getGlobalInfFn(const double &r)
Returns the value of influence function.
Definition material.h:40

References d_horizon.

Referenced by getBondEF().

Here is the caller graph for this function:

◆ getMoment()

double material::PmbMaterial::getMoment ( const size_t &  i) const
inlineoverridevirtual

Returns the moment of influence function.

If \( J(r) \) is the influence function for \( r\in [0,1)\) then \( i^{th}\) moment is given by

\[ M_i = \int_0^1 J(r) r^i dr. \]

Parameters
iith moment
Returns
value Moment

Implements material::Material.

Definition at line 828 of file material.h.

828 {
829 return getGlobalMoment(i);
830 };
double getGlobalMoment(const size_t &i)
Returns the moment of influence function.
Definition material.h:53

◆ getS()

double material::PmbMaterial::getS ( const util::Point dx,
const util::Point du 
) const
inlineoverridevirtual

Returns the bond strain.

Parameters
dxReference bond vector
duDifference of displacement
Returns
strain Bond strain \( S = \frac{du \cdot dx}{|dx|^2} \)

Implements material::Material.

Definition at line 791 of file material.h.

791 {
792 return ((dx + du).length() - dx.length()) / dx.length();
793 };
double length() const
Computes the Euclidean length of the vector.
Definition point.h:118

References util::Point::length().

Here is the call graph for this function:

◆ getSc()

double material::PmbMaterial::getSc ( const double &  r) const
inlineoverridevirtual

Returns critical bond strain.

Parameters
rReference length of bond
Returns
strain Critical strain

Implements material::Material.

Definition at line 801 of file material.h.

801{ return d_s0; };

References d_s0.

◆ isStateActive()

bool material::PmbMaterial::isStateActive ( ) const
inlineoverridevirtual

Returns true if state-based potential is active.

Returns
bool True/false

Implements material::Material.

Definition at line 721 of file material.h.

721{ return false; };

◆ print() [1/2]

void material::PmbMaterial::print ( ) const
inlineoverridevirtual

Prints the information about the object.

Reimplemented from material::Material.

Definition at line 910 of file material.h.

910{ print(0, 0); };
void print() const override
Prints the information about the object.
Definition material.h:910

References print().

Referenced by print().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ print() [2/2]

void material::PmbMaterial::print ( int  nt,
int  lvl 
) const
inlineoverridevirtual

Prints the information about the object.

Parameters
ntNumber of tabs to append before printing
lvlInformation level (higher means more information)

Reimplemented from material::Material.

Definition at line 905 of file material.h.

905 {
906 std::cout << printStr(nt, lvl);
907 };
std::string printStr(int nt, int lvl) const override
Returns the string containing printable information about the object.
Definition material.h:880

References printStr().

Here is the call graph for this function:

◆ printStr()

std::string material::PmbMaterial::printStr ( int  nt,
int  lvl 
) const
inlineoverridevirtual

Returns the string containing printable information about the object.

Parameters
ntNumber of tabs to append before printing
lvlInformation level (higher means more information)
Returns
string String containing printable information about the object

Reimplemented from material::Material.

Definition at line 880 of file material.h.

880 {
881
882 auto tabS = util::io::getTabS(nt);
883 std::ostringstream oss;
884 oss << tabS << "------- particle::PmbMaterial --------" << std::endl
885 << std::endl;
886 oss << tabS << "State active = " << 0 << std::endl;
887 oss << tabS << "Horizon = " << d_horizon << std::endl;
888 oss << tabS << "Influence fn address = " << influence_fn.get() << std::endl;
889 oss << tabS << "Influence fn info: " << std::endl;
890 oss << influence_fn->printStr(nt + 1, lvl);
891 oss << tabS << "Peridynamic parameters: " << std::endl;
892 oss << tabS << " c = " << d_c << std::endl;
893 oss << tabS << " s0 = " << d_s0 << std::endl;
894 oss << tabS << std::endl;
895
896 return oss.str();
897 };
std::string getTabS(int nt)
Returns tab spaces of given size.
Definition io.h:40

References d_c, d_horizon, d_s0, and util::io::getTabS().

Referenced by print().

Here is the call graph for this function:
Here is the caller graph for this function:

Field Documentation

◆ d_c

double material::PmbMaterial::d_c
private

Parameter C.

Definition at line 1013 of file material.h.

Referenced by computeMaterialProperties(), computeParameters(), getBondEF(), PmbMaterial(), and printStr().

◆ d_density

double material::PmbMaterial::d_density
private

Density.

Definition at line 1005 of file material.h.

Referenced by getDensity().

◆ d_horizon

double material::PmbMaterial::d_horizon
private

Horizon.

Definition at line 1002 of file material.h.

Referenced by computeMaterialProperties(), computeParameters(), getHorizon(), getInfFn(), and printStr().

◆ d_s0

double material::PmbMaterial::d_s0
private

Parameter \( \beta \).

Definition at line 1016 of file material.h.

Referenced by computeMaterialProperties(), computeParameters(), getBondEF(), getSc(), PmbMaterial(), and printStr().


The documentation for this class was generated from the following file: