-
Notifications
You must be signed in to change notification settings - Fork 1
/
Cartesian2d.h
373 lines (284 loc) · 10.3 KB
/
Cartesian2d.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
/*
(c - MIT) T.W.J. de Geus (Tom) | www.geus.me | github.com/tdegeus/GMatElastoPlasticQPot
*/
#ifndef GMATELASTOPLASTICQPOT_CARTESIAN2D_H
#define GMATELASTOPLASTICQPOT_CARTESIAN2D_H
#include "config.h"
namespace GMatElastoPlasticQPot {
namespace Cartesian2d {
// Alias
#if defined(_WIN32) || defined(_WIN64)
using Tensor2 = xt::xtensor<double, 2>;
using Tensor4 = xt::xtensor<double, 4>;
#else
using Tensor2 = xt::xtensor_fixed<double, xt::xshape<2, 2>>;
using Tensor4 = xt::xtensor_fixed<double, xt::xshape<2, 2, 2, 2>>;
#endif
// Unit tensors
inline Tensor2 I2();
inline Tensor4 II();
inline Tensor4 I4();
inline Tensor4 I4rt();
inline Tensor4 I4s();
inline Tensor4 I4d();
// Tensor decomposition
template <class T, class U>
inline void hydrostatic(const T& A, U& B);
template <class T>
inline auto Hydrostatic(const T& A);
template <class T, class U>
inline void deviatoric(const T& A, U& B);
template <class T>
inline auto Deviatoric(const T& A);
// Equivalent strain
template <class T, class U>
inline void epsd(const T& A, U& B);
template <class T>
inline auto Epsd(const T& A);
// Equivalent stress
template <class T, class U>
inline void sigd(const T& A, U& B);
template <class T>
inline auto Sigd(const T& A);
// Material point
class Elastic
{
public:
// Constructors
Elastic() = default;
Elastic(double K, double G);
// Parameters
double K() const;
double G() const;
// Set strain
template <class T>
void setStrain(const T& Eps);
template <class T>
void setStrainIterator(const T& begin); // presumes: contiguous + row-major & symmetric
// Stress (no allocation, overwrites "Sig" / writes to "begin")
template <class T>
void stress(T& Sig) const;
template <class T>
void stressIterator(const T& begin) const; // presumes: contiguous + row-major
// Tangent (no allocation, overwrites "C")
template <class T>
void tangent(T& C) const;
// Auto-allocation
Tensor2 Stress() const;
Tensor4 Tangent() const;
// Return current state
double energy() const; // potential energy
private:
double m_K; // bulk modulus
double m_G; // shear modulus
std::array<double, 4> m_Eps; // strain tensor [xx, xy, yx, yy]
std::array<double, 4> m_Sig; // stress tensor [xx, xy, yx, yy]
};
// Material point
class Cusp
{
public:
// Constructors
Cusp() = default;
Cusp(double K, double G, const xt::xtensor<double, 1>& epsy, bool init_elastic = true);
// Parameters
double K() const;
double G() const;
xt::xtensor<double, 1> epsy() const;
// Set strain
template <class T>
void setStrain(const T& Eps);
template <class T>
void setStrainIterator(const T& begin); // presumes: contiguous + row-major & symmetric
// Stress (no allocation, overwrites "Sig" / writes to "begin")
template <class T>
void stress(T& Sig) const;
template <class T>
void stressIterator(const T& begin) const; // presumes: contiguous + row-major
// Tangent (no allocation, overwrites "C")
template <class T>
void tangent(T& C) const;
// Auto-allocation
Tensor2 Stress() const;
Tensor4 Tangent() const;
// Return current state
size_t currentIndex() const; // yield index
double currentYieldLeft() const; // yield strain left epsy[index]
double currentYieldRight() const; // yield strain right epsy[index + 1]
double epsp() const; // "plastic strain" (mean of currentYieldLeft and currentYieldRight)
double energy() const; // potential energy
// Check that 'the particle' is at least "n" wells from the far-left/right
bool checkYieldBoundLeft(size_t n = 0) const;
bool checkYieldBoundRight(size_t n = 0) const;
private:
double m_K; // bulk modulus
double m_G; // shear modulus
QPot::Static m_yield; // potential energy landscape
std::array<double, 4> m_Eps; // strain tensor [xx, xy, yx, yy]
std::array<double, 4> m_Sig; // stress tensor [xx, xy, yx, yy]
};
// Material point
class Smooth
{
public:
// Constructors
Smooth() = default;
Smooth(double K, double G, const xt::xtensor<double, 1>& epsy, bool init_elastic = true);
// Parameters
double K() const;
double G() const;
xt::xtensor<double, 1> epsy() const;
// Set strain
template <class T>
void setStrain(const T& Eps);
template <class T>
void setStrainIterator(const T& begin); // presumes: contiguous + row-major & symmetric
// Stress (no allocation, overwrites "Sig" / writes to "begin")
template <class T>
void stress(T& Sig) const;
template <class T>
void stressIterator(const T& begin) const; // presumes: contiguous + row-major
// Tangent (no allocation, overwrites "C")
template <class T>
void tangent(T& C) const;
// Auto-allocation
Tensor2 Stress() const;
Tensor4 Tangent() const;
// Return current state
size_t currentIndex() const; // yield index
double currentYieldLeft() const; // yield strain left epsy[index]
double currentYieldRight() const; // yield strain right epsy[index + 1]
double epsp() const; // "plastic strain" (mean of currentYieldLeft and currentYieldRight)
double energy() const; // potential energy
// Check that 'the particle' is at least "n" wells from the far-left/right
bool checkYieldBoundLeft(size_t n = 0) const;
bool checkYieldBoundRight(size_t n = 0) const;
private:
double m_K; // bulk modulus
double m_G; // shear modulus
QPot::Static m_yield; // potential energy landscape
std::array<double, 4> m_Eps; // strain tensor [xx, xy, yx, yy]
std::array<double, 4> m_Sig; // stress tensor [xx, xy, yx, yy]
};
// Material identifier
struct Type {
enum Value {
Unset,
Elastic,
Cusp,
Smooth,
};
};
// Array of material points
template <size_t rank>
class Array
{
public:
// Constructors
Array() = default;
Array(const std::array<size_t, rank>& shape);
// Shape
std::array<size_t, rank> shape() const;
// Type
xt::xtensor<size_t, rank> type() const;
xt::xtensor<size_t, rank> isElastic() const;
xt::xtensor<size_t, rank> isPlastic() const;
xt::xtensor<size_t, rank> isCusp() const;
xt::xtensor<size_t, rank> isSmooth() const;
// Parameters
xt::xtensor<double, rank> K() const;
xt::xtensor<double, rank> G() const;
// Array of unit tensors
xt::xtensor<double, rank + 2> I2() const;
xt::xtensor<double, rank + 4> II() const;
xt::xtensor<double, rank + 4> I4() const;
xt::xtensor<double, rank + 4> I4rt() const;
xt::xtensor<double, rank + 4> I4s() const;
xt::xtensor<double, rank + 4> I4d() const;
// Check that a type has been set everywhere (throws if unset points are found)
void check() const;
// Set parameters for a batch of points
// (uniform for all points specified: that have "I(i,j) == 1")
void setElastic(
const xt::xtensor<size_t, rank>& I,
double K,
double G);
void setCusp(
const xt::xtensor<size_t, rank>& I,
double K,
double G,
const xt::xtensor<double, 1>& epsy,
bool init_elastic = true);
void setSmooth(
const xt::xtensor<size_t, rank>& I,
double K,
double G,
const xt::xtensor<double, 1>& epsy,
bool init_elastic = true);
// Set parameters for a batch of points:
// each to the same material, but with different parameters:
// the matrix "idx" refers to a which entry to use: "K[idx]", "G[idx]", or "epsy[idx,:]"
void setElastic(
const xt::xtensor<size_t, rank>& I,
const xt::xtensor<size_t, rank>& idx,
const xt::xtensor<double, 1>& K,
const xt::xtensor<double, 1>& G);
void setCusp(
const xt::xtensor<size_t, rank>& I,
const xt::xtensor<size_t, rank>& idx,
const xt::xtensor<double, 1>& K,
const xt::xtensor<double, 1>& G,
const xt::xtensor<double, 2>& epsy,
bool init_elastic = true);
void setSmooth(
const xt::xtensor<size_t, rank>& I,
const xt::xtensor<size_t, rank>& idx,
const xt::xtensor<double, 1>& K,
const xt::xtensor<double, 1>& G,
const xt::xtensor<double, 2>& epsy,
bool init_elastic = true);
// Set strain tensor, get the response
void setStrain(const xt::xtensor<double, rank + 2>& Eps);
void stress(xt::xtensor<double, rank + 2>& Sig) const;
void tangent(xt::xtensor<double, rank + 4>& C) const;
void currentIndex(xt::xtensor<size_t, rank>& arg) const;
void currentYieldLeft(xt::xtensor<double, rank>& arg) const;
void currentYieldRight(xt::xtensor<double, rank>& arg) const;
bool checkYieldBoundLeft(size_t n = 0) const;
bool checkYieldBoundRight(size_t n = 0) const;
void epsp(xt::xtensor<double, rank>& arg) const;
void energy(xt::xtensor<double, rank>& arg) const;
// Auto-allocation of the functions above
xt::xtensor<double, rank + 2> Stress() const;
xt::xtensor<double, rank + 4> Tangent() const;
xt::xtensor<size_t, rank> CurrentIndex() const;
xt::xtensor<double, rank> CurrentYieldLeft() const;
xt::xtensor<double, rank> CurrentYieldRight() const;
xt::xtensor<double, rank> Epsp() const;
xt::xtensor<double, rank> Energy() const;
private:
// Material vectors
std::vector<Elastic> m_Elastic;
std::vector<Cusp> m_Cusp;
std::vector<Smooth> m_Smooth;
// Identifiers for each matrix entry
xt::xtensor<size_t, rank> m_type; // type (e.g. "Type::Elastic")
xt::xtensor<size_t, rank> m_index; // index from the relevant material vector (e.g. "m_Elastic")
// Shape
static const size_t m_ndim = 2;
size_t m_size;
std::array<size_t, rank> m_shape;
std::array<size_t, rank + 2> m_shape_tensor2;
std::array<size_t, rank + 4> m_shape_tensor4;
// Internal check
bool m_allSet = false; // true if all points have a material assigned
void checkAllSet(); // check if all points have a material assigned (modifies "m_allSet")
};
} // namespace Cartesian2d
} // namespace GMatElastoPlasticQPot
#include "Cartesian2d.hpp"
#include "Cartesian2d_Array.hpp"
#include "Cartesian2d_Cusp.hpp"
#include "Cartesian2d_Elastic.hpp"
#include "Cartesian2d_Smooth.hpp"
#endif