PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
Loading...
Searching...
No Matches
materialDeck.h
Go to the documentation of this file.
1/*
2 * -------------------------------------------
3 * Copyright (c) 2021 - 2024 Prashant K. Jha
4 * -------------------------------------------
5 * PeriDEM https://github.com/prashjha/PeriDEM
6 *
7 * Distributed under the Boost Software License, Version 1.0. (See accompanying
8 * file LICENSE)
9 */
10
11#ifndef INP_MATERIALDECK_H
12#define INP_MATERIALDECK_H
13
14#include "util/io.h"
15#include <cmath>
16#include <string>
17#include <vector>
18
19namespace inp {
20
22struct MatData {
23
30 double d_E;
31
33 double d_G;
34
36 double d_K;
37
39 double d_nu;
40
42 double d_lambda;
43
45 double d_mu;
46
55 double d_KIc;
56
58 double d_Gc;
59
66 : d_E(-1.), d_G(-1.), d_K(-1.), d_nu(-1.), d_lambda(-1.), d_mu(-1.),
67 d_KIc(-1.), d_Gc(-1.){};
68
74 MatData(const MatData &md)
75 : d_E(md.d_E), d_G(md.d_G), d_K(md.d_K), d_nu(md.d_nu),
76 d_lambda(md.d_lambda), d_mu(md.d_mu), d_KIc(md.d_KIc), d_Gc(md.d_Gc){};
77
85 std::string printStr(int nt = 0, int lvl = 0) const {
86
87 auto tabS = util::io::getTabS(nt);
88 std::ostringstream oss;
89 oss << tabS << "------- MatData --------" << std::endl << std::endl;
90 oss << tabS << "Young's modulus = " << d_E << std::endl;
91 oss << tabS << "Shear modulus = " << d_G << std::endl;
92 oss << tabS << "Bulk modulus = " << d_K << std::endl;
93 oss << tabS << "Poisson ratio = " << d_nu << std::endl;
94 oss << tabS << "Lame parameter Lambda = " << d_lambda << std::endl;
95 oss << tabS << "Lame parameter Mu = " << d_mu << std::endl;
96 oss << tabS << "Critical stress intensity factor = " << d_KIc << std::endl;
97 oss << tabS << "Critical energy release rate = " << d_Gc << std::endl;
98 oss << tabS << std::endl;
99
100 return oss.str();
101 }
102
109 void print(int nt = 0, int lvl = 0) const { std::cout << printStr(nt, lvl); }
110
122 double toNu(double lambda, double mu) { return lambda * 0.5 / (lambda + mu); }
123
130 double toNuEG(double E, double G) { return E * 0.5 / G - 1.; }
131
138 double toE(double K, double nu) { return K * (3. * (1. - 2. * nu)); }
139
146 double toK(double E, double nu) { return E / (3. * (1. - 2. * nu)); }
147
155 double toLambdaE(double E, double nu) {
156 return E * nu / ((1. + nu) * (1. - 2. * nu));
157 }
158
166 double toLambdaK(double K, double nu) { return 3. * K * nu / (1. + nu); }
167
174 double toGE(double E, double nu) { return E / (2. * (1. + nu)); }
175
182 double toGK(double K, double nu) {
183 return 3. * K * (1. - 2. * nu) / (2. * (1. + nu));
184 }
185
193 double toELambda(double lambda, double nu) {
194 return lambda * (1. + nu) * (1. - 2. * nu) / nu;
195 }
196
210 double toGc(double KIc, double nu, double E) { return KIc * KIc / E; }
211
225 double toKIc(double Gc, double nu, double E) { return std::sqrt(Gc * E); }
226
228};
229
237
243
245 std::string d_materialType;
246
249
252
255
257 std::vector<double> d_bondPotentialParams;
258
260 std::vector<double> d_statePotentialParams;
261
263 std::vector<double> d_influenceFnParams;
264
271
274
277
280
283
285 double d_density;
286
288 double d_horizon;
289
292
302
323
331 std::string printStr(int nt = 0, int lvl = 0) const {
332
333 auto tabS = util::io::getTabS(nt);
334 std::ostringstream oss;
335 oss << tabS << "------- MaterialDeck --------" << std::endl << std::endl;
336 oss << tabS << "Is plain strain = " << d_isPlaneStrain << std::endl;
337 oss << tabS << "Material type = " << d_materialType << std::endl;
338 oss << tabS << "Bond potential type = " << d_bondPotentialType << std::endl;
339 oss << tabS << "Bond potential params = ["
340 << util::io::printStr<double>(d_bondPotentialParams, 0) << "]"
341 << std::endl;
342 oss << tabS << "State potential type = " << d_statePotentialType
343 << std::endl;
344 oss << tabS << "State potential params = ["
345 << util::io::printStr<double>(d_statePotentialParams, 0) << "]"
346 << std::endl;
347 oss << tabS << "Influence function type = " << d_influenceFnType
348 << std::endl;
349 oss << tabS << "Influence function params = ["
350 << util::io::printStr<double>(d_influenceFnParams, 0) << "]"
351 << std::endl;
352 oss << tabS
353 << "Irreversible bond breaking enabled = " << d_irreversibleBondBreak
354 << std::endl;
355 oss << tabS << "State contribution from broken bond enabled = "
356 << d_stateContributionFromBrokenBond << std::endl;
357 oss << tabS << "Check Sc factor = " << d_checkScFactor << std::endl;
358 oss << tabS << "Compute parameters from elastic properties = "
359 << d_computeParamsFromElastic << std::endl;
360 oss << d_matData.printStr(nt + 1, lvl);
361 oss << tabS << "Density = " << d_density << std::endl;
362 oss << tabS << "Horizon = " << d_horizon << std::endl;
363 oss << tabS << "Horizon to mesh ratio = " << d_horizonMeshRatio
364 << std::endl;
365 oss << tabS << std::endl;
366
367 return oss.str();
368 }
369
376 void print(int nt = 0, int lvl = 0) const { std::cout << printStr(nt, lvl); }
377};
378
381} // namespace inp
382
383#endif // INP_MATERIALDECK_H
Collection of methods and database related to input.
Definition mesh.h:20
std::string getTabS(int nt)
Returns tab spaces of given size.
Definition io.h:40
Structure for elastic properties and fracture properties.
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.
void print(int nt=0, int lvl=0) const
Prints the information about the object.
double d_mu
Lame second parameter.
double toGK(double K, double nu)
Compute shear modulus from Bulk modulus K and Poisson's ratio nu.
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.
MatData(const MatData &md)
Copy constructor.
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 toNu(double lambda, double mu)
Compute Poisson's ratio from Lame parameters.
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 toELambda(double lambda, double nu)
Compute Young's modulus E from Lame first parameter lambda and Poisson's ratio nu.
double toLambdaK(double K, double nu)
Compute Lame first parameter lambda from Bulk modulus K and Poisson's ratio nu.
double toNuEG(double E, double G)
Compute Poisson's ratio from Young's and Shear modulus.
double d_E
Young's elastic modulus.
MatData()
Constructor.
std::string printStr(int nt=0, int lvl=0) const
Returns the string containing printable information about the object.
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.
Structure to read and store material related data.
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_bondPotentialType
Type of pairwise (bond-based) potential.
double d_horizonMeshRatio
Horizon to mesh ratio.
size_t d_influenceFnType
Type of influence function.
std::vector< double > d_bondPotentialParams
List of parameters for pairwise potential.
MaterialDeck()
Constructor.
std::string d_materialType
Material type.
double d_horizon
Horizon for peridynamic interaction.
bool d_stateContributionFromBrokenBond
Flag for contribution to hydrostatic force from the broken bond.
bool d_computeParamsFromElastic
Compute Peridynamic material properties from elastic properties.
inp::MatData d_matData
List of elastic and fracture properties.
std::vector< double > d_influenceFnParams
List of parameters for influence function.
size_t d_statePotentialType
Type of hydrostatic (state-based) potential.
std::vector< double > d_statePotentialParams
List of parameters for hydrostatic potential.
bool d_irreversibleBondBreak
Flag for irreversible breaking of bonds.
double d_checkScFactor
Factor to check if bond is broken.
std::string printStr(int nt=0, int lvl=0) const
Returns the string containing printable information about the object.
void print(int nt=0, int lvl=0) const
Prints the information about the object.
MaterialDeck(const MaterialDeck &md)
Copy constructor.