NLMech  0.1.0
ElasticState.h
1 // Copyright (c) 2019 Patrick Diehl
2 //
3 // Distributed under the GNU GENERAL PUBLIC LICENSE, Version 3.0.
4 // (See accompanying file LICENSE.txt)
5 
6 #ifndef MATERIAL_PD_ELASTIC_STATEMATERIAL_H
7 #define MATERIAL_PD_ELASTIC_STATEMATERIAL_H
8 
9 #include "baseMaterial.h"
10 #include "util/point.h"
11 #include "util/matrixBlaze.h"
12 #include "geometry/neighbor.h"
13 #include "fe/mesh.h"
14 #include "geometry/volumeCorrection.h"
15 
16 #include <vector>
17 
18 #include <hpx/include/parallel_algorithm.hpp>
19 
20 // forward declaration
21 namespace inp {
22 struct MaterialDeck;
23 struct OutputDeck;
24 }
25 
26 // forward declaration
27 namespace data {
28 class DataManager;
29 }
30 
31 // forward declaration of class
32 namespace fe {
33 class Mesh;
34 } // namespace fe
35 
36 namespace geometry {
37 class VolumeCorrection;
38 } // namespace geometry
39 
40 namespace util {
41 class StateBasedHelperFunctions;
42 } // namespace material
43 
44 namespace material {
45 
46 namespace pd {
47 
56 class ElasticState : public BaseMaterial {
57 
58  public:
65 
73  std::pair<util::Point3, double> getBondEF(size_t i, size_t j);
74 
82  double getSc(size_t i, size_t j);
83 
89  util::Matrix33 getStress(size_t i);
90 
96  util::Matrix33 getStrain(size_t i);
97 
101  double getFactor2D();
102 
103  void update();
104 
111  const util::Point3 &du) const override {
112  return dx / dx.length();
113  }
114 
121  double getInfFn(const double &r) const {return 1.; };
122 
123  private:
133  void computeParameters(inp::MaterialDeck *deck, size_t dim);
134 
147 
155  util::Point3 X_vector_state(size_t i, size_t j);
156 
164  util::Point3 Y_vector_state(size_t i, size_t j);
165 
172  util::Matrix33 K_shape_tensor(size_t i);
173 
183  util::Matrix33 K_modulus_tensor(size_t i, size_t j, size_t k, size_t m);
184 
196  double dirac_delta(util::Point3 x, size_t i, size_t j, size_t m);
197 
203  double d_factor2D;
204 
206  size_t dim;
207 
208  double horizon;
209 
212 
213 
221 
223  //const std::vector<util::Point3> *d_nodes;
224 
226  //const std::vector<double> *d_weightedVolume;
227 
229  //const std::vector<std::vector<double>> *d_volumeCorrection;
230 
232  //const std::vector<double> * d_volumes;
233 
235  //geometry::Neighbor *d_neighbors;
236 
238 
242 };
243 
244 } // namespace pd
245 
246 } // namespace material
247 
248 #endif // MATERIAL_PD_ELASTIC_STATEMATERIAL_H
Data manager to collect the global simulation data.
Definition: DataManager.h:58
A base class providing methods to compute energy density and force.
Definition: baseMaterial.h:26
A class implementing the state-based elastic material model.
Definition: ElasticState.h:56
util::Matrix33 K_modulus_tensor(size_t i, size_t j, size_t k, size_t m)
Computes the K modulus tensor.
Definition: ElasticState.cpp:359
double getInfFn(const double &r) const
Returns the value of influence function.
Definition: ElasticState.h:121
util::Point3 X_vector_state(size_t i, size_t j)
Computes the x vector state between node i and node j.
Definition: ElasticState.cpp:282
util::Matrix33 getStress(size_t i)
Computes the stress tensor for one node.
Definition: ElasticState.cpp:432
util::Matrix33 getStrain(size_t i)
Computes the strain vector for one node.
Definition: ElasticState.cpp:336
void computeParameters(inp::MaterialDeck *deck, size_t dim)
Computes elastic state-based material parameters from elastic constants.
Definition: ElasticState.cpp:47
const inp::MaterialDeck * d_deck
Pointer to the material deck.
Definition: ElasticState.h:220
ElasticState(inp::MaterialDeck *deck, data::DataManager *dataManager)
Constructor.
Definition: ElasticState.cpp:17
util::Matrix33 deformation_gradient(size_t i)
Computes the deformation gradient for node i.
Definition: ElasticState.cpp:304
util::Point3 getBondForceDirection(const util::Point3 &dx, const util::Point3 &du) const override
Get direction of bond force.
Definition: ElasticState.h:110
util::Matrix33 K_shape_tensor(size_t i)
Computes the K shape tensor for node i.
Definition: ElasticState.cpp:287
double getSc(size_t i, size_t j)
Returns critical bond strain between node i and node j.
Definition: ElasticState.cpp:273
double getFactor2D()
return the factor for two-dimensional problems
Definition: ElasticState.cpp:473
data::DataManager * d_dataManager_p
Pointer to the nodes.
Definition: ElasticState.h:237
size_t dim
Dimension of the problem.
Definition: ElasticState.h:206
util::Point3 Y_vector_state(size_t i, size_t j)
Computes the y vector state between node i and node j.
Definition: ElasticState.cpp:275
void update()
Let the material class in the quasi-static case know that there is a new loading step.
Definition: ElasticState.cpp:475
bool strainEnergy
Compute strain energy.
Definition: ElasticState.h:211
std::pair< util::Point3, double > getBondEF(size_t i, size_t j)
Returns energy and force state between node i and node j.
Definition: ElasticState.cpp:160
double dirac_delta(util::Point3 x, size_t i, size_t j, size_t m)
Provide the image of x under the Dirac Delta Function.
Definition: ElasticState.cpp:347
double d_factor2D
Correction factor for using plain stress.
Definition: ElasticState.h:203
Data mamanger to share the global simulation data between the classes.
Definition: DataManager.h:55
Collection of methods and data related to finite element and mesh.
Definition: baseElem.h:15
Collection of methods and data related to geometry.
Definition: DataManager.h:23
Collection of methods and database related to input.
Definition: main.cpp:21
Collection of methods useful in simulation.
Definition: DataManager.h:50
blaze::StaticMatrix< double, 3UL, 3UL > Matrix33
Blaze: Definition of 3 x 3 matrix.
Definition: matrixBlaze.h:23
Structure to read and store material related data.
Definition: materialDeck.h:226
A structure to represent 3d vectors.
Definition: point.h:29
double length() const
Computes the Euclidean length of the vector.
Definition: point.h:119