NLMech  0.1.0
materialDeck.h
1 // Copyright (c) 2019 Prashant K. Jha
3 // Copyright (c) 2019 Patrick Diehl
4 //
5 // Distributed under the Boost Software License, Version 1.0. (See accompanying
6 // file LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
8 
9 #ifndef INP_MATERIALDECK_H
10 #define INP_MATERIALDECK_H
11 
12 #include "util/utilIO.h"
13 #include <cmath>
14 #include <string>
15 #include <vector>
16 
17 namespace inp {
18 
20 struct MatData {
21 
28  double d_E;
29 
31  double d_G;
32 
34  double d_K;
35 
37  double d_nu;
38 
40  double d_lambda;
41 
43  double d_mu;
44 
53  double d_KIc;
54 
56  double d_Gc;
57 
64  : d_E(-1.), d_G(-1.), d_K(-1.), d_nu(-1.), d_lambda(-1.), d_mu(-1.),
65  d_KIc(-1.), d_Gc(-1.){};
66 
74  std::string printStr(int nt = 0, int lvl = 0) const {
75 
76  auto tabS = util::io::getTabS(nt);
77  std::ostringstream oss;
78  oss << tabS << "------- MatData --------" << std::endl << std::endl;
79  oss << tabS << "Young's modulus = " << d_E << std::endl;
80  oss << tabS << "Shear modulus = " << d_G << std::endl;
81  oss << tabS << "Bulk modulus = " << d_K << std::endl;
82  oss << tabS << "Poisson ratio = " << d_nu << std::endl;
83  oss << tabS << "Lame parameter Lambda = " << d_lambda << std::endl;
84  oss << tabS << "Lame parameter Mu = " << d_mu << std::endl;
85  oss << tabS << "Critical stress intensity factor = " << d_KIc << std::endl;
86  oss << tabS << "Critical energy release rate = " << d_Gc << std::endl;
87  oss << tabS << std::endl;
88 
89  return oss.str();
90  }
91 
99  void print(int nt = 0, int lvl = 0) const { std::cout << printStr(nt, lvl); }
100 
112  double toNu(double lambda, double mu) { return lambda * 0.5 / (lambda + mu); }
113 
120  double toNuClassical(double K, double G) { return (3*K - 2*G) / (2*(3*K+G)); }
121 
128  double toE(double K, double nu) { return K * (3. * (1. - 2. * nu)); }
129 
136  double toK(double E, double nu) { return E / (3. * (1. - 2. * nu)); }
137 
145  double toLambdaE(double E, double nu) {
146  return E * nu / ((1. + nu) * (1. - 2. * nu));
147  }
148 
156  double toLambdaK(double K, double nu) { return 3. * K * nu / (1. + nu); }
157 
164  double toGE(double E, double nu) { return E / (2. * (1. + nu)); }
165 
172  double toGK(double K, double nu) {
173  return 3. * K * (1. - 2. * nu) / (2. * (1. + nu));
174  }
175 
183  double toELambda(double lambda, double nu) {
184  return lambda * (1. + nu) * (1. - 2. * nu) / nu;
185  }
186 
200  double toGc(double KIc, double nu, double E) { return KIc * KIc / E; }
201 
215  double toKIc(double Gc, double nu, double E) { return std::sqrt(Gc * E); }
216 
218 };
219 
226 struct MaterialDeck {
227 
233 
235  std::string d_materialType;
236 
239 
242 
245 
247  std::vector<double> d_bondPotentialParams;
248 
250  std::vector<double> d_statePotentialParams;
251 
253  std::vector<double> d_influenceFnParams;
254 
261 
264 
267 
270 
273 
275  double d_density;
276 
279 
288  d_density(1.), d_applyContact(false){};
289 
297  std::string printStr(int nt = 0, int lvl = 0) const {
298 
299  auto tabS = util::io::getTabS(nt);
300  std::ostringstream oss;
301  oss << tabS << "------- MaterialDeck --------" << std::endl << std::endl;
302  oss << tabS << "Is plain strain = " << d_isPlaneStrain << std::endl;
303  oss << tabS << "Material type = " << d_materialType << std::endl;
304  oss << tabS << "Bond potential type = " << d_bondPotentialType << std::endl;
305  oss << tabS << "Bond potential params = ["
306  << util::io::printStr<double>(d_bondPotentialParams, 0) << "]"
307  << std::endl;
308  oss << tabS << "State potential type = " << d_statePotentialType
309  << std::endl;
310  oss << tabS << "State potential params = ["
311  << util::io::printStr<double>(d_statePotentialParams, 0) << "]"
312  << std::endl;
313  oss << tabS << "Influence function type = " << d_influenceFnType
314  << std::endl;
315  oss << tabS << "Influence function params = ["
316  << util::io::printStr<double>(d_influenceFnParams, 0) << "]"
317  << std::endl;
318  oss << tabS
319  << "Irreversible bond breaking enabled = " << d_irreversibleBondBreak
320  << std::endl;
321  oss << tabS << "State contribution from broken bond enabled = "
322  << d_stateContributionFromBrokenBond << std::endl;
323  oss << tabS << "Check Sc factor = " << d_checkScFactor << std::endl;
324  oss << tabS << "Compute parameters from elastic properties = "
325  << d_computeParamsFromElastic << std::endl;
326  oss << d_matData.printStr(nt + 1, lvl);
327  oss << tabS << "Density = " << d_density << std::endl;
328  oss << tabS << "Apply contact on broken bonds = " << d_applyContact << std::endl;
329  oss << tabS << std::endl;
330 
331  return oss.str();
332  }
333 
340  void print(int nt = 0, int lvl = 0) const { std::cout << printStr(nt, lvl); }
341 };
342 
345 } // namespace inp
346 
347 #endif // INP_MATERIALDECK_H
Collection of methods and database related to input.
Definition: main.cpp:21
std::string getTabS(int nt)
Generate a string contaning nt tabs.
Definition: utilIO.h:26
Structure for elastic properties and fracture properties.
Definition: materialDeck.h:20
double toGc(double KIc, double nu, double E)
Compute critical energy release rate Gc from critical stress-intensity factor KIc,...
Definition: materialDeck.h:200
double toLambdaE(double E, double nu)
Compute Lame first parameter lambda from Young's modulus E and Poisson's ratio nu.
Definition: materialDeck.h:145
void print(int nt=0, int lvl=0) const
Prints the information.
Definition: materialDeck.h:99
double d_mu
Lame second parameter.
Definition: materialDeck.h:43
double toGK(double K, double nu)
Compute shear modulus from Bulk modulus K and Poisson's ratio nu.
Definition: materialDeck.h:172
double toKIc(double Gc, double nu, double E)
Compute critical stress-intensity factor KIc from critical energy release rate Gc,...
Definition: materialDeck.h:215
double d_KIc
Critical stress intensity factor.
Definition: materialDeck.h:53
double toK(double E, double nu)
Compute Bulk modulus K from Young's modulus K and Poisson's ratio nu.
Definition: materialDeck.h:136
double d_K
Bulk modulus.
Definition: materialDeck.h:34
double d_lambda
Lame first parameter.
Definition: materialDeck.h:40
double toNu(double lambda, double mu)
Compute Poisson's ratio from Lame parameters.
Definition: materialDeck.h:112
double toGE(double E, double nu)
Compute shear modulus from Young's modulus E and Poisson's ratio nu.
Definition: materialDeck.h:164
double d_nu
Poisson's ratio.
Definition: materialDeck.h:37
double d_G
Shear modulus or Lame second parameter.
Definition: materialDeck.h:31
double toELambda(double lambda, double nu)
Compute Young's modulus E from Lame first parameter lambda and Poisson's ratio nu.
Definition: materialDeck.h:183
double toLambdaK(double K, double nu)
Compute Lame first parameter lambda from Bulk modulus K and Poisson's ratio nu.
Definition: materialDeck.h:156
double d_E
Young's elastic modulus.
Definition: materialDeck.h:28
MatData()
Constructor.
Definition: materialDeck.h:63
std::string printStr(int nt=0, int lvl=0) const
Prints the information.
Definition: materialDeck.h:74
double d_Gc
Critical energy release rate.
Definition: materialDeck.h:56
double toNuClassical(double K, double G)
Compute Poisson's ratio from Bulk modulus and Shear modulus.
Definition: materialDeck.h:120
double toE(double K, double nu)
Compute Young's modulus E from Bulk modulus K and Poisson's ratio nu.
Definition: materialDeck.h:128
Structure to read and store material related data.
Definition: materialDeck.h:226
bool d_isPlaneStrain
Indicates if the 2-d simulation is of plane-strain type (thick material) or plane-stress type (thin m...
Definition: materialDeck.h:232
double d_density
Density of material.
Definition: materialDeck.h:275
size_t d_bondPotentialType
Type of pairwise (bond-based) potential.
Definition: materialDeck.h:238
size_t d_influenceFnType
Type of influence function.
Definition: materialDeck.h:244
std::vector< double > d_bondPotentialParams
List of parameters for pairwise potential.
Definition: materialDeck.h:247
MaterialDeck()
Constructor.
Definition: materialDeck.h:283
std::string d_materialType
Material type.
Definition: materialDeck.h:235
bool d_stateContributionFromBrokenBond
Flag for contribution to hydrostatic force from the broken bond.
Definition: materialDeck.h:263
bool d_computeParamsFromElastic
Compute Peridynamic material properties from elastic properties.
Definition: materialDeck.h:269
inp::MatData d_matData
List of elastic and fracture properties.
Definition: materialDeck.h:272
std::vector< double > d_influenceFnParams
List of parameters for influence function.
Definition: materialDeck.h:253
size_t d_statePotentialType
Type of hydrostatic (state-based) potential.
Definition: materialDeck.h:241
std::vector< double > d_statePotentialParams
List of parameters for hydrostatic potential.
Definition: materialDeck.h:250
bool d_irreversibleBondBreak
Flag for irreversible breaking of bonds.
Definition: materialDeck.h:260
double d_checkScFactor
Factor to check if bond is broken.
Definition: materialDeck.h:266
std::string printStr(int nt=0, int lvl=0) const
Prints the information.
Definition: materialDeck.h:297
void print(int nt=0, int lvl=0) const
Prints the information.
Definition: materialDeck.h:340
bool d_applyContact
Enable non-penetration of broken bonds.
Definition: materialDeck.h:278