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

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

#include <material.h>

Inheritance diagram for material::PdState:
Collaboration diagram for material::PdState:

Public Member Functions

 PdState (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_K
 Bulk modulus.
 
double d_G
 Shear modulus.
 
double d_s0
 Critical stretch.
 

Detailed Description

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

Definition at line 1330 of file material.h.

Constructor & Destructor Documentation

◆ PdState()

material::PdState::PdState ( 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 1339 of file material.h.

1340 : Material("PDState"), d_horizon(horizon), d_density(deck.d_density),
1341 d_K(0.), d_G(0.), d_s0(0.) {
1342
1343 // set global fields
1344 if (dimension != dim)
1345 dimension = dim;
1346
1347 if (is_plane_strain != deck.d_isPlaneStrain)
1349
1350 // create influence function
1351 if (deck.d_influenceFnType == 0) {
1352 if (influence_fn == nullptr)
1353 influence_fn = std::make_shared<material::ConstInfluenceFn>(
1354 deck.d_influenceFnParams, dim);
1355 }
1356 else if (deck.d_influenceFnType == 1) {
1357 if (influence_fn == nullptr)
1358 influence_fn = std::make_shared<material::LinearInfluenceFn>(
1359 deck.d_influenceFnParams, dim);
1360 }
1361 else if (deck.d_influenceFnType == 2) {
1362 if (influence_fn == nullptr)
1363 influence_fn = std::make_shared<material::GaussianInfluenceFn>(
1364 deck.d_influenceFnParams, dim);
1365 }
1366 else {
1367 std::cerr << "Error: Influence function type = "
1368 << deck.d_influenceFnType
1369 << " is invalid.\n";
1370 exit(1);
1371 }
1372
1373 // check if we need to compute the material parameters
1375 computeParameters(deck, dim);
1376 else {
1377 d_K = deck.d_bondPotentialParams[0];
1378 d_G = deck.d_bondPotentialParams[1];
1379 d_s0 = deck.d_bondPotentialParams[2];
1380 }
1381 };
Material(std::string name="")
Constructor.
Definition material.h:83
double d_G
Shear modulus.
Definition material.h:1718
void computeParameters(inp::MaterialDeck &deck, const size_t &dim)
Compute material model parameters.
Definition material.h:1579
double d_horizon
Horizon.
Definition material.h:1704
double d_K
Bulk modulus.
Definition material.h:1715
double d_density
Density.
Definition material.h:1707
double d_s0
Critical stretch.
Definition material.h:1721
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, inp::MaterialDeck::d_computeParamsFromElastic, d_G, inp::MaterialDeck::d_influenceFnParams, inp::MaterialDeck::d_influenceFnType, inp::MaterialDeck::d_isPlaneStrain, d_K, and d_s0.

Here is the call graph for this function:

Member Function Documentation

◆ computeMaterialProperties()

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

Computes elastic and fracture material properties and returns the data.

Parameters
dimDimension of the problem
Returns
Data Material data

Implements material::Material.

Definition at line 1504 of file material.h.

1504 {
1505
1506 auto data = inp::MatData();
1507
1508 // we already have G and K
1509 data.d_G = d_G;
1510 data.d_K = d_K;
1511
1512 // get Poisson ratio and Young's modulus
1513 data.d_nu = (3. * d_K - 2. * d_G) / (2. * (3. * d_K + d_G));
1514 data.d_E = data.toE(d_K, data.d_nu);
1515
1516 // get lame parameters
1517 data.d_lambda = data.toLambdaE(data.d_E, data.d_nu);
1518 data.d_mu = d_G;
1519
1520 // get Gc
1521 double d =
1522 (3. * d_G + std::pow(3. / 4., 4) * (d_K - 5. * d_G / 3.)) * d_horizon;
1523 data.d_Gc = d_s0 * d_s0 * d;
1524
1525 // KIc
1526 data.d_KIc = data.toKIc(
1527 data.d_Gc, data.d_nu, data.d_E);
1528
1529 return data;
1530 };
float d
initial distance between particle and wall
Structure for elastic properties and fracture properties.

References d_G, d_horizon, d_K, and d_s0.

◆ computeParameters()

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

Compute material model parameters.

Parameters
deckMaterialDeck
dimDimension of the domain

Definition at line 1579 of file material.h.

1579 {
1580 //
1581 // Need following elastic and fracture properties
1582 // 1. E or K
1583 // 2. Poisson ratio or shear modulus
1584 // 3. Gc or KIc
1585 //
1586 bool found_E = false;
1587 bool found_K = false;
1588 bool found_G = false;
1589 bool found_nu = false;
1590 size_t num_props = 0;
1591
1592 found_E = util::isGreater(deck.d_matData.d_E, 0.);
1593 if (found_E)
1594 num_props++;
1595
1596 found_K = util::isGreater(deck.d_matData.d_K, 0.);
1597 if (found_K)
1598 num_props++;
1599
1600 found_G = util::isGreater(deck.d_matData.d_G, 0.);
1601 if (found_G)
1602 num_props++;
1603
1604 found_nu = util::isGreater(deck.d_matData.d_nu, 0.);
1605 if (found_nu)
1606 num_props++;
1607
1608 if (num_props != 2) {
1609 std::ostringstream oss;
1610 oss << "Error: Require two different elastic properties for the "
1611 "PdState material. Pairs supported are (E, K), (E, G), "
1612 "(E, nu), (K, G).\n";
1613 oss << deck.printStr(0, 0) << "\n";
1614 std::cout << oss.str();
1615 exit(1);
1616 }
1617
1618 if (util::isLess(deck.d_matData.d_Gc, 0.) &&
1619 util::isLess(deck.d_matData.d_KIc, 0.)) {
1620 std::cerr << "Error: Require either critical energy release rate Gc or "
1621 "critical stress intensity factor KIc to compute the RNP "
1622 "bond-based peridynamic parameters.\n";
1623 exit(1);
1624 } else if (util::isGreater(deck.d_matData.d_Gc, 0.) &&
1625 util::isGreater(deck.d_matData.d_KIc, 0.)) {
1626 std::cout << "Warning: Both critical energy release rate Gc and critical "
1627 "stress intensity factor KIc are provided.\n";
1628 std::cout << "Warning: To compute the RNP bond-based peridynamic "
1629 "parameters, we only require one of those.\n";
1630 std::cout << "Warning: Selecting critical energy release rate Gc to "
1631 "compute parameters.\n";
1632 }
1633
1634 // compute nu if not provided
1635 if (!found_nu) {
1636 if (found_E and found_G)
1637 deck.d_matData.d_nu =
1638 (0.5 * deck.d_matData.d_E / deck.d_matData.d_G) - 1.;
1639
1640 if (found_E and found_K)
1641 deck.d_matData.d_nu =
1642 (3. * deck.d_matData.d_K - deck.d_matData.d_E) /
1643 (6. * deck.d_matData.d_K);
1644
1645 if (found_G and found_K)
1646 deck.d_matData.d_nu =
1647 (3. * deck.d_matData.d_K - 2. * deck.d_matData.d_G) /
1648 (2. * (3. * deck.d_matData.d_K + deck.d_matData.d_G));
1649 }
1650 found_nu = true;
1651
1652 // compute E if not provided
1653 if (!found_E) {
1654 if (found_K)
1655 deck.d_matData.d_E =
1656 deck.d_matData.toE(deck.d_matData.d_K, deck.d_matData.d_nu);
1657
1658 if (found_G)
1659 deck.d_matData.d_E =
1660 2. * deck.d_matData.d_G * (1. + deck.d_matData.d_nu);
1661 }
1662 found_E = true;
1663
1664 // compute K if not provided
1665 if (!found_K)
1666 deck.d_matData.d_K =
1667 deck.d_matData.toK(deck.d_matData.d_E, deck.d_matData.d_nu);
1668
1669 found_K = true;
1670
1671
1672 // compute G if not provided
1673 if (!found_G)
1674 deck.d_matData.d_G =
1675 deck.d_matData.toGE(deck.d_matData.d_E, deck.d_matData.d_nu);
1676
1677 found_G = true;
1678
1679 // compute Gc (if not provided) and KIc (if not provided)
1680 if (deck.d_matData.d_Gc > 0.)
1681 deck.d_matData.d_KIc = deck.d_matData.toKIc(
1682 deck.d_matData.d_Gc, deck.d_matData.d_nu, deck.d_matData.d_E);
1683
1684 if (deck.d_matData.d_KIc > 0. && deck.d_matData.d_Gc < 0.)
1685 deck.d_matData.d_Gc = deck.d_matData.toGc(
1686 deck.d_matData.d_KIc, deck.d_matData.d_nu, deck.d_matData.d_E);
1687
1688 // compute lame parameter
1689 deck.d_matData.d_lambda =
1690 deck.d_matData.toLambdaE(deck.d_matData.d_E, deck.d_matData.d_nu);
1691 deck.d_matData.d_mu = deck.d_matData.d_G;
1692
1693 // compute peridynamic parameters
1694 d_K = deck.d_matData.d_K;
1695 d_G = deck.d_matData.d_G;
1696
1697 double d =
1698 (3. * d_G + std::pow(3. / 4., 4) * (d_K - 5. * d_G / 3.)) * d_horizon;
1699 d_s0 = std::sqrt(deck.d_matData.d_Gc / d);
1700 };
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.
std::string printStr(int nt=0, int lvl=0) const
Returns the string containing printable information about the object.

References inp::MatData::d_E, inp::MatData::d_G, d_G, inp::MatData::d_Gc, d_horizon, inp::MatData::d_K, 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::MaterialDeck::printStr(), inp::MatData::toE(), inp::MatData::toGc(), inp::MatData::toGE(), inp::MatData::toK(), inp::MatData::toKIc(), and inp::MatData::toLambdaE().

Referenced by PdState().

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

◆ getBondEF() [1/2]

std::pair< double, double > material::PdState::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 1398 of file material.h.

1400 {
1401
1402 return {0., 0.};
1403 };

◆ getBondEF() [2/2]

std::pair< double, double > material::PdState::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 1416 of file material.h.

1417 {
1418
1419 if (fs)
1420 return {0., 0.};
1421
1422 double J = getInfFn(r);
1423 double change_length = s * r;
1424
1425 double alpha = 15. * d_G / mx;
1426 double factor = (3. * d_K / mx) - alpha / 3.;
1427
1428 return {0., J * (r * thetax * factor + change_length * alpha)};
1429 };
double getInfFn(const double &r) const override
Returns the value of influence function.
Definition material.h:1473

References d_G, d_K, and getInfFn().

Here is the call graph for this function:

◆ getBondForceDirection()

util::Point material::PdState::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 1438 of file material.h.

1439 {
1440 return (dx + du) / (dx + du).length();
1441 };

◆ getDensity()

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

Returns the density of the material.

Returns
density Density of the material

Implements material::Material.

Definition at line 1465 of file material.h.

1465{ return d_density; };

References d_density.

◆ getHorizon()

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

Returns horizon.

Returns
horizon Horizon

Implements material::Material.

Definition at line 1495 of file material.h.

1495{ return d_horizon; };

References d_horizon.

◆ getInfFn()

double material::PdState::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 1473 of file material.h.

1473 {
1474 return getGlobalInfFn(r / d_horizon);
1475 };
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::PdState::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 1486 of file material.h.

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

◆ getS()

double material::PdState::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 1449 of file material.h.

1449 {
1450 return ((dx + du).length() - dx.length()) / dx.length();
1451 };
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::PdState::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 1459 of file material.h.

1459{ return d_s0; };

References d_s0.

◆ isStateActive()

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

Returns true if state-based potential is active.

Returns
bool True/false

Implements material::Material.

Definition at line 1387 of file material.h.

1387{ return true; };

◆ print() [1/2]

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

Prints the information about the object.

Reimplemented from material::Material.

Definition at line 1570 of file material.h.

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

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::PdState::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 1565 of file material.h.

1565 {
1566 std::cout << printStr(nt, lvl);
1567 };
std::string printStr(int nt, int lvl) const override
Returns the string containing printable information about the object.
Definition material.h:1539

References printStr().

Here is the call graph for this function:

◆ printStr()

std::string material::PdState::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 1539 of file material.h.

1539 {
1540
1541 auto tabS = util::io::getTabS(nt);
1542 std::ostringstream oss;
1543 oss << tabS << "------- particle::PdState --------" << std::endl
1544 << std::endl;
1545 oss << tabS << "State active = " << 1 << std::endl;
1546 oss << tabS << "Horizon = " << d_horizon << std::endl;
1547 oss << tabS << "Influence fn address = " << influence_fn.get() << std::endl;
1548 oss << tabS << "Influence fn info: " << std::endl;
1549 oss << influence_fn->printStr(nt + 1, lvl);
1550 oss << tabS << "Peridynamic parameters: " << std::endl;
1551 oss << tabS << " K = " << d_K << std::endl;
1552 oss << tabS << " G = " << d_G << std::endl;
1553 oss << tabS << " s0 = " << d_s0 << std::endl;
1554 oss << tabS << std::endl;
1555
1556 return oss.str();
1557 };
std::string getTabS(int nt)
Returns tab spaces of given size.
Definition io.h:40

References d_G, d_horizon, d_K, 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_density

double material::PdState::d_density
private

Density.

Definition at line 1707 of file material.h.

Referenced by getDensity().

◆ d_G

double material::PdState::d_G
private

Shear modulus.

Definition at line 1718 of file material.h.

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

◆ d_horizon

double material::PdState::d_horizon
private

Horizon.

Definition at line 1704 of file material.h.

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

◆ d_K

double material::PdState::d_K
private

Bulk modulus.

Definition at line 1715 of file material.h.

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

◆ d_s0

double material::PdState::d_s0
private

Critical stretch.

Definition at line 1721 of file material.h.

Referenced by computeMaterialProperties(), computeParameters(), getSc(), PdState(), and printStr().


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