NLMech  0.1.0
fractureDeck.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_FRACTUREDECK_H
10 #define INP_FRACTUREDECK_H
11 
12 #include "util/compare.h" // compare utility
13 #include "util/point.h" // definition of Point3
14 #include "util/transfomation.h" // rotation of point
15 #include <cstdlib> // definition of siz_t
16 #include <vector>
17 #include "util/utilIO.h"
18 
19 namespace inp {
20 
22 struct EdgeCrack {
23 
31  int d_o;
32 
40  double d_theta;
41 
43  double d_l;
44 
46  double d_lt;
47 
49  double d_lb;
50 
53 
56 
58  int d_it;
59 
61  int d_ib;
62 
65 
68 
71 
74 
77 
80 
82  bool d_trackt;
83 
85  bool d_trackb;
86 
90 
93 
98  : d_o(1), d_theta(0.), d_l(0.), d_lt(0.), d_lb(0.), d_it(-1), d_ib(-1),
99  d_vt(util::Point3()), d_vb(util::Point3()),
100  d_pt(util::Point3()), d_pb(util::Point3()),
101  d_oldPt(util::Point3()), d_oldPb(util::Point3()),
102  d_initPt(util::Point3()), d_initPb(util::Point3()),
103  d_trackt(false), d_trackb(false), d_activationTime(-1.),
104  d_crackAcrivated(false) {};
105 
116  double theta = 0.0) {
117 
118  if (o == -1) {
119 
120  // straight crack along y-axis
121  return util::compare::definitelyLessThan(p.d_y, pb.d_y) or
123 
124  } else if (o == 1) {
125 
126  // straight crack along x-axis ===> pb == pl, pt == pr
127  return util::compare::definitelyLessThan(p.d_x, pb.d_x) or
129 
130  } else if (o == 0) {
131 
132  // straight crack at theta angle with x-axis ===> pb == pl, pt == pr
133 
134  // 1. pb ==> (0,0), pt ==> pt - pb, p ==> p - pb
135  // 2. Apply CW rotation to new pt and new p so that crack line
136  // after rotation is simply along x-axis with left point at
137  // origin and right point at transformed new pt
138 
139  util::Point3 pmap = util::transformation::rotateCW2D(p - pb, theta);
140  util::Point3 ptmap = util::transformation::rotateCW2D(pt - pb, theta);
141 
142  return util::compare::definitelyLessThan(pmap.d_x, 0.0) or
144  }
145 
146  return true;
147  }
148 
157 
158  // algorithm is same for any orientation of crack
159  // pt
160  // o
161  // /
162  // p /
163  // o /
164  // /
165  // /
166  // /
167  // o
168  // pb
169  double a = 0.5 * ((pt.d_x - pb.d_x) * (p.d_y - pb.d_y) -
170  (p.d_x - pb.d_x) * (pt.d_y - pb.d_y));
171 
172  // crack is closer to left nodes
173  return !util::compare::definitelyLessThan(a, 0.0);
174  }
175 
184 
185  return !this->ptLeftside(p, pb, pt);
186  }
187 
197  std::string printStr(int nt = 0, int lvl = 0) const {
198  auto tabS = util::io::getTabS(nt);
199  std::ostringstream oss;
200  oss << tabS << "------- EdgeCrack --------" << std::endl << std::endl;
201  oss << tabS << "Orientation = " << d_o << std::endl;
202  oss << tabS << "Point 1 = " << d_pb << std::endl;
203  oss << tabS << "Point 2 = " << d_pt << std::endl;
204  oss << tabS << "Length = " << d_l << std::endl;
205  oss << tabS << "Angle = " << d_theta << std::endl;
206  oss << tabS << std::endl;
207 
208  return oss.str();
209  };
210 
218  void print(int nt = 0, int lvl = 0) const { std::cout << printStr(nt, lvl); };
219 };
220 
227 struct FractureDeck {
228 
230  std::vector<inp::EdgeCrack> d_cracks;
231 
236 
246  std::string printStr(int nt = 0, int lvl = 0) const {
247  auto tabS = util::io::getTabS(nt);
248  std::ostringstream oss;
249  oss << tabS << "------- FractureDeck --------" << std::endl << std::endl;
250  oss << tabS << "Number of cracks = " << d_cracks.size() << std::endl;
251  for (size_t i=0; i<d_cracks.size(); i++) {
252  oss << tabS << "Crack " << i + 1 << " information" << std::endl;
253  oss << d_cracks[i].printStr(nt + 1, lvl);
254  oss << std::endl;
255  }
256  oss << tabS << std::endl;
257 
258  return oss.str();
259  };
260 
268  void print(int nt = 0, int lvl = 0) const { std::cout << printStr(nt, lvl); };
269 };
270 
273 } // namespace inp
274 
275 #endif // INP_FRACTUREDECK_H
Collection of methods and database related to input.
Definition: main.cpp:21
bool definitelyGreaterThan(const double &a, const double &b)
Compares if a > to b.
Definition: compare.cpp:25
bool definitelyLessThan(const double &a, const double &b)
Compares if a is < to b.
Definition: compare.cpp:30
std::string getTabS(int nt)
Generate a string contaning nt tabs.
Definition: utilIO.h:26
std::vector< double > rotateCW2D(const std::vector< double > &x, const double &theta)
Rotates a vector in xy-plane in clockwise direction.
Definition: transformation.cpp:13
Collection of methods useful in simulation.
Definition: DataManager.h:50
A structure to edge crack of any orientation.
Definition: fractureDeck.h:22
std::string printStr(int nt=0, int lvl=0) const
Returns the string containing information about the instance of the object.
Definition: fractureDeck.h:197
void print(int nt=0, int lvl=0) const
Prints the information about the instance of the object.
Definition: fractureDeck.h:218
double d_lb
Current length of crack on bottom side (left side)
Definition: fractureDeck.h:49
double d_lt
Current length of crack on top side (right side)
Definition: fractureDeck.h:46
bool ptRightside(util::Point3 p, util::Point3 pb, util::Point3 pt)
Checks if point lies on left(top) or right(bottom) of crack.
Definition: fractureDeck.h:183
util::Point3 d_vb
Velocity of bottom (left) crack tip.
Definition: fractureDeck.h:55
util::Point3 d_initPt
Top (right) crack tip location (initial)
Definition: fractureDeck.h:76
bool d_crackAcrivated
Is crack activated.
Definition: fractureDeck.h:92
double d_l
Total current length of crack.
Definition: fractureDeck.h:43
bool d_trackb
Track bottom point (left) of crack.
Definition: fractureDeck.h:85
bool d_trackt
Track top point (right) of crack.
Definition: fractureDeck.h:82
double d_theta
Angle the edge crack makes with the horizontal axis.
Definition: fractureDeck.h:40
util::Point3 d_oldPt
Top (right) crack tip location (old)
Definition: fractureDeck.h:70
int d_it
Closest id of node to top (right) crack tip.
Definition: fractureDeck.h:58
util::Point3 d_pb
Bottom (left) crack tip location.
Definition: fractureDeck.h:67
EdgeCrack()
Constructor.
Definition: fractureDeck.h:97
int d_ib
Closest id of node to bottom (left) crack tip.
Definition: fractureDeck.h:61
util::Point3 d_initPb
Bottom (left) crack tip location (initial)
Definition: fractureDeck.h:79
bool ptLeftside(util::Point3 p, util::Point3 pb, util::Point3 pt)
Checks if point lies on left(top) or right(bottom) of crack.
Definition: fractureDeck.h:156
util::Point3 d_oldPb
Bottom (left) crack tip location (old)
Definition: fractureDeck.h:73
bool ptOutside(util::Point3 p, int o, util::Point3 pb, util::Point3 pt, double theta=0.0)
Checks if point lies outside the crack.
Definition: fractureDeck.h:115
int d_o
Orientation of crack.
Definition: fractureDeck.h:31
util::Point3 d_vt
Velocity of top (right) crack tip.
Definition: fractureDeck.h:52
double d_activationTime
Activation time (for example when we abruptly slit the domain at specified time)
Definition: fractureDeck.h:89
util::Point3 d_pt
Top (right) crack tip location.
Definition: fractureDeck.h:64
Structure to read and store fracture related input data.
Definition: fractureDeck.h:227
std::vector< inp::EdgeCrack > d_cracks
Vector of pre-crack data.
Definition: fractureDeck.h:230
std::string printStr(int nt=0, int lvl=0) const
Returns the string containing information about the instance of the object.
Definition: fractureDeck.h:246
FractureDeck()
Constructor.
Definition: fractureDeck.h:235
void print(int nt=0, int lvl=0) const
Prints the information about the instance of the object.
Definition: fractureDeck.h:268
A structure to represent 3d vectors.
Definition: point.h:29
double d_x
the x coordinate
Definition: point.h:32
double d_y
the y coordinate
Definition: point.h:35