NLMech  0.1.0
DataManager.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 SRC_DATA_DATAMANAGER_H_
10 #define SRC_DATA_DATAMANAGER_H_
11 
12 #include "util/point.h"
13 
14 #include <vector>
15 
16 // forward declaration of class
17 namespace fe {
18 class MassMatrix;
19 class Mesh;
20 // class Quadrature;
21 }// namespace fe
22 
23 namespace geometry {
24 class Fracture;
25 class InteriorFlags;
26 class Neighbor;
27 class VolumeCorrection;
28 } // namespace geometry
29 
30 namespace inp {
31 struct ModelDeck;
32 struct RestartDeck;
33 struct OutputDeck;
34 class Input;
35 class Policy;
36 } // namespace inp
37 
38 namespace loading {
39 class InitialCondition;
40 class ULoading;
41 class FLoading;
42 } // namespace loading
43 
44 namespace material {
45 namespace pd {
46 class BaseMaterial;
47 }
48 } // namespace material
49 
50 namespace util {
52 } // namespace util
53 
55 namespace data {
56 
58 class DataManager {
59 
60 public:
61 
63  DataManager();
64 
74  void setModelDeckP(inp::ModelDeck *pointer);
75 
81 
85  void setOutputDeckP(inp::OutputDeck *pointer);
86 
92 
104  void setBodyForceP(std::vector<util::Point3> *pointer);
105 
109  std::vector<util::Point3>* getBodyForceP();
110 
114  void setForceP(std::vector<util::Point3> *pointer);
115 
119  std::vector<util::Point3>* getForceP();
120 
124  void setReactionForceP(std::vector<util::Point3> *pointer);
125 
129  std::vector<util::Point3>* getReactionForceP();
130 
134  void setVelocityP(std::vector<util::Point3> *pointer);
135 
139  std::vector<util::Point3>* getVelocityP();
140 
144  void setDisplacementP(std::vector<util::Point3> *pointer);
145 
146 
150  std::vector<util::Point3>* getDisplacementP();
151 
155  void setMeshP(fe::Mesh *d_mesh_p);
156 
157 
161  fe::Mesh* getMeshP();
162 
167 
172 
177 
178 
183 
187  void setFractureP(geometry::Fracture *pointer);
188 
189 
194 
199 
200 
205 
210 
215 
220 
221 
226 
230  void setForceLoadingP(loading::FLoading * pointer);
231 
236 
240  void setExtensionP(std::vector<std::vector<double>>* pointer);
241 
242 
246  std::vector<std::vector<double>>* getExtensionP();
247 
251  void setStressTensorP(std::vector<util::Matrix33>* pointer);
252 
253 
257  std::vector<util::Matrix33>* getStressTensorP();
258 
262  void setStrainTensorP(std::vector<util::Matrix33>* pointer);
263 
267  std::vector<util::Matrix33>* getStrainTensorP();
268 
272  void setDilatationP(std::vector<double> * pointer);
273 
277  std::vector<double> * getDilatationP();
278 
282  void setWeightedVolP(std::vector<double> * pointer);
283 
287  std::vector<double> * getWeightedVolP();
288 
292  void setTotalReactionForceP(std::vector<double>* pointer);
293 
294 
298  std::vector<double>* getTotalReactionForceP();
299 
303  void setStrainEnergyP(std::vector<float>* pointer);
304 
308  std::vector<float>* getStrainEnergyP();
309 
313  void setWorkDoneP(std::vector<float>* pointer);
314 
318  std::vector<float>* getWorkDoneP();
319 
323  void setPhiP(std::vector<float>* pointer);
324 
325 
326 
330  std::vector<float>* getPhiP();
331 
335  void setDamageFunctionP(std::vector<float>* pointer);
336 
337 
341  std::vector<float>* getDamageFunctionP();
342 
346  void setBBFractureEnergyP(std::vector<float>*pointer);
347 
348 
352  std::vector<float>* getBBFractureEnergyP();
353 
357  void setFractureEnergyP(std::vector<float>*pointer);
358 
362  std::vector<float>* getFractureEnergyP();
363 
367  void setKineticEnergyP(std::vector<float>* pointer);
368 
372  std::vector<float>* getKineticEnergyP();
373 
385  std::string printStr(int nt = 0, int lvl = 0) const;
386 
394  void print(int nt = 0, int lvl = 0) const { std::cout << printStr(nt, lvl); };
395 
396 private:
397 
406 
409 
419  std::vector<util::Point3> *d_b_p = nullptr;
420 
422  std::vector<util::Point3> *d_u_p = nullptr;
423 
425  std::vector<util::Point3> *d_v_p = nullptr;
426 
428  std::vector<util::Point3> *d_f_p = nullptr;
429 
431  std::vector<util::Point3> *d_reaction_force_p = nullptr;
432 
434  std::vector<std::vector<double>>* d_extension_p = nullptr;
435 
437  std::vector<double> *d_dilatation_p = nullptr;
438 
440  std::vector<double> *d_weighted_vol_p = nullptr;
441 
447  std::vector<double> d_thetaX;
448 
453  std::vector<double> d_mX;
454 
456  std::vector<float>* d_e_p = nullptr;
457 
459  std::vector<double>* d_total_reaction_force_p = nullptr;
460 
462  std::vector<float>* d_w_p = nullptr;
463 
465  std::vector<float>* d_phi_p = nullptr;
466 
468  std::vector<float>* d_Z_p = nullptr;
469 
471  std::vector<float>* d_eFB_p = nullptr;
472 
474  std::vector<float>* d_eF_p = nullptr;
475 
477  std::vector<float>* d_ke_p = nullptr;
478 
480  std::vector<util::Matrix33>* d_strain_p = nullptr;
481 
483  std::vector<util::Matrix33>* d_stress_p = nullptr;
484 
490  fe::Mesh *d_mesh_p = nullptr;
491 
494 
497 
500 
503 
506 
507 
510 
513 
516 };
517 
518 }
519 
520 #endif /* SRC_DATA_DATAMANAGER_H_ */
Data manager to collect the global simulation data.
Definition: DataManager.h:58
std::vector< float > * d_eF_p
Pointer to fracture energy of the nodes.
Definition: DataManager.h:474
loading::ULoading * getDisplacementLoadingP()
Definition: DataManager.cpp:122
std::vector< double > * getWeightedVolP()
Definition: DataManager.cpp:177
std::vector< util::Point3 > * d_v_p
Pointer to the velocity vector.
Definition: DataManager.h:425
void setFractureEnergyP(std::vector< float > *pointer)
Definition: DataManager.cpp:223
void setStateBasedHelperFunctionsP(util::StateBasedHelperFunctions *pointer)
Definition: DataManager.cpp:108
std::vector< double > * getTotalReactionForceP()
Definition: DataManager.cpp:193
geometry::InteriorFlags * getInteriorFlagsP()
Definition: DataManager.cpp:100
std::vector< util::Matrix33 > * d_stress_p
Pointer to the stress tensor vector.
Definition: DataManager.h:483
void setStrainTensorP(std::vector< util::Matrix33 > *pointer)
Definition: DataManager.cpp:157
std::vector< float > * d_ke_p
Pointer to the kinetic energy of the nodes.
Definition: DataManager.h:477
std::vector< util::Matrix33 > * d_strain_p
Pointer to the strain tensor vector.
Definition: DataManager.h:480
void setVolumeCorrectionP(geometry::VolumeCorrection *pointer)
Definition: DataManager.cpp:85
std::vector< std::vector< double > > * d_extension_p
Extension for the neighborhood of each node.
Definition: DataManager.h:434
std::vector< float > * d_eFB_p
Pointer to bond-based fracture energy of the nodes.
Definition: DataManager.h:471
void setExtensionP(std::vector< std::vector< double >> *pointer)
Definition: DataManager.cpp:134
void setNeighborP(geometry::Neighbor *d_neighbor_p)
Definition: DataManager.cpp:79
std::vector< std::vector< double > > * getExtensionP()
Definition: DataManager.cpp:139
std::string printStr(int nt=0, int lvl=0) const
Returns the string containing information about the instance of the object.
Definition: DataManager.cpp:235
std::vector< float > * getWorkDoneP()
Definition: DataManager.cpp:201
void setInteriorFlagsP(geometry::InteriorFlags *pointer)
Definition: DataManager.cpp:104
util::StateBasedHelperFunctions * d_sbhelper_p
Pointer to StatebasedHelper function object.
Definition: DataManager.h:505
inp::OutputDeck * getOutputDeckP()
Definition: DataManager.cpp:47
std::vector< util::Point3 > * d_u_p
Pointer to the displacement vector.
Definition: DataManager.h:422
void setVelocityP(std::vector< util::Point3 > *pointer)
Definition: DataManager.cpp:61
loading::FLoading * d_fLoading_p
Pointer to force Loading object.
Definition: DataManager.h:512
std::vector< util::Point3 > * d_f_p
Pointer to the force vector.
Definition: DataManager.h:428
void setDamageFunctionP(std::vector< float > *pointer)
Definition: DataManager.cpp:209
void setWorkDoneP(std::vector< float > *pointer)
Definition: DataManager.cpp:197
void setReactionForceP(std::vector< util::Point3 > *pointer)
Definition: DataManager.cpp:181
void print(int nt=0, int lvl=0) const
Prints the information about the instance of the object.
Definition: DataManager.h:394
std::vector< float > * getKineticEnergyP()
Definition: DataManager.cpp:233
loading::ULoading * d_uLoading_p
Pointer to displacement Loading object.
Definition: DataManager.h:509
std::vector< float > * getFractureEnergyP()
Definition: DataManager.cpp:227
geometry::VolumeCorrection * d_volumeCorrection_p
Pointer to Volume Correction object.
Definition: DataManager.h:496
fe::Mesh * d_mesh_p
Pointer to Mesh object.
Definition: DataManager.h:490
std::vector< float > * getBBFractureEnergyP()
Definition: DataManager.cpp:219
void setTotalReactionForceP(std::vector< double > *pointer)
Definition: DataManager.cpp:189
std::vector< util::Point3 > * getDisplacementP()
Definition: DataManager.cpp:71
void setStressTensorP(std::vector< util::Matrix33 > *pointer)
Definition: DataManager.cpp:149
std::vector< double > * d_total_reaction_force_p
Pointer to the total reaction force vector.
Definition: DataManager.h:459
std::vector< float > * getPhiP()
Definition: DataManager.cpp:207
geometry::InteriorFlags * d_interiorFlags_p
Pointer to InteriorFlags object.
Definition: DataManager.h:502
void setForceLoadingP(loading::FLoading *pointer)
Definition: DataManager.cpp:126
std::vector< double > * d_dilatation_p
Dilatation of nodes.
Definition: DataManager.h:437
void setDisplacementP(std::vector< util::Point3 > *pointer)
Definition: DataManager.cpp:67
std::vector< util::Point3 > * getForceP()
Definition: DataManager.cpp:59
void setWeightedVolP(std::vector< double > *pointer)
Definition: DataManager.cpp:173
std::vector< double > * getDilatationP()
Definition: DataManager.cpp:169
std::vector< float > * d_e_p
Pointer to the strain energy vector.
Definition: DataManager.h:456
std::vector< float > * d_Z_p
Pointer to damage function at the nodes.
Definition: DataManager.h:468
void setBBFractureEnergyP(std::vector< float > *pointer)
Definition: DataManager.cpp:215
void setModelDeckP(inp::ModelDeck *pointer)
Definition: DataManager.cpp:37
void setDisplacementLoadingP(loading::ULoading *pointer)
Definition: DataManager.cpp:118
std::vector< util::Matrix33 > * getStrainTensorP()
Definition: DataManager.cpp:161
DataManager()
Constructor.
Definition: DataManager.cpp:30
std::vector< double > d_thetaX
Dilation.
Definition: DataManager.h:447
std::vector< double > d_mX
Weighted volume.
Definition: DataManager.h:453
std::vector< util::Matrix33 > * getStressTensorP()
Definition: DataManager.cpp:153
std::vector< float > * getStrainEnergyP()
Definition: DataManager.cpp:147
geometry::Fracture * d_fracture_p
Pointer to Fracture object.
Definition: DataManager.h:499
void setPhiP(std::vector< float > *pointer)
Definition: DataManager.cpp:203
geometry::VolumeCorrection * getVolumeCorrectionP()
Definition: DataManager.cpp:90
geometry::Neighbor * getNeighborP()
Definition: DataManager.cpp:83
std::vector< util::Point3 > * getVelocityP()
Definition: DataManager.cpp:65
std::vector< float > * d_w_p
Pointer to the Work done on each of the nodes.
Definition: DataManager.h:462
geometry::Fracture * getFractureP()
Definition: DataManager.cpp:94
void setBodyForceP(std::vector< util::Point3 > *pointer)
Definition: DataManager.cpp:49
void setFractureP(geometry::Fracture *pointer)
Definition: DataManager.cpp:96
void setDilatationP(std::vector< double > *pointer)
Definition: DataManager.cpp:165
inp::ModelDeck * getModelDeckP()
Definition: DataManager.cpp:41
void setKineticEnergyP(std::vector< float > *pointer)
Definition: DataManager.cpp:229
inp::ModelDeck * d_modelDeck_p
Model deck.
Definition: DataManager.h:405
std::vector< double > * d_weighted_vol_p
Weighted volume.
Definition: DataManager.h:440
void setStrainEnergyP(std::vector< float > *pointer)
Definition: DataManager.cpp:143
util::StateBasedHelperFunctions * getStateBasedHelperFunctionsP()
Definition: DataManager.cpp:114
void setOutputDeckP(inp::OutputDeck *pointer)
Definition: DataManager.cpp:43
std::vector< util::Point3 > * getBodyForceP()
Definition: DataManager.cpp:53
geometry::Neighbor * d_neighbor_p
Pointer to Neighbor object.
Definition: DataManager.h:493
loading::FLoading * getForceLoadingP()
Definition: DataManager.cpp:130
std::vector< util::Point3 > * getReactionForceP()
Definition: DataManager.cpp:185
std::vector< util::Point3 > * d_b_p
Pointer to the body force vector.
Definition: DataManager.h:419
void setForceP(std::vector< util::Point3 > *pointer)
Definition: DataManager.cpp:55
fe::Mesh * getMeshP()
Definition: DataManager.cpp:77
std::vector< float > * getDamageFunctionP()
Definition: DataManager.cpp:213
void setMeshP(fe::Mesh *d_mesh_p)
Definition: DataManager.cpp:75
std::vector< util::Point3 > * d_reaction_force_p
Pointer to the reaction force.
Definition: DataManager.h:431
std::vector< float > * d_phi_p
Pointer to damage function at the nodes.
Definition: DataManager.h:465
inp::OutputDeck * d_outputDeck_p
Output deck.
Definition: DataManager.h:408
A class for mesh data.
Definition: mesh.h:49
A class for fracture state of bonds.
Definition: fracture.h:31
A class to store interior/exterior flags of node.
Definition: interiorFlags.h:198
A class to store neighbor list and provide access to the list.
Definition: neighbor.h:32
A class to to compute the volume correction for state-based models.
Definition: volumeCorrection.h:44
A class to apply force boundary condition.
Definition: fLoading.h:25
A class to apply initial condition.
Definition: initialCondition.h:34
A class to apply displacement boundary condition.
Definition: uLoading.h:25
A class for global methods of state-based material models.
Definition: stateBasedHelperFunctions.h:50
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 and database related to loading.
Definition: DataManager.h:38
Collection of methods useful in simulation.
Definition: DataManager.h:50
Structure to read and store model related input data.
Definition: modelDeck.h:23
Structure to read input data for performing simulation output.
Definition: outputDeck.h:24