High Performance Plasticity  0.5.0
casesUtils.h
Go to the documentation of this file.
1 
7 #ifndef HPP_CASESUTILS_H
8 #define HPP_CASESUTILS_H
9 
10 #include <hpp/config.h>
11 #include <hpp/tensor.h>
12 #include <hpp/rotation.h>
13 #include <functional>
14 #include <stdexcept>
15 #include <memory>
16 
17 namespace hpp
18 {
19 
29 template <typename U>
31  Tensor2<U> F = identityTensor2<U>(3);
32  F(0,1) = shear_rate*t;
33  return F;
34 }
35 
45 template <typename U>
47  Tensor2<U> L(3,3);
48  L(0,1) = shear_rate;
49  return L;
50 }
51 
63 template <typename U>
65  Tensor2<U> F(3,3);
66  F(0,0) = std::exp(-0.5*comp_rate*t);
67  F(1,1) = std::exp(-0.5*comp_rate*t);
68  F(2,2) = std::exp(1.0*comp_rate*t);
69  return F;
70 }
71 
83 template <typename U>
85  Tensor2<U> L(3,3);
86  L(0,0) = -0.5*comp_rate;
87  L(1,1) = -0.5*comp_rate;
88  L(2,2) = 1.0*comp_rate;
89  return L;
90 }
91 
103 template <typename U>
105  Tensor2<U> F(3,3);
106  F(0,0) = std::exp(1.0*comp_rate*t);
107  F(1,1) = 1.0;
108  F(2,2) = std::exp(-1.0*comp_rate*t);
109  return F;
110 }
111 
122 template <typename U>
124  Tensor2<U> L(3,3);
125  L(0,0) = 1.0*comp_rate;
126  L(2,2) = -1.0*comp_rate;
127  return L;
128 }
129 
137 template <typename U>
139  Tensor2<U> F = identityTensor2<U>(3);
140  return F;
141 }
142 
150 template <typename U>
152  Tensor2<U> L(3,3);
153  return L;
154 }
155 
164 template <typename U>
166  public:
172  virtual void generateNext(Tensor2<U>& rotMatrix) = 0;
178  virtual void generateNext(EulerAngles<U>& angles) {
179  Tensor2<U> rotMatrix(3,3);
180  this->generateNext(rotMatrix);
181  angles = getEulerZXZAngles(rotMatrix);
182  }
186  virtual ~OrientationGenerator() {};
187 };
188 
197 template <typename U>
199  public:
209  virtual void generateNext(Tensor2<U>& rotMatrix) {
210  randomRotationTensorInPlace(3, rotMatrix, true);
211  }
212 };
213 
224 template <typename U>
226  public:
232  GridOrientationGenerator(int nTheta=8, int nPhi=4);
238  virtual void generateNext(EulerAngles<U>& angles) {
239  // Get angles
240  U theta = iTheta*dTheta;
241  U phi = iPhi*dPhi;
242  angles = polarToEuler(theta, phi);
243 
244  // Prepare next angle
245  iPhi++;
246  if (iPhi == nPhi) {
247  iPhi = 0;
248  iTheta++;
249  iTheta = iTheta % nTheta;
250  }
251  };
257  virtual void generateNext(Tensor2<U>& rotMatrix) {
258  EulerAngles<U> angles;
259  this->generateNext(angles);
260  rotMatrix = EulerZXZRotationMatrix(angles);
261  };
262 
263  private:
264  int nTheta;
265  int nPhi;
266  int iTheta = 0;
267  int iPhi = 0;
268  U dTheta;
269  U dPhi;
270 };
271 
279 template <typename U>
280 class Experiment {
281 
282  public:
287  explicit Experiment(std::string experimentName);
288  U tStart;
289  U tEnd;
291  std::function<Tensor2<U>(U)> F_of_t;
296  std::function<Tensor2<U>(U)> L_of_t;
298  std::shared_ptr<OrientationGenerator<U>> orientationGenerator;
299 
300  // Generating next orientation
301  void generateNextOrientationInPlace(Tensor2<U>& rotMatrix) {this->orientationGenerator->generateNext(rotMatrix);}
302  void generateNextOrientationInPlace(EulerAngles<U>& angles) {this->orientationGenerator->generateNext(angles);}
304  EulerAngles<U> angles;
305  this->orientationGenerator->generateNext(angles);
306  return angles;
307  }
308 
309  // Getters/setters (mainly intended for Python interface)
310  U getTStart() const {return tStart;}
311  U getTEnd() const {return tEnd;}
312  U getStrainRate() const {return strainRate;}
313  std::function<Tensor2<U>(U)> getF_of_t() const {return F_of_t;}
314  std::function<Tensor2<U>(U)> getL_of_t() const {return L_of_t;}
315 };
316 
317 } // END NAMESPACE hpp
318 
319 #endif /* HPP_CASESUTILS_H */
Tensor2< U > staticDeformationGradient(U t)
The deformation gradient for no deformation.
Definition: casesUtils.h:138
virtual void generateNext(Tensor2< U > &rotMatrix)
Generate the next orientation in the sequence.
Definition: casesUtils.h:257
U getTEnd() const
Definition: casesUtils.h:311
Tensor2< U > simpleCompressionDeformationGradient(U t, U comp_rate)
The deformation gradient from a simple compression.
Definition: casesUtils.h:64
An abstract class for reproducibly generating a sequence of orientations.
Definition: casesUtils.h:165
virtual void generateNext(Tensor2< U > &rotMatrix)=0
Generate the next orientation in the sequence.
U getStrainRate() const
Definition: casesUtils.h:312
RandomOrientationGenerator()
Default constructor.
Definition: casesUtils.h:203
Definition: casesUtils.cpp:4
Tensor2< T > EulerZXZRotationMatrix(T alpha, T beta, T gamma)
Definition: rotation.h:118
virtual void generateNext(EulerAngles< U > &angles)
Generate the next orientation in the sequence.
Definition: casesUtils.h:238
U dTheta
the spacing between azimuthal angles on the grid
Definition: casesUtils.h:268
Generates orientations on a fixed grid over an azimuthal angle of and a zenithal angle of ...
Definition: casesUtils.h:225
EulerAngles< U > generateNextOrientationAngles()
Definition: casesUtils.h:303
Tensor2< U > planeStrainCompressionDeformationGradient(U t, U comp_rate)
The deformation gradient from a plane strain compression.
Definition: casesUtils.h:104
Header file for rotation classes and functions.
int nPhi
the number of grid points for the zenithal angle
Definition: casesUtils.h:265
U strainRate
Definition: casesUtils.h:290
Header file for tensor classes.
std::function< Tensor2< U >U)> getF_of_t() const
Definition: casesUtils.h:313
A class for reproducibly generating a sequence of random orientations.
Definition: casesUtils.h:198
virtual void generateNext(EulerAngles< U > &angles)
Generate the next orientation in the sequence.
Definition: casesUtils.h:178
U tEnd
the end time of the experiment
Definition: casesUtils.h:289
Tensor2< U > simpleCompressionVelocityGradient(U t, U comp_rate)
The velocity gradient from a simple compression.
Definition: casesUtils.h:84
Tensor2< U > planeStrainCompressionVelocityGradient(U t, U comp_rate)
The velocity gradient from a plane strain compression.
Definition: casesUtils.h:123
void generateNextOrientationInPlace(Tensor2< U > &rotMatrix)
Definition: casesUtils.h:301
EulerAngles< T > polarToEuler(T theta, T phi)
Convert polar angles to Euler angles.
Definition: rotation.h:76
void generateNextOrientationInPlace(EulerAngles< U > &angles)
Definition: casesUtils.h:302
U tStart
the start time of the experiment
Definition: casesUtils.h:288
Tensor2< U > simpleShearDeformationGradient(U t, U shear_rate)
The deformation gradient from a simple shear.
Definition: casesUtils.h:30
virtual void generateNext(Tensor2< U > &rotMatrix)
Generate the next orientation in the sequence.
Definition: casesUtils.h:209
std::shared_ptr< OrientationGenerator< U > > orientationGenerator
the generator of crystal orientations for the experiment
Definition: casesUtils.h:298
int nTheta
the number of grid points for the azimuthal angle
Definition: casesUtils.h:261
Tensor2< U > simpleShearVelocityGradient(U t, U shear_rate)
The velocity gradient from a simple shear.
Definition: casesUtils.h:46
std::function< Tensor2< U >U)> getL_of_t() const
Definition: casesUtils.h:314
U getTStart() const
Definition: casesUtils.h:310
Tensor2< U > staticVelocityGradient(U t)
The velocity gradient for no deformation.
Definition: casesUtils.h:151
U dPhi
the spacing between zenithal angles on the grid
Definition: casesUtils.h:269
EulerAngles< T > getEulerZXZAngles(Tensor2< T > R)
Get Euler angles from rotation matrix.
Definition: rotation.h:149
Definition: casesUtils.h:280
void randomRotationTensorInPlace(unsigned int dim, Tensor2< T > &A, bool defaultSeed=false)
Generate a random rotation tensor.
Definition: rotation.h:321
Definition: rotation.h:43
virtual ~OrientationGenerator()
Default destructor.
Definition: casesUtils.h:186