-
Notifications
You must be signed in to change notification settings - Fork 1
/
friction-optimized.cpp
447 lines (334 loc) · 13.3 KB
/
friction-optimized.cpp
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
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
#include <GMatElastoPlasticQPot/Cartesian2d.h>
#include <GooseFEM/GooseFEM.h>
#include <GooseFEM/Matrix.h>
#include <xtensor/xtensor.hpp>
#include <xtensor/xrandom.hpp>
#include <xtensor/xnorm.hpp>
#include <xtensor/xio.hpp>
#include <xtensor/xcsv.hpp>
#include <fstream>
namespace GF = GooseFEM;
namespace QD = GooseFEM::Element::Quad4;
namespace GM = GMatElastoPlasticQPot::Cartesian2d;
// #define CHECK_MYSHORTCUT
class System {
private:
// mesh parameters
xt::xtensor<size_t, 2> m_conn;
xt::xtensor<size_t, 2> m_conn_elastic;
xt::xtensor<size_t, 2> m_conn_plastic;
xt::xtensor<double, 2> m_coor;
xt::xtensor<size_t, 2> m_dofs;
xt::xtensor<size_t, 1> m_iip;
// mesh dimensions
size_t m_N; // == nelem_plas
size_t m_nelem;
size_t m_nelem_elas;
size_t m_nelem_plas;
size_t m_nne;
size_t m_nnode;
size_t m_ndim;
size_t m_nip;
// element sets
xt::xtensor<size_t, 1> m_elastic;
xt::xtensor<size_t, 1> m_plastic;
// numerical quadrature
QD::Quadrature m_quad;
QD::Quadrature m_quad_elas;
QD::Quadrature m_quad_plas;
// convert vectors between 'nodevec', 'elemvec', ...
GF::VectorPartitioned m_vector;
GF::VectorPartitioned m_vector_elas;
GF::VectorPartitioned m_vector_plas;
// mass matrix
GF::MatrixDiagonalPartitioned m_M;
// damping matrix
GF::MatrixDiagonal m_D;
// material definition
GM::Array<2> m_material;
GM::Array<2> m_material_elas;
GM::Array<2> m_material_plas;
// convergence check
GF::Iterate::StopList m_stop = GF::Iterate::StopList(20);
// time evolution
double m_t = 0.0; // current time
double m_dt; // time step
// nodal displacements, velocities, and accelerations (current and last time-step)
xt::xtensor<double, 2> m_u;
xt::xtensor<double, 2> m_v;
xt::xtensor<double, 2> m_a;
xt::xtensor<double, 2> m_v_n;
xt::xtensor<double, 2> m_a_n;
// element vectors
xt::xtensor<double, 3> m_ue;
xt::xtensor<double, 3> m_fe;
xt::xtensor<double, 3> m_ue_plas;
xt::xtensor<double, 3> m_fe_plas;
// nodal forces
#ifdef CHECK_MYSHORTCUT
xt::xtensor<double, 2> m_fmaterial;
#endif
xt::xtensor<double, 2> m_felas;
xt::xtensor<double, 2> m_fplas;
xt::xtensor<double, 2> m_fdamp;
xt::xtensor<double, 2> m_fint;
xt::xtensor<double, 2> m_fext;
xt::xtensor<double, 2> m_fres;
// integration point tensors
xt::xtensor<double, 4> m_Eps;
xt::xtensor<double, 4> m_Eps_elas;
xt::xtensor<double, 4> m_Eps_plas;
xt::xtensor<double, 4> m_Sig;
xt::xtensor<double, 4> m_Sig_elas;
xt::xtensor<double, 4> m_Sig_plas;
// stiffness matrix
GF::Matrix m_K_elas;
public:
System(size_t N) : m_N(N)
{
// ----
// mesh
// ----
double h = xt::numeric_constants<double>::PI;
double L = h * static_cast<double>(m_N);
GF::Mesh::Quad4::FineLayer mesh(m_N, m_N, h);
m_coor = mesh.coor();
m_conn = mesh.conn();
m_dofs = mesh.dofs();
m_nnode = m_coor.shape(0);
m_ndim = m_coor.shape(1);
m_nelem = m_conn.shape(0);
m_nne = m_conn.shape(1);
m_plastic = mesh.elementsMiddleLayer();
m_elastic = xt::setdiff1d(xt::arange(m_nelem), m_plastic);
m_nelem_elas = m_elastic.size();
m_nelem_plas = m_plastic.size();
m_conn_elastic = xt::view(m_conn, xt::keep(m_elastic), xt::all());
m_conn_plastic = xt::view(m_conn, xt::keep(m_plastic), xt::all());
auto left = mesh.nodesLeftOpenEdge();
auto right = mesh.nodesRightOpenEdge();
xt::view(m_dofs, xt::keep(right), 0) = xt::view(m_dofs, xt::keep(left), 0);
xt::view(m_dofs, xt::keep(right), 1) = xt::view(m_dofs, xt::keep(left), 1);
auto top = mesh.nodesTopEdge();
auto bottom = mesh.nodesBottomEdge();
size_t nfix = top.size();
m_iip = xt::empty<decltype(m_iip)::value_type>({2 * m_ndim * nfix});
xt::view(m_iip, xt::range(0 * nfix, 1 * nfix)) = xt::view(m_dofs, xt::keep(bottom), 0);
xt::view(m_iip, xt::range(1 * nfix, 2 * nfix)) = xt::view(m_dofs, xt::keep(bottom), 1);
xt::view(m_iip, xt::range(2 * nfix, 3 * nfix)) = xt::view(m_dofs, xt::keep(top), 0);
xt::view(m_iip, xt::range(3 * nfix, 4 * nfix)) = xt::view(m_dofs, xt::keep(top), 1);
m_vector = GF::VectorPartitioned(m_conn, m_dofs, m_iip);
m_vector_elas = GF::VectorPartitioned(m_conn_elastic, m_dofs, m_iip);
m_vector_plas = GF::VectorPartitioned(m_conn_plastic, m_dofs, m_iip);
m_quad = QD::Quadrature(m_vector.AsElement(m_coor));
m_quad_elas = QD::Quadrature(m_vector_elas.AsElement(m_coor));
m_quad_plas = QD::Quadrature(m_vector_plas.AsElement(m_coor));
m_nip = m_quad.nip();
m_u = m_vector.AllocateNodevec(0.0);
m_v = m_vector.AllocateNodevec(0.0);
m_a = m_vector.AllocateNodevec(0.0);
m_v_n = m_vector.AllocateNodevec(0.0);
m_a_n = m_vector.AllocateNodevec(0.0);
m_ue = m_vector.AllocateElemvec(0.0);
m_fe = m_vector.AllocateElemvec(0.0);
m_ue_plas = m_vector_plas.AllocateElemvec(0.0);;
m_fe_plas = m_vector_plas.AllocateElemvec(0.0);;
#ifdef CHECK_MYSHORTCUT
m_fmaterial = m_vector.AllocateNodevec(0.0);
#endif
m_felas = m_vector.AllocateNodevec(0.0);
m_fplas = m_vector.AllocateNodevec(0.0);
m_fdamp = m_vector.AllocateNodevec(0.0);
m_fint = m_vector.AllocateNodevec(0.0);
m_fext = m_vector.AllocateNodevec(0.0);
m_fres = m_vector.AllocateNodevec(0.0);
m_Eps = m_quad.AllocateQtensor<2>(0.0);
m_Sig = m_quad.AllocateQtensor<2>(0.0);
m_Eps_elas = m_quad_elas.AllocateQtensor<2>(0.0);
m_Sig_elas = m_quad_elas.AllocateQtensor<2>(0.0);
m_Eps_plas = m_quad_plas.AllocateQtensor<2>(0.0);
m_Sig_plas = m_quad_plas.AllocateQtensor<2>(0.0);
// --------
// material
// --------
// material parameters
double c = 1.0;
double G = 1.0;
double K = 10.0 * G;
double rho = G / std::pow(c, 2.0);
double qL = 2.0 * xt::numeric_constants<double>::PI / L;
double qh = 2.0 * xt::numeric_constants<double>::PI / h;
double alpha = std::sqrt(2.0) * qL * c * rho;
m_dt = 1.0 / (c * qh) / 10.0;
// material definition
m_material = GM::Array<2>({m_nelem, m_nip});
m_material_elas = GM::Array<2>({m_nelem_elas, m_nip});
m_material_plas = GM::Array<2>({m_nelem_plas, m_nip});
// assign elastic material points
{
xt::xtensor<size_t, 2> I = xt::zeros<size_t>({m_nelem, m_nip});
xt::view(I, xt::keep(m_elastic), xt::all()) = 1ul;
m_material.setElastic(I, K, G);
I = xt::ones<size_t>({m_nelem_elas, m_nip});
m_material_elas.setElastic(I, K, G);
}
// assign plastic material points
{
double k = 2.0;
xt::xtensor<double, 2> epsy = 1e-5 + 1e-3 * xt::random::weibull<double>(std::array<size_t, 2>{m_N, 1000}, k, 1.0);
xt::view(epsy, xt::all(), 0) = 1e-5 + 1e-3 * xt::random::rand<double>({m_N});
epsy = xt::cumsum(epsy, 1);
xt::xtensor<size_t, 2> I = xt::zeros<size_t>({m_nelem, m_nip});
xt::xtensor<size_t, 2> idx = xt::zeros<size_t>({m_nelem, m_nip});
xt::view(I, xt::keep(m_plastic), xt::all()) = 1ul;
xt::view(idx, xt::keep(m_plastic), xt::all()) = xt::arange<size_t>(m_N).reshape({-1, 1});
xt::xtensor<double, 1> unit = xt::ones<double>({m_N});
m_material.setCusp(I, idx, K * unit, G * unit, epsy);
I = xt::ones<size_t>({m_N, m_nip});
idx = xt::zeros<size_t>({m_N, m_nip});
xt::view(idx, xt::all(), xt::all()) = xt::arange<size_t>(m_N).reshape({-1, 1});
m_material_plas.setCusp(I, idx, K * unit, G * unit, epsy);
}
m_material.check();
m_material_elas.check();
m_material_plas.check();
m_material.setStrain(m_Eps);
m_material_elas.setStrain(m_Eps_elas);
m_material_plas.setStrain(m_Eps_plas);
// -----------
// mass matrix
// -----------
{
auto x = m_vector.AsElement(m_coor);
QD::Quadrature nodalQuad(x, QD::Nodal::xi(), QD::Nodal::w());
xt::xtensor<double, 2> val_quad = rho * xt::ones<double>({m_nelem, nodalQuad.nip()});
m_M = GF::MatrixDiagonalPartitioned(m_conn, m_dofs, m_iip);
m_M.assemble(nodalQuad.Int_N_scalar_NT_dV(val_quad));
}
// --------------
// damping matrix
// --------------
{
auto x = m_vector.AsElement(m_coor);
QD::Quadrature nodalQuad(x, QD::Nodal::xi(), QD::Nodal::w());
xt::xtensor<double, 2> val_quad = alpha * xt::ones<double>({m_nelem, nodalQuad.nip()});
m_D = GF::MatrixDiagonal(m_conn, m_dofs);
m_D.assemble(nodalQuad.Int_N_scalar_NT_dV(val_quad));
}
// ------------------------
// partial stiffness matrix
// ------------------------
m_K_elas = GF::Matrix(m_conn_elastic, m_dofs);
m_K_elas.assemble(m_quad_elas.Int_gradN_dot_tensor4_dot_gradNT_dV(m_material_elas.Tangent()));
}
public:
void computeStrainStress()
{
m_vector.asElement(m_u, m_ue);
m_quad.symGradN_vector(m_ue, m_Eps);
m_material.setStrain(m_Eps);
m_material.stress(m_Sig);
}
public:
void computeStrainStressForcesWeakLayer()
{
m_vector_plas.asElement(m_u, m_ue_plas);
m_quad_plas.symGradN_vector(m_ue_plas, m_Eps_plas);
m_material_plas.setStrain(m_Eps_plas);
m_material_plas.stress(m_Sig_plas);
m_quad_plas.int_gradN_dot_tensor2_dV(m_Sig_plas, m_fe_plas);
m_vector_plas.assembleNode(m_fe_plas, m_fplas);
m_K_elas.dot(m_u, m_felas);
}
public:
void timeStep()
{
// history
m_t += m_dt;
xt::noalias(m_v_n) = m_v;
xt::noalias(m_a_n) = m_a;
// new displacement
xt::noalias(m_u) = m_u + m_dt * m_v + 0.5 * std::pow(m_dt, 2.) * m_a;
// compute strain/strain, and corresponding force
computeStrainStressForcesWeakLayer();
#ifdef CHECK_MYSHORTCUT
computeStrainStress();
m_quad.int_gradN_dot_tensor2_dV(m_Sig, m_fe);
m_vector.assembleNode(m_fe, m_fmaterial);
if (!xt::allclose(m_fmaterial, m_felas + m_fplas)) {
throw std::runtime_error("Forces");
}
#endif
// estimate new velocity, update corresponding force
xt::noalias(m_v) = m_v_n + m_dt * m_a_n;
m_D.dot(m_v, m_fdamp);
// compute residual force & solve
xt::noalias(m_fint) = m_felas + m_fplas + m_fdamp;
m_vector.copy_p(m_fint, m_fext);
xt::noalias(m_fres) = m_fext - m_fint;
m_M.solve(m_fres, m_a);
// re-estimate new velocity, update corresponding force
xt::noalias(m_v) = m_v_n + .5 * m_dt * (m_a_n + m_a);
m_D.dot(m_v, m_fdamp);
// compute residual force & solve
xt::noalias(m_fint) = m_felas + m_fplas + m_fdamp;
m_vector.copy_p(m_fint, m_fext);
xt::noalias(m_fres) = m_fext - m_fint;
m_M.solve(m_fres, m_a);
// new velocity, update corresponding force
xt::noalias(m_v) = m_v_n + .5 * m_dt * (m_a_n + m_a);
m_D.dot(m_v, m_fdamp);
// compute residual force & solve
xt::noalias(m_fint) = m_felas + m_fplas + m_fdamp;
m_vector.copy_p(m_fint, m_fext);
xt::noalias(m_fres) = m_fext - m_fint;
m_M.solve(m_fres, m_a);
}
public:
xt::xtensor<double, 2> run()
{
xt::xtensor<double, 3> dF = xt::zeros<double>({1001, 2, 2});
xt::view(dF, xt::range(1, dF.shape(0)), 0, 1) = 0.004 / 1000.0;
xt::xtensor<double, 2> ret = xt::zeros<double>(std::array<size_t, 2>{dF.shape(0), 2});
auto dV = m_quad.AsTensor<2>(m_quad.dV());
for (size_t inc = 0 ; inc < dF.shape(0); ++inc) {
for (size_t i = 0; i < m_nnode; ++i) {
for (size_t j = 0; j < m_ndim; ++j) {
for (size_t k = 0; k < m_ndim; ++k) {
m_u(i, j) += dF(inc, j, k) * (m_coor(i, k) - m_coor(0, k));
}
}
}
computeStrainStressForcesWeakLayer();
for (size_t iiter = 0; iiter < 99999 ; ++iiter) {
timeStep();
if (m_stop.stop(xt::norm_l2(m_fres)() / xt::norm_l2(m_fext)(), 1e-5)) {
std::cout << inc << ", " << iiter << std::endl;
break;
}
}
m_v.fill(0.0);
m_a.fill(0.0);
m_stop.reset();
computeStrainStress();
xt::xtensor<double, 2> Epsbar = xt::average(m_Eps, dV, {0, 1});
xt::xtensor<double, 2> Sigbar = xt::average(m_Sig, dV, {0, 1});
ret(inc, 0) = GM::Epsd(Epsbar)();
ret(inc, 1) = GM::Epsd(Sigbar)();
}
return ret;
}
};
int main(void)
{
size_t N = std::pow(3, 4);
System sys(N);
auto ret = sys.run();
std::ofstream outfile("friction.txt");
xt::dump_csv(outfile, ret);
outfile.close();
#ifdef CHECK_MYSHORTCUT
std::cout << "Shortcut successfully checked" << std::endl;
#endif
return 0;
}