NLMech  0.1.0
compute.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 TOOLS_PP_COMPUTE_H
10 #define TOOLS_PP_COMPUTE_H
11 
12 #include "util.h"
13 #include "util/matrix.h" // definition of SymMatrix3
14 
15 // forward declarations of decks
16 namespace inp {
17 class Input;
18 struct MaterialDeck;
19 struct FractureDeck;
20 struct ModelDeck;
21 struct OutputDeck;
22 } // namespace inp
23 
24 namespace fe {
25 class Mesh;
26 }
27 
28 namespace geometry {
29 class Fracture;
30 class Neighbor;
31 }
32 
33 namespace material {
34 namespace pd {
35 class BaseMaterial;
36 }
37 }
38 
39 namespace rw {
40 namespace writer {
41 class Writer;
42 }
43 }
44 
45 namespace data {
46 class DataManager;
47 }
48 
49 namespace tools {
50 
58 namespace pp {
59 
75 class Compute {
76 
77 public:
83  explicit Compute(const std::string &filename);
84 
88  void init();
89 
93  void finalize();
94 
101  void readComputeInstruction(const std::string &set,
103 
111  void readCrackTipData(const std::string &filename, int crack_id,
112  std::vector<tools::pp::CrackTipData> *data, bool
113  is_inclined = false);
114 
121  void initWriter(rw::writer::Writer *writer,
122  std::vector<util::Point3> *u);
123 
137  void transformU(rw::writer::Writer *writer);
138 
148  void transformV(rw::writer::Writer *writer);
149 
195  void computeStrain(rw::writer::Writer *writer);
196 
210  void computeDamage(rw::writer::Writer *writer,
211  std::vector<double> *Z, bool perf_out = false);
212 
230  void findCrackTip(std::vector<double> *Z, rw::writer::Writer
231  *writer);
232 
356  void computeJIntegral();
357 
363 
366 private:
378  void addUniqueToList(size_t i, std::vector<size_t> *list);
379 
392  size_t findNode(const util::Point3 &x, const std::vector<util::Point3> *nodes,
393  const std::vector<util::Point3> *u = nullptr);
394 
412  void
413  listElemsAndNodesInDomain(const std::pair<util::Point3, util::Point3> &cd,
414  const double &tol, const double &tol_elem,
415  std::vector<size_t> *nodes,
416  std::vector<size_t> *elements,
417  bool calc_in_ref);
418 
429  void decomposeSearchNodes(const std::pair<util::Point3, util::Point3> &cd,
430  std::vector<size_t> *nodes,
431  std::vector<size_t> *nodes_new);
432 
455  util::Point3 &vp,
456  const std::vector<size_t> &ids,
457  bool calc_in_ref,
458  bool check_only = false);
459 
489  void interpolateUV(const util::Point3 &p, util::Point3 &up, util::Point3 &vp,
490  const std::vector<size_t> *nodes,
491  const std::vector<size_t> *elements, bool calc_in_ref);
492 
493  size_t interpolateUVNodes(const util::Point3 &p, util::Point3 &up,
494  util::Point3 &vp,
495  const std::vector<size_t> *nodes, bool calc_in_ref);
496 
515  void getContourContribJInt(const util::Point3 &p,
516  const std::vector<size_t> *nodes,
517  const std::vector<size_t> *elements,
518  const util::Point3 &normal,
519  double &pd_energy,
520  double &kinetic_energy,
521  double &elastic_energy,
522  util::Point3 &dot_u,
523  bool calc_in_ref,
524  const tools::pp::CrackTipData &ctip);
525 
573  void updateCrack(const double &time, const std::vector<double> *Z);
574 
576  void crackOutput();
577 
595  inp::EdgeCrack &crack,
596  std::vector<std::pair<util::Point3, util::Point3>> &rects_t,
597  std::vector<std::pair<util::Point3, util::Point3>> &rects_b,
598  std::vector<std::vector<size_t>> &nodes_t,
599  std::vector<std::vector<size_t>> &nodes_b,
600  std::vector<std::vector<double>> &Z_t,
601  std::vector<std::vector<double>> &Z_b, const std::vector<double> *Z);
602 
611  void addNewCrackTip(inp::EdgeCrack &crack, util::Point3 pnew, double time,
612  bool is_top);
613 
629  inp::EdgeCrack &crack, const double &max_Z,
630  const std::vector<std::pair<util::Point3, util::Point3>> &rects,
631  const std::vector<std::vector<size_t>> &nodes,
632  const std::vector<std::vector<double>> &Zs, const std::vector<double> *Z,
633  bool is_top);
634 
641  void safeExit(const std::string &err_message, rw::writer::Writer *writer = nullptr);
642 
646  std::string d_inpFilename;
647 
649  std::vector<tools::pp::InstructionData> d_computeData;
650 
652  std::string d_fileZ;
653 
655  std::string d_tagZ;
656 
658  int d_nOut;
659 
661  size_t d_nC;
662 
664  double d_time;
665 
668 
670  std::string d_outFilename;
671 
673  std::string d_outPath;
674 
676  std::string d_outPreTag;
677 
682  std::string d_simInpFilename;
683 
688  std::string d_simOutFilename;
689 
703 
706 
708  bool d_uPlus;
709 
711  size_t d_dtN;
712 
714  size_t d_dtStart;
715 
720  size_t d_dtEnd;
721 
726  std::string d_fnErrMsg;
727 
730 
733 
735  std::vector<util::Point3> d_u;
736 
738  std::vector<util::Point3> d_v;
739 
741  std::vector<double> d_Z;
742 
748  std::vector<double> d_thetaX;
749 
754  std::vector<double> d_mX;
755 
757  //inp::ModelDeck *d_modelDeck_p;
758 
761 
764 
767 
769  //geometry::Fracture *d_fracture_p;
770 
773 
776 
779 
780 };
781 
782 } // namespace pp
783 
784 } // namespace tools
785 
786 #endif // TOOLS_PP_COMPUTE_H
Data manager to collect the global simulation data.
Definition: DataManager.h:58
A class to read input file.
Definition: input.h:63
A base class providing methods to compute energy density and force.
Definition: baseMaterial.h:26
A interface class writing data.
Definition: writer.h:41
Processes simulation results and computes post-processing quantities.
Definition: compute.h:75
void computeStrain(rw::writer::Writer *writer)
Compute strain and stress.
Definition: compute.cpp:1019
bool d_needDamageZ
Need damage at nodes.
Definition: compute.h:729
inp::MaterialDeck * d_matDeck_p
Material deck.
Definition: compute.h:766
void finalize()
Finalize.
Definition: compute.cpp:485
inp::FractureDeck * d_fractureDeck_p
Fracture deck.
Definition: compute.h:763
void readComputeInstruction(const std::string &set, tools::pp::InstructionData *data)
Read compute instruction.
Definition: compute.cpp:515
void computeJIntegral()
Compute J integral.
Definition: compute.cpp:1230
size_t d_dtOutChange
Time step at which output interval changes.
Definition: compute.h:702
int d_nOut
Current active simulation output step.
Definition: compute.h:658
Compute(const std::string &filename)
Constructor.
Definition: compute.cpp:133
std::string d_tagZ
Tag of damage function Z in vtu file.
Definition: compute.h:655
size_t findNode(const util::Point3 &x, const std::vector< util::Point3 > *nodes, const std::vector< util::Point3 > *u=nullptr)
Find node closest to given point.
Definition: compute.cpp:1893
util::Point3 findTipInRects(inp::EdgeCrack &crack, const double &max_Z, const std::vector< std::pair< util::Point3, util::Point3 >> &rects, const std::vector< std::vector< size_t >> &nodes, const std::vector< std::vector< double >> &Zs, const std::vector< double > *Z, bool is_top)
Finds new crack tip location from the rectangle data.
Definition: compute.cpp:2657
void computeDamage(rw::writer::Writer *writer, std::vector< double > *Z, bool perf_out=false)
Compute damage at nodes.
Definition: compute.cpp:1170
void initWriter(rw::writer::Writer *writer, std::vector< util::Point3 > *u)
Create new writer object if it is not already done so.
Definition: compute.cpp:867
bool d_uPlus
Specify if we consider u(n+1) or u(n)
Definition: compute.h:708
material::pd::BaseMaterial * d_material_p
Pointer to Material object.
Definition: compute.h:775
void readCrackTipData(const std::string &filename, int crack_id, std::vector< tools::pp::CrackTipData > *data, bool is_inclined=false)
Read crack tip data.
Definition: compute.cpp:830
void decomposeSearchNodes(const std::pair< util::Point3, util::Point3 > &cd, std::vector< size_t > *nodes, std::vector< size_t > *nodes_new)
Decompose node list into two lists.
Definition: compute.cpp:1968
bool d_writerReady
State of writer class.
Definition: compute.h:705
double d_time
Current time.
Definition: compute.h:664
size_t d_dtStart
Global output start step. Default is 1.
Definition: compute.h:714
void addNewCrackTip(inp::EdgeCrack &crack, util::Point3 pnew, double time, bool is_top)
Adds new tip to crack, computes velocity and crack length.
Definition: compute.cpp:2629
inp::OutputDeck * d_outputDeck_p
Model deck.
Definition: compute.h:760
std::string d_simInpFilename
Simulation input filename.
Definition: compute.h:682
size_t d_dtEnd
Global output start step. Default is maximum output step for the simulation, i.e. d_dtN.
Definition: compute.h:720
std::string d_simOutFilename
Mesh filename to get displacement and velocity (results of simulation)
Definition: compute.h:688
void findCrackTip(std::vector< double > *Z, rw::writer::Writer *writer)
Find crack tip location and velocity.
Definition: compute.cpp:1218
void safeExit(const std::string &err_message, rw::writer::Writer *writer=nullptr)
Exits code safely in the event there is error.
Definition: compute.cpp:501
void transformU(rw::writer::Writer *writer)
Transform displacement.
Definition: compute.cpp:892
std::vector< util::Point3 > d_u
Displacement of nodes.
Definition: compute.h:735
std::string d_fileZ
Filename (if any) to read damage at nodes.
Definition: compute.h:652
bool triCheckAndInterpolateUV(const util::Point3 &p, util::Point3 &up, util::Point3 &vp, const std::vector< size_t > &ids, bool calc_in_ref, bool check_only=false)
Interpolates displacement and velocity at given point in triangle.
Definition: compute.cpp:1985
void getRectsAndNodesForCrackTip(inp::EdgeCrack &crack, std::vector< std::pair< util::Point3, util::Point3 >> &rects_t, std::vector< std::pair< util::Point3, util::Point3 >> &rects_b, std::vector< std::vector< size_t >> &nodes_t, std::vector< std::vector< size_t >> &nodes_b, std::vector< std::vector< double >> &Z_t, std::vector< std::vector< double >> &Z_b, const std::vector< double > *Z)
Creates sequence of rectangles for crack tip search, and finds the ids of nodes in each rectangle and...
Definition: compute.cpp:2445
std::string d_inpFilename
Input filename for compute instruction.
Definition: compute.h:646
std::vector< tools::pp::InstructionData > d_computeData
List of compute sets.
Definition: compute.h:649
std::string d_outFilename
Current active output filename to save postprocessed data.
Definition: compute.h:670
void listElemsAndNodesInDomain(const std::pair< util::Point3, util::Point3 > &cd, const double &tol, const double &tol_elem, std::vector< size_t > *nodes, std::vector< size_t > *elements, bool calc_in_ref)
Find node within contour and elements intersecting contour.
Definition: compute.cpp:1917
std::vector< double > d_Z
Damage of nodes.
Definition: compute.h:741
void updateCrack(const double &time, const std::vector< double > *Z)
Updates crack tip location and crack velocity.
Definition: compute.cpp:2353
void crackOutput()
Performs output of crack tip data.
Definition: compute.cpp:2400
void transformV(rw::writer::Writer *writer)
Transform velocity.
Definition: compute.cpp:921
size_t d_dtN
Total number of output files to process.
Definition: compute.h:711
bool d_fnSuccess
Save success/failur of function.
Definition: compute.h:725
std::vector< util::Point3 > d_v
Velocity of nodes.
Definition: compute.h:738
tools::pp::InstructionData * d_currentData
Current active compute data.
Definition: compute.h:667
void getContourContribJInt(const util::Point3 &p, const std::vector< size_t > *nodes, const std::vector< size_t > *elements, const util::Point3 &normal, double &pd_energy, double &kinetic_energy, double &elastic_energy, util::Point3 &dot_u, bool calc_in_ref, const tools::pp::CrackTipData &ctip)
Computes contribution to energy into crack from the quadrature point on contour.
Definition: compute.cpp:2138
size_t d_nC
Current active compute set.
Definition: compute.h:661
std::string d_outPath
Output path where postprocessed files are created.
Definition: compute.h:673
std::vector< double > d_mX
Weighted volume.
Definition: compute.h:754
void addUniqueToList(size_t i, std::vector< size_t > *list)
Inserts element to list if not present.
Definition: compute.cpp:1886
std::string d_outPreTag
Initial tag for each post-processing output.
Definition: compute.h:676
inp::Input * d_input_p
Pointer to Fracture object.
Definition: compute.h:772
void init()
Initialize the class based on input data.
Definition: compute.cpp:307
std::vector< double > d_thetaX
Dilation.
Definition: compute.h:748
void computeJIntegralAngledCrack()
Compute J integral where crack line is angled as opposed to either horizontal or verticle line.
Definition: compute.cpp:1878
data::DataManager * d_dataManager_p
Pointer to the Data Mananger.
Definition: compute.h:778
bool d_needNeighborList
Need neighbor list.
Definition: compute.h:732
void interpolateUV(const util::Point3 &p, util::Point3 &up, util::Point3 &vp, const std::vector< size_t > *nodes, const std::vector< size_t > *elements, bool calc_in_ref)
Interpolates displacement and velocity at given point.
Definition: compute.cpp:2025
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 reading and writing.
Definition: QuasiStaticModel.h:60
Namespace consisting of libraries which are used either in pre-process or in post-process.
Definition: compute.h:49
A structure to edge crack of any orientation.
Definition: fractureDeck.h:22
Structure to read and store fracture related input data.
Definition: fractureDeck.h:227
Structure to read and store material related data.
Definition: materialDeck.h:226
Structure to read input data for performing simulation output.
Definition: outputDeck.h:24
Datatype to hold crack tip data.
Definition: util.h:112
Datatype to hold instructions for post-processing operation.
Definition: util.h:325
A structure to represent 3d vectors.
Definition: point.h:29