-
Notifications
You must be signed in to change notification settings - Fork 36
/
nonlinear_elasticity_global_rom.cpp
2727 lines (2290 loc) · 89.6 KB
/
nonlinear_elasticity_global_rom.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
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
// libROM MFEM Example: parametric ROM for nonlinear elasticity problem (adapted from ex10p.cpp)
//
// Compile with: make nonlinear_elasticity_global_rom
//
// Description: This examples solves a time dependent nonlinear elasticity
// problem of the form dv/dt = H(x) + S v, dx/dt = v, where H is a
// hyperelastic model and S is a viscosity operator of Laplacian
// type. The geometry of the domain is assumed to be as follows:
//
// +---------------------+
// boundary --->| |
// attribute 1 | |
// (fixed) +---------------------+
//
// The example demonstrates the use of hyper reduction to solve a
// nonlinear elasticity problem. Time integration is done with various
// explicit time integrator solvers. The basis for the velocity field
// is either constructed using a separate velocity basis or using the
// displacement basis. It is possible to set the initial condition in
// terms of either velocity or deformation. The velocity initial condition
// works better when both velocity and displacement bases are used. The
// deformation initial condition is better when only the displacement
// basis is used. The input flag -sc controls the scaling of the initial
// condition applied. This is what parameterizes the ROM. If the scaling
// factor is within the range +-10%, the results are generally accurate.
// =================================================================================
//
// Sample runs and results for parametric ROM using displacement basis, velocity basis
// and nonlinear term basis, with velocity initial condition:
//
// Offline phase:
// ./nonlinear_elasticity_global_rom -offline -dt 0.01 -tf 5.0 -s 14 -vs 100 -sc 0.9 -id 0
// ./nonlinear_elasticity_global_rom -offline -dt 0.01 -tf 5.0 -s 14 -vs 100 -sc 1.1 -id 1
//
// Merge phase:
// ./nonlinear_elasticity_global_rom -merge -ns 2 -dt 0.01 -tf 5.0
//
// Create FOM comparison data:
// ./nonlinear_elasticity_global_rom -offline -dt 0.01 -tf 5.0 -s 14 -vs 100 -sc 1.0 -id 2
// Output message:
// Elapsed time for time integration loop 9.17153
// Elapsed time for entire simulation 10.0001
//
// Online phase with full sampling:
// ./nonlinear_elasticity_global_rom -online -dt 0.01 -tf 5.0 -s 14 -vs 100 -hyp -hrtype gnat -rvdim 40 -rxdim 10 -hdim 71 -nsr 1170 -sc 1.0
// Output message:
// Elapsed time for time integration loop 7.03877
// Relative error of ROM position (x) at t_final: 5 is 0.000231698
// Relative error of ROM velocity (v) at t_final: 5 is 0.466941
// Elapsed time for entire simulation 8.98334
//
// Online phase with strong hyper-reduction, using GNAT (over-sampled DEIM):
// ./nonlinear_elasticity_global_rom -online -dt 0.01 -tf 5.0 -s 14 -vs 100 -hyp -hrtype gnat -rvdim 40 -rxdim 10 -hdim 71 -nsr 100 -sc 1.0
// Output message:
// Elapsed time for time integration loop 4.18114
// Relative error of ROM position (x) at t_final: 5 is 0.00209877
// Relative error of ROM velocity (v) at t_final: 5 is 1.39472
// Elapsed time for entire simulation 4.52802
//
// Online phase with strong hyper-reduction, using QDEIM:
// ./nonlinear_elasticity_global_rom -online -dt 0.01 -tf 5.0 -s 14 -vs 100 -hyp -hrtype qdeim -rvdim 40 -rxdim 10 -hdim 71 -nsr 100 -sc 1.0
// Output message:
// Elapsed time for time integration loop 4.2972
// Relative error of ROM position (x) at t_final: 5 is 0.00188458
// Relative error of ROM velocity (v) at t_final: 5 is 0.978726
// Elapsed time for entire simulation 5.33154
//
// Online phase with EQP hyper-reduction
// ./nonlinear_elasticity_global_rom -online -dt 0.01 -tf 5.0 -s 14 -vs 100 -eqp -ns 2 -ntw 50 -rvdim 40 -rxdim 10 -hdim 1 -sc 1.0
// Output message:
// Elapsed time for time integration loop 2.70006
// Relative error of ROM position (x) at t_final: 5 is 0.00169666
// Relative error of ROM velocity (v) at t_final: 5 is 17.71
// Elapsed time for entire simulation 831.216
//
// Online phase with EQP hyper-reduction and subsampling of snapshots
// ./nonlinear_elasticity_global_rom -online -dt 0.01 -tf 5.0 -s 14 -vs 100 -eqp -ns 2 -ntw 50 -rvdim 40 -rxdim 10 -hdim 1 -sc 1.0 -sbs 10
// Output message:
// Elapsed time for time integration loop 0.552546
// Relative error of ROM position (x) at t_final: 5 is 0.00247479
// Relative error of ROM velocity (v) at t_final: 5 is 1.33349
// Elapsed time for entire simulation 11.469
//
// =================================================================================
//
// Sample runs and results for parametric ROM using only displacement basis
// and nonlinear term basis:
// Offline phase:
// ./nonlinear_elasticity_global_rom -offline -dt 0.01 -tf 5.0 -s 14 -vs 100 -sc 0.9 -xbo -def-ic -id 0
// ./nonlinear_elasticity_global_rom -offline -dt 0.01 -tf 5.0 -s 14 -vs 100 -sc 1.1 -xbo -def-ic -id 1
//
// Merge phase:
// ./nonlinear_elasticity_global_rom -merge -ns 2 -dt 0.01 -tf 5.0 -xbo
//
// Create FOM comparison data:
// ./nonlinear_elasticity_global_rom -offline -dt 0.01 -tf 5.0 -s 14 -vs 100 -sc 1.0 -xbo -def-ic -id 2
// Output message:
// Elapsed time for time integration loop 9.61062
// Elapsed time for entire simulation 10.3806
//
// Online phase with full sampling:
// ./nonlinear_elasticity_global_rom -online -dt 0.01 -tf 5.0 -s 14 -vs 100 -hyp -hrtype gnat -rxdim 57 -hdim 183 -nsr 1170 -sc 1.0 -xbo -def-ic
// Output message:
// Elapsed time for time integration loop 8.43858
// Relative error of ROM position (x) at t_final: 5 is 7.08272e-05
// Relative error of ROM velocity (v) at t_final: 5 is 0.00387647
// Elapsed time for entire simulation 24.7245
//
// Online phase with strong hyper reduction, using GNAT (over-sampled DEIM):
// ./nonlinear_elasticity_global_rom -online -dt 0.01 -tf 5.0 -s 14 -vs 100 -hyp -hrtype gnat -rxdim 2 -hdim 4 -nsr 10 -sc 1.0 -xbo -def-ic
// Output message:
// Elapsed time for time integration loop 0.507339
// Relative error of ROM position (x) at t_final: 5 is 0.0130818
// Relative error of ROM velocity (v) at t_final: 5 is 0.979978
// Elapsed time for entire simulation 0.583618
//
// Online phase with strong hyper reduction, using QDEIM:
// ./nonlinear_elasticity_global_rom -online -dt 0.01 -tf 5.0 -s 14 -vs 100 -hyp -hrtype qdeim -rxdim 2 -hdim 4 -nsr 10 -sc 1.0 -xbo -def-ic
// Output message:
// Elapsed time for time integration loop 0.504657
// Relative error of ROM position (x) at t_final: 5 is 0.00883079
// Relative error of ROM velocity (v) at t_final: 5 is 1.26721
// Elapsed time for entire simulation 0.581446
//
// Online phase with EQP hyper reduction:
// ./nonlinear_elasticity_global_rom -online -dt 0.01 -tf 5.0 -s 14 -vs 100 -eqp -ns 2 -rxdim 2 -hdim 1 -ntw 25 -sc 1.00 -xbo -def-ic
// Elapsed time for time integration loop 0.144368
// Relative error of ROM position (x) at t_final: 5 is 0.0189361
// Relative error of ROM velocity (v) at t_final: 5 is 0.830724
// Elapsed time for entire simulation 4.97996
//
// Online phase with EQP hyper reduction and subsampling of snapshots:
// ./nonlinear_elasticity_global_rom -online -dt 0.01 -tf 5.0 -s 14 -vs 100 -eqp -ns 2 -rxdim 2 -hdim 1 -ntw 25 -sc 1.00 -xbo -def-ic -sbs 10
// Elapsed time for time integration loop 0.0141128
// Relative error of ROM position (x) at t_final: 5 is 0.0243952
// Relative error of ROM velocity (v) at t_final: 5 is 0.986434
// Elapsed time for entire simulation 0.782108
//
// This example runs in parallel with MPI, by using the same number of MPI ranks
// in all phases (offline, merge, online).
#include "mfem.hpp"
#include "linalg/Vector.h"
#include "linalg/BasisGenerator.h"
#include "linalg/BasisReader.h"
#include "linalg/NNLS.h"
#include "hyperreduction/Hyperreduction.h"
#include "mfem/SampleMesh.hpp"
#include "mfem/Utilities.hpp"
#include <memory>
#include <cmath>
#include <limits>
using namespace std;
using namespace mfem;
#include "nonlinear_elasticity_global_rom_eqp.hpp"
class RomOperator : public TimeDependentOperator
{
private:
int rxdim, rvdim, hdim;
int nsamp_H;
double current_dt;
bool oversampling;
CAROM::Matrix *V_v_sp, *V_x_sp;
CAROM::Vector *V_xTv_0_sp;
CAROM::Matrix *V_xTV_v_sp;
CAROM::Matrix *V_v_sp_dist;
CAROM::Vector *psp_librom, *psp_v_librom;
CAROM::Vector *Vx_librom_temp;
Vector *psp;
Vector *psp_x;
Vector *psp_v;
Vector *Vx_temp;
mutable Vector zH;
mutable CAROM::Vector zX;
mutable Vector zX_MFEM;
mutable CAROM::Vector zN;
const CAROM::Matrix *Hsinv;
mutable CAROM::Vector *z_librom;
mutable Vector z;
mutable Vector z_x;
mutable Vector z_v;
CAROM::Vector *v0_fom_librom, *v0_librom, *V_xTv_0;
CAROM::Matrix *V_xTV_v;
bool hyperreduce;
bool x_base_only;
CAROM::Vector *pfom_librom, *pfom_v_librom;
Vector *pfom;
Vector *pfom_x;
Vector *pfom_v;
mutable Vector *zfom_x;
mutable Vector *zfom_v;
CAROM::Vector *zfom_x_librom;
CAROM::SampleMeshManager *smm;
CAROM::Vector *z_v_librom;
CAROM::Vector *z_x_librom;
// Data for EQP
bool eqp;
const IntegrationRule *ir_eqp;
std::vector<double> eqp_rw;
std::vector<int> eqp_qp;
Vector eqp_coef;
Vector eqp_DS_coef;
const bool fastIntegration = true;
ElemMatrices *em;
CAROM::Matrix eqp_lifting;
std::vector<int> eqp_liftDOFs;
mutable CAROM::Vector eqp_lifted;
int rank;
NeoHookeanModel *model;
protected:
CAROM::Matrix *S_hat;
CAROM::Vector *S_hat_v0;
Vector *S_hat_v0_temp;
CAROM::Vector *S_hat_v0_temp_librom;
CAROM::Matrix *M_hat;
CAROM::Matrix *M_hat_inv;
const CAROM::Matrix *U_H;
HyperelasticOperator *fomSp;
CGSolver M_hat_solver; // Krylov solver for inverting the reduced mass matrix M_hat
HypreSmoother M_hat_prec; // Preconditioner for the reduced mass matrix M_hat
public:
HyperelasticOperator *fom;
RomOperator(HyperelasticOperator *fom_,
HyperelasticOperator *fomSp_, const int rvdim_, const int rxdim_,
const int hdim_, CAROM::SampleMeshManager *smm_, const Vector *v0_,
const Vector *x0_, const Vector v0_fom_, const CAROM::Matrix *V_v_,
const CAROM::Matrix *V_x_, const CAROM::Matrix *U_H_,
const CAROM::Matrix *Hsinv_, const int myid, const bool oversampling_,
const bool hyperreduce_, const bool x_base_only_, const bool use_eqp,
CAROM::Vector *eqpSol,
const IntegrationRule *ir_eqp_, NeoHookeanModel *model_);
virtual void Mult(const Vector &y, Vector &dy_dt) const;
void Mult_Hyperreduced(const Vector &y, Vector &dy_dt) const;
void Mult_FullOrder(const Vector &y, Vector &dy_dt) const;
void SetEQP(CAROM::Vector *eqpSol);
CAROM::Matrix V_v, V_x, V_vTU_H;
const Vector *x0;
const Vector *v0;
Vector v0_fom;
virtual ~RomOperator();
};
/** Function representing the elastic energy density for the given hyperelastic
model+deformation. Used in HyperelasticOperator::GetElasticEnergyDensity. */
class ElasticEnergyCoefficient : public Coefficient
{
private:
HyperelasticModel &model;
const ParGridFunction &x;
DenseMatrix J;
public:
ElasticEnergyCoefficient(HyperelasticModel &m, const ParGridFunction &x_)
: model(m), x(x_) {}
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip);
virtual ~ElasticEnergyCoefficient() {}
};
void InitialDeformationIC1(const Vector &x, Vector &y);
void InitialVelocityIC1(const Vector &x, Vector &v);
void InitialDeformationIC2(const Vector &x, Vector &y);
void InitialVelocityIC2(const Vector &x, Vector &v);
void visualize(ostream &out, ParMesh *mesh, ParGridFunction *deformed_nodes,
ParGridFunction *field, const char *field_name = NULL,
bool init_vis = false);
// TODO: move this to the library?
CAROM::Matrix *GetFirstColumns(const int N, const CAROM::Matrix *A)
{
CAROM::Matrix *S = new CAROM::Matrix(A->numRows(), std::min(N, A->numColumns()),
A->distributed());
for (int i = 0; i < S->numRows(); ++i)
{
for (int j = 0; j < S->numColumns(); ++j)
(*S)(i, j) = (*A)(i, j);
}
// delete A; // TODO: find a good solution for this.
return S;
}
// TODO: move this to the library?
void MergeBasis(const int dimFOM, const int nparam, const int max_num_snapshots,
std::string name)
{
MFEM_VERIFY(nparam > 0, "Must specify a positive number of parameter sets");
bool update_right_SV = false;
bool isIncremental = false;
CAROM::Options options(dimFOM, nparam * max_num_snapshots, 1, update_right_SV);
CAROM::BasisGenerator generator(options, isIncremental, "basis" + name);
for (int paramID = 0; paramID < nparam; ++paramID)
{
std::string snapshot_filename = "basis" + std::to_string(
paramID) + "_" + name + "_snapshot";
generator.loadSamples(snapshot_filename, "snapshot");
}
generator.endSamples(); // save the merged basis file
int cutoff = 0;
generator.finalSummary(1e-7, cutoff, "mergedSV_" + name + ".txt");
}
const CAROM::Matrix *GetSnapshotMatrix(const int dimFOM, const int nparam,
const int max_num_snapshots, std::string name)
{
MFEM_VERIFY(nparam > 0, "Must specify a positive number of parameter sets");
bool update_right_SV = false;
bool isIncremental = false;
CAROM::Options options(dimFOM, nparam * max_num_snapshots, 1, update_right_SV);
CAROM::BasisGenerator generator(options, isIncremental, "basis" + name);
for (int paramID = 0; paramID < nparam; ++paramID)
{
std::string snapshot_filename = "basis" + std::to_string(
paramID) + "_" + name + "_snapshot";
generator.loadSamples(snapshot_filename, "snapshot");
}
// TODO: this deep copy is inefficient, just to get around generator owning the matrix.
CAROM::Matrix *s = new CAROM::Matrix(*generator.getSnapshotMatrix());
return s;
// return generator.getSnapshotMatrix(); // BUG: the matrix is deleted when generator goes out of scope.
}
// TODO: remove this by making online computation serial?
void BroadcastUndistributedRomVector(CAROM::Vector *v)
{
const int N = v->dim();
MFEM_VERIFY(N > 0, "");
double *d = new double[N];
MFEM_VERIFY(d != 0, "");
for (int i = 0; i < N; ++i)
d[i] = (*v)(i);
MPI_Bcast(d, N, MPI_DOUBLE, 0, MPI_COMM_WORLD);
for (int i = 0; i < N; ++i)
(*v)(i) = d[i];
delete[] d;
}
// Scaling factor for parameterization
double s = 1.0;
int main(int argc, char *argv[])
{
// 1. Initialize MPI.
int num_procs, myid;
MPI_Init(&argc, &argv);
MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
MPI_Comm_rank(MPI_COMM_WORLD, &myid);
// 2. Parse command-line options.
const char *mesh_file = "../data/beam-quad.mesh";
int ser_ref_levels = 2;
int par_ref_levels = 0;
int order = 2;
int ode_solver_type = 14; // RK4
int vis_steps = 1;
double t_final = 15.0;
double dt = 0.03;
double visc = 1e-2;
double mu = 0.25;
double K = 5.0;
bool adaptive_lin_rtol = true;
bool visualization = true;
bool visit = false;
bool def_ic = false;
// ROM parameters
bool offline = false;
bool merge = false;
bool online = false;
bool use_eqp = false;
bool writeSampleMesh = false;
bool hyperreduce = true;
bool x_base_only = false;
int num_samples_req = -1;
const char *samplingType = "gnat";
int nsets = 0;
int id_param = 0;
// Number of basis vectors to use
int rxdim = -1;
int rvdim = -1;
int hdim = -1;
bool preconditionNNLS = false;
double tolNNLS = 1.0e-14;
int maxNNLSnnz = 0;
// Number of time windows
int n_windows = 0;
// Subsampling steps
int snap_step = 1;
OptionsParser args(argc, argv);
args.AddOption(&mesh_file, "-m", "--mesh",
"Mesh file to use.");
args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
"Number of times to refine the mesh uniformly in serial.");
args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
"Number of times to refine the mesh uniformly in parallel.");
args.AddOption(&order, "-o", "--order",
"Order (degree) of the finite elements.");
args.AddOption(&ode_solver_type, "-s", "--ode-solver",
"ODE solver: 1 - Backward Euler, 2 - SDIRK2, 3 - SDIRK3,\n\t"
" 11 - Forward Euler, 12 - RK2,\n\t"
" 13 - RK3 SSP, 14 - RK4."
" 22 - Implicit Midpoint Method,\n\t"
" 23 - SDIRK23 (A-stable), 24 - SDIRK34");
args.AddOption(&t_final, "-tf", "--t-final",
"Final time; start time is 0.");
args.AddOption(&dt, "-dt", "--time-step",
"Time step.");
args.AddOption(&visc, "-v", "--viscosity",
"Viscosity coefficient.");
args.AddOption(&mu, "-mu", "--shear-modulus",
"Shear modulus in the Neo-Hookean hyperelastic model.");
args.AddOption(&K, "-K", "--bulk-modulus",
"Bulk modulus in the Neo-Hookean hyperelastic model.");
args.AddOption(&adaptive_lin_rtol, "-alrtol", "--adaptive-lin-rtol",
"-no-alrtol", "--no-adaptive-lin-rtol",
"Enable or disable adaptive linear solver rtol.");
args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
"--no-visualization",
"Enable or disable GLVis visualization.");
args.AddOption(&visit, "-visit", "--visit-datafiles", "-no-visit",
"--no-visit-datafiles",
"Save data files for VisIt (visit.llnl.gov) visualization.");
args.AddOption(&vis_steps, "-vs", "--visualization-steps",
"Visualize every n-th timestep.");
args.AddOption(&nsets, "-ns", "--nset", "Number of parametric snapshot sets");
args.AddOption(&offline, "-offline", "--offline", "-no-offline", "--no-offline",
"Enable or disable the offline phase.");
args.AddOption(&online, "-online", "--online", "-no-online", "--no-online",
"Enable or disable the online phase.");
args.AddOption(&merge, "-merge", "--merge", "-no-merge", "--no-merge",
"Enable or disable the merge phase.");
args.AddOption(&samplingType, "-hrtype", "--hrsamplingtype",
"Sampling type for hyperreduction.");
args.AddOption(&num_samples_req, "-nsr", "--nsr",
"number of samples we want to select for the sampling algorithm.");
args.AddOption(&rxdim, "-rxdim", "--rxdim",
"Basis dimension for displacement solution space.");
args.AddOption(&rvdim, "-rvdim", "--rvdim",
"Basis dimension for velocity solution space.");
args.AddOption(&hdim, "-hdim", "--hdim",
"Basis dimension for the nonlinear term.");
args.AddOption(&id_param, "-id", "--id", "Parametric index");
args.AddOption(&hyperreduce, "-hyp", "--hyperreduce", "-no-hyp",
"--no-hyperreduce", "Hyper reduce nonlinear term");
args.AddOption(&x_base_only, "-xbo", "--xbase-only", "-no-xbo",
"--not-xbase-only", "Use the displacement (X) basis to approximate velocity.");
args.AddOption(&def_ic, "-def-ic", "--deformation-ic", "-vel-ic",
"--velocity-ic",
"Use a deformation-, or velocity initial condition. Default is velocity IC.");
args.AddOption(&s, "-sc", "--scaling", "Scaling factor for initial condition.");
args.AddOption(&use_eqp, "-eqp", "--eqp", "-no-eqp", "--no-eqp",
"Use EQP instead of DEIM for the hyperreduction.");
args.AddOption(&writeSampleMesh, "-smesh", "--sample-mesh", "-no-smesh",
"--no-sample-mesh", "Write the sample mesh to file.");
args.AddOption(&preconditionNNLS, "-preceqp", "--preceqp", "-no-preceqp",
"--no-preceqp", "Precondition the NNLS system for EQP.");
args.AddOption(&tolNNLS, "-tolnnls", "--tol-nnls",
"Tolerance for NNLS solver.");
args.AddOption(&maxNNLSnnz, "-maxnnls", "--max-nnls",
"Maximum nnz for NNLS");
args.AddOption(&n_windows, "-ntw", "--n-time-windows",
"Number of time windows. Default is 0");
args.AddOption(&snap_step, "-sbs", "--subsampling-step",
"Subsamble every n-th snapshot. Default is 1, meaning no subsampling");
args.Parse();
if (!args.Good())
{
if (myid == 0)
{
args.PrintUsage(cout);
}
MPI_Finalize();
return 1;
}
if (myid == 0)
{
args.PrintOptions(cout);
}
const bool check = (offline && !merge && !online) || (!offline && merge
&& !online) || (!offline && !merge && online);
MFEM_VERIFY(check, "only one of offline, merge, or online must be true!");
StopWatch solveTimer, totalTimer;
totalTimer.Start();
// 3. Read the serial mesh from the given mesh file on all processors. We can
// handle triangular, quadrilateral, tetrahedral and hexahedral meshes
// with the same code.
Mesh *mesh = new Mesh(mesh_file, 1, 1);
int dim = mesh->Dimension();
// 4. Define the ODE solver used for time integration. Several implicit
// singly diagonal implicit Runge-Kutta (SDIRK) methods, as well as
// explicit Runge-Kutta methods are available.
ODESolver *ode_solver;
switch (ode_solver_type)
{
// Implicit L-stable methods
// To be added...
// Explicit methods
case 11:
ode_solver = new ForwardEulerSolver;
break;
case 12:
ode_solver = new RK2Solver(0.5);
break; // midpoint method
case 13:
ode_solver = new RK3SSPSolver;
break;
case 14:
ode_solver = new RK4Solver;
break;
case 15:
ode_solver = new GeneralizedAlphaSolver(0.5);
break;
// Implicit A-stable methods (not L-stable)
// To be added...
default:
if (myid == 0)
{
cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
}
delete mesh;
MPI_Finalize();
return 3;
}
// 5. Refine the mesh in serial to increase the resolution. In this example
// we do 'ser_ref_levels' of uniform refinement, where 'ser_ref_levels' is
// a command-line parameter.
for (int lev = 0; lev < ser_ref_levels; lev++)
{
mesh->UniformRefinement();
}
// 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
// this mesh further in parallel to increase the resolution. Once the
// parallel mesh is defined, the serial mesh can be deleted.
ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
delete mesh;
for (int lev = 0; lev < par_ref_levels; lev++)
{
pmesh->UniformRefinement();
}
// 7. Define the parallel vector finite element spaces representing the mesh
// deformation x_gf, the velocity v_gf, and the initial configuration,
// x_ref. Define also the elastic energy density, w_gf, which is in a
// discontinuous higher-order space. Since x and v are integrated in time
// as a system, we group them together in block vector vx, on the unique
// parallel degrees of freedom, with offsets given by array true_offset.
H1_FECollection fe_coll(order, dim);
ParFiniteElementSpace fespace(pmesh, &fe_coll, dim);
HYPRE_BigInt glob_size = fespace.GlobalTrueVSize();
if (myid == 0)
{
cout << "Number of velocity/deformation unknowns: " << glob_size << endl;
}
int true_size = fespace.TrueVSize();
Array<int> true_offset(3);
true_offset[0] = 0;
true_offset[1] = true_size;
true_offset[2] = 2 * true_size;
BlockVector vx(true_offset);
ParGridFunction v_gf, x_gf;
// Associate a FiniteElementSpace and true-dof data with the GridFunctions.
v_gf.MakeTRef(&fespace, vx, true_offset[0]);
x_gf.MakeTRef(&fespace, vx, true_offset[1]);
ParGridFunction x_ref(&fespace);
pmesh->GetNodes(x_ref);
L2_FECollection w_fec(order + 1, dim);
ParFiniteElementSpace w_fespace(pmesh, &w_fec);
ParGridFunction w_gf(&w_fespace);
// Basis params
bool update_right_SV = false;
bool isIncremental = false;
const std::string basisFileName = "basis" + std::to_string(id_param);
int max_num_snapshots = int(t_final / dt) + 2;
// The merge phase
if (merge)
{
totalTimer.Clear();
totalTimer.Start();
// Merge bases
if (x_base_only == false)
{
MergeBasis(true_size, nsets, max_num_snapshots, "V");
}
MergeBasis(true_size, nsets, max_num_snapshots, "X");
MergeBasis(true_size, nsets, max_num_snapshots, "H");
totalTimer.Stop();
if (myid == 0)
{
printf("Elapsed time for merging and building ROM basis: %e second\n",
totalTimer.RealTime());
}
MPI_Finalize();
return 0;
}
// 8. Set the initial conditions for v_gf, x_gf and vx, and define the
// boundary conditions on a beam-like mesh (see description above).
VectorFunctionCoefficient *velo = 0;
VectorFunctionCoefficient *deform = 0;
if (def_ic)
{
velo = new VectorFunctionCoefficient(dim, InitialVelocityIC2);
}
else
{
velo = new VectorFunctionCoefficient(dim, InitialVelocityIC1);
}
v_gf.ProjectCoefficient(*velo);
v_gf.SetTrueVector();
if (def_ic)
{
deform = new VectorFunctionCoefficient(dim, InitialDeformationIC2);
}
else
{
deform = new VectorFunctionCoefficient(dim, InitialDeformationIC1);
}
x_gf.ProjectCoefficient(*deform);
x_gf.SetTrueVector();
v_gf.SetFromTrueVector();
x_gf.SetFromTrueVector();
Array<int> ess_bdr(fespace.GetMesh()->bdr_attributes.Max());
ess_bdr = 0;
ess_bdr[0] = 1; // boundary attribute 1 (index 0) is fixed
Array<int> ess_tdof_list;
fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
// Store initial vx
BlockVector vx0 = BlockVector(vx);
BlockVector vx_diff = BlockVector(vx);
// Reduced order solution
Vector *wMFEM = 0;
CAROM::Vector *w = 0;
CAROM::Vector *w_v = 0;
CAROM::Vector *w_x = 0;
// Initialize reconstructed solution
Vector *v_rec = new Vector(v_gf.GetTrueVector());
Vector *x_rec = new Vector(x_gf.GetTrueVector());
CAROM::Vector *v_rec_librom = new CAROM::Vector(v_rec->GetData(), v_rec->Size(),
true, false);
CAROM::Vector *x_rec_librom = new CAROM::Vector(x_rec->GetData(), x_rec->Size(),
true, false);
// 9. Initialize the hyperelastic operator, the GLVis visualization and print
// the initial energies.
HyperelasticOperator oper(fespace, ess_tdof_list, visc, mu, K);
NeoHookeanModel *model = new NeoHookeanModel(mu, K);
HyperelasticOperator *soper = 0;
// Fill dvdt and dxdt
Vector dvxdt(true_size * 2);
Vector dvdt(dvxdt.GetData() + 0, true_size);
Vector dxdt(dvxdt.GetData() + true_size, true_size);
socketstream vis_v, vis_w;
if (visualization)
{
char vishost[] = "localhost";
int visport = 19916;
vis_v.open(vishost, visport);
vis_v.precision(8);
visualize(vis_v, pmesh, &x_gf, &v_gf, "Velocity", true);
// Make sure all ranks have sent their 'v' solution before initiating
// another set of GLVis connections (one from each rank)
MPI_Barrier(pmesh->GetComm());
vis_w.open(vishost, visport);
if (vis_w)
{
oper.GetElasticEnergyDensity(x_gf, w_gf);
vis_w.precision(8);
visualize(vis_w, pmesh, &x_gf, &w_gf, "Elastic energy density", true);
}
}
// Create data collection for solution output: either VisItDataCollection for
// ascii data files, or SidreDataCollection for binary data files.
DataCollection *dc = NULL;
if (visit)
{
if (offline)
dc = new VisItDataCollection("nlelast-fom", pmesh);
else
dc = new VisItDataCollection("nlelast-rom", pmesh);
dc->SetPrecision(8);
// To save the mesh using MFEM's parallel mesh format:
// dc->SetFormat(DataCollection::PARALLEL_FORMAT);
dc->RegisterField("x", &x_gf);
dc->RegisterField("v", &v_gf);
dc->SetCycle(0);
dc->SetTime(0.0);
dc->Save();
}
double ee0 = oper.ElasticEnergy(x_gf);
double ke0 = oper.KineticEnergy(v_gf);
if (myid == 0)
{
cout << "initial elastic energy (EE) = " << ee0 << endl;
cout << "initial kinetic energy (KE) = " << ke0 << endl;
cout << "initial total energy (TE) = " << (ee0 + ke0) << endl;
}
// 10. Create pROM object.
CAROM::BasisGenerator *basis_generator_v = 0;
CAROM::BasisGenerator *basis_generator_x = 0;
CAROM::BasisGenerator *basis_generator_H = 0;
if (offline)
{
CAROM::Options options(fespace.GetTrueVSize(), max_num_snapshots, 1,
update_right_SV);
if (x_base_only == false)
{
basis_generator_v = new CAROM::BasisGenerator(options, isIncremental,
basisFileName + "_V");
}
basis_generator_x = new CAROM::BasisGenerator(options, isIncremental,
basisFileName + "_X");
basis_generator_H = new CAROM::BasisGenerator(options, isIncremental,
basisFileName + "_H");
}
RomOperator *romop = 0;
const CAROM::Matrix *BV_librom = 0;
const CAROM::Matrix *BX_librom = 0;
const CAROM::Matrix *H_librom = 0;
const CAROM::Matrix *Hsinv = 0;
int nsamp_H = -1;
CAROM::SampleMeshManager *smm = nullptr;
CAROM::Vector *eqpSol = nullptr;
CAROM::Vector *eqpSol_S = nullptr;
CAROM::Vector *window_ids = nullptr;
CAROM::Vector *load_eqpsol = new CAROM::Vector(1,
false); // Will be resized later
// 11. Initialize ROM operator
// I guess most of this should be done on id =0
if (online)
{
// Read bases
CAROM::BasisReader *readerV = 0;
if (x_base_only)
{
readerV = new
CAROM::BasisReader("basisX"); // The basis for v uses the x basis instead.
rvdim = rxdim;
}
else
{
readerV = new CAROM::BasisReader("basisV");
}
BV_librom = readerV->getSpatialBasis();
if (rvdim == -1) // Change rvdim
rvdim = BV_librom->numColumns();
else
BV_librom = GetFirstColumns(rvdim,
BV_librom);
MFEM_VERIFY(BV_librom->numRows() == true_size, "");
if (myid == 0)
printf("reduced V dim = %d\n", rvdim);
CAROM::BasisReader readerX("basisX");
BX_librom = readerX.getSpatialBasis();
if (rxdim == -1) // Change rxdim
rxdim = BX_librom->numColumns();
else
BX_librom = GetFirstColumns(rxdim,
BX_librom);
MFEM_VERIFY(BX_librom->numRows() == true_size, "");
if (myid == 0)
printf("reduced X dim = %d\n", rxdim);
// Hyper reduce H
CAROM::BasisReader readerH("basisH");
H_librom = readerH.getSpatialBasis();
// Compute sample points
if (hdim == -1)
{
hdim = H_librom->numColumns();
}
CAROM::Matrix *Hsinv = new CAROM::Matrix(hdim, hdim, false);
MFEM_VERIFY(H_librom->numColumns() >= hdim, "");
if (H_librom->numColumns() > hdim)
H_librom = GetFirstColumns(hdim, H_librom);
if (myid == 0)
printf("reduced H dim = %d\n", hdim);
// Setup hyperreduction, using either EQP or sampled DOFs and a sample mesh.
ParFiniteElementSpace *sp_XV_space;
const IntegrationRule *ir0 = NULL;
if (ir0 == NULL)
{
const FiniteElement &fe = *fespace.GetFE(0);
ir0 = &(IntRules.Get(fe.GetGeomType(), 2 * fe.GetOrder() + 3));
}
// Timewindowing setup for EQP
int n_step = int(t_final / dt);
if (n_windows > 1)
{
window_ids = new CAROM::Vector(n_windows + 1, false);
get_window_ids(n_step, n_windows, window_ids);
}
if (use_eqp)
{
// EQP setup
eqpSol = new CAROM::Vector(ir0->GetNPoints() * fespace.GetNE(), true);
SetupEQP_snapshots(ir0, myid, &fespace, nsets, BV_librom,
GetSnapshotMatrix(fespace.GetTrueVSize(), nsets, max_num_snapshots, "X"),
vx0.GetBlock(0),
preconditionNNLS, tolNNLS, maxNNLSnnz, model, *eqpSol, window_ids, snap_step);
if (writeSampleMesh)
WriteMeshEQP(pmesh, myid, ir0->GetNPoints(), *eqpSol);
}
else
{
vector<int> num_sample_dofs_per_proc(num_procs);
if (num_samples_req != -1)
{
nsamp_H = num_samples_req;
}
else
{
nsamp_H = hdim;
}
Hsinv->setSize(nsamp_H, hdim);
vector<int> sample_dofs(nsamp_H);
// Setup hyperreduction using DEIM, GNAT, or S-OPT
CAROM::Hyperreduction hr(samplingType);
hr.ComputeSamples(H_librom,
hdim,
sample_dofs,
num_sample_dofs_per_proc,
*Hsinv,
myid,
num_procs,
nsamp_H);
// Construct sample mesh
const int nspaces = 1;
std::vector<ParFiniteElementSpace *> spfespace(nspaces);
spfespace[0] = &fespace;
smm = new CAROM::SampleMeshManager(spfespace);
vector<int> sample_dofs_empty;
vector<int> num_sample_dofs_per_proc_empty;
num_sample_dofs_per_proc_empty.assign(num_procs, 0);
smm->RegisterSampledVariable("V", 0, sample_dofs,
num_sample_dofs_per_proc);
smm->RegisterSampledVariable("X", 0, sample_dofs,
num_sample_dofs_per_proc);
smm->RegisterSampledVariable("H", 0, sample_dofs,
num_sample_dofs_per_proc);
smm->ConstructSampleMesh();
}
w = new CAROM::Vector(rxdim + rvdim, false);
w_v = new CAROM::Vector(rvdim, false);
w_x = new CAROM::Vector(rxdim, false);
*w = 0.0;
// Note that some of this could be done only on the ROM solver process,
// but it is tricky, since RomOperator assembles Bsp in parallel.
wMFEM = new Vector(&((*w)(0)), rxdim + rvdim);
// Initial conditions
Vector *w_v0 = 0;
Vector *w_x0 = 0;
int sp_size = 0;
Array<int> ess_tdof_list_sp; // Initialize sample essential boundary conditions
if (myid == 0 && !use_eqp)
{
sp_XV_space = smm->GetSampleFESpace(0);
sp_size = sp_XV_space->TrueVSize();
Array<int> sp_offset(3);
sp_offset[0] = 0;
sp_offset[1] = sp_size;
sp_offset[2] = 2 * sp_size;
// Initialize sp_p with initial conditions.
BlockVector sp_vx(sp_offset);
ParGridFunction sp_v_gf, sp_x_gf;
// 12. Set the initial conditions for v_gf, x_gf and vx, and define the
// boundary conditions on a beam-like mesh (see description above).
// Associate a FiniteElementSpace and true-dof data with the GridFunctions.
sp_v_gf.MakeTRef(sp_XV_space, sp_vx, sp_offset[0]);
sp_x_gf.MakeTRef(sp_XV_space, sp_vx, sp_offset[1]);