-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathmhd3d_stationary.py
1018 lines (877 loc) · 35.8 KB
/
mhd3d_stationary.py
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
# -*- coding: utf-8 -*-
from firedrake import *
import alfi
from alfi.stabilisation import *
from alfi.transfer import *
from pyop2.datatypes import IntType
import ufl.algorithms
from datetime import datetime
import argparse
import numpy
import sys
import os
from mpi4py import MPI
import petsc4py
petsc4py.PETSc.Sys.popErrorHandler()
# Parallel distribution parameters
distribution_parameters = {"partition": True, "overlap_type": (DistributedMeshOverlapType.VERTEX, 2)}
# Definition of Burman Stabilisation term to avoid oscillations in velocity u
def BurmanStab(B, C, wind, stab_weight, mesh):
n = FacetNormal(mesh)
h = FacetArea(mesh)
beta = avg(facet_avg(sqrt(inner(wind, wind)+1e-10)))
gamma1 = stab_weight # as chosen in doi:10.1016/j.apnum.2007.11.001
stabilisation_form = 0.5 * gamma1 * avg(h)**2 * beta * dot(jump(grad(B), n), jump(grad(C), n))*dS
return stabilisation_form
# Definition of outer Schur complement for the order ((E,B),(u,p)), which is Navier-Stokes block + Lorentz force
class SchurPC(AuxiliaryOperatorPC):
def form(self, pc, V, U):
mesh = V.function_space().mesh()
[u, p] = split(U)
[v, q] = split(V)
Z = V.function_space()
U_ = Function(Z)
state = self.get_appctx(pc)['state']
[u_n, p_n, B_n, E_n] = split(state)
eps = lambda x: sym(grad(x))
# Weak form of NS + Lorentz force
A = (
2/Re * inner(eps(u), eps(v))*dx
+ gamma * inner(div(u), div(v)) * dx
- inner(p, div(v)) * dx
- inner(div(u), q) * dx
+ S * inner(cross(B_n, cross(u, B_n)), v) * dx
)
# For H(div)-L2 discretization we have to add a DG-formulation of the advection and diffusion term
if discr in ["rt", "bdm"]:
h = CellVolume(mesh)/FacetArea(mesh)
uflux_int = 0.5*(dot(u, n) + theta*abs(dot(u, n)))*u
uflux_ext_1 = 0.5*(inner(u, n) + theta*abs(inner(u, n)))*u
uflux_ext_2 = 0.5*(inner(u, n) - theta*abs(inner(u, n)))*u_ex
A_DG = (
- 1/Re * inner(avg(2*sym(grad(u))), 2*avg(outer(v, n))) * dS
- 1/Re * inner(avg(2*sym(grad(v))), 2*avg(outer(u, n))) * dS
+ 1/Re * sigma/avg(h) * inner(2*avg(outer(u, n)), 2*avg(outer(v, n))) * dS
- inner(outer(v, n), 2/Re*sym(grad(u))) * ds
- inner(outer(u-u_ex, n), 2/Re*sym(grad(v))) * ds(bcs_ids_apply)
+ 1/Re*(sigma/h)*inner(v, u-u_ex) * ds(bcs_ids_apply)
- advect * dot(u, div(outer(v, u))) * dx
+ advect * dot(v('+')-v('-'), uflux_int('+')-uflux_int('-'))*dS
+ advect * dot(v, uflux_ext_1) * ds
+ advect * dot(v, uflux_ext_2) * ds(bcs_ids_apply)
)
if bcs_ids_dont_apply is not None:
A_DG += (
- inner(outer(u, n), 2/Re*sym(grad(v))) * ds(bcs_ids_dont_apply)
+ 1/Re*(sigma/h)*inner(v, u) * ds(bcs_ids_dont_apply)
)
# Linearize A_DG
A_DG = action(A_DG, U_)
A_DG_linear = derivative(A_DG, U_, U)
A_DG_linear = replace(A_DG_linear, {split(U_)[0]: u_n, split(U_)[1]: p_n})
A = A + A_DG_linear
elif discr in ["cg"]:
A = A + advect * inner(dot(grad(u_n), u), v) * dx + advect * inner(dot(grad(u), u_n), v) * dx
if stab:
stabilisation2 = BurmanStabilisation(V.function_space().sub(0), state=z_last_u, h=FacetArea(mesh), weight=stab_weight)
stabilisation_form_2 = stabilisation2.form(u, v)
A = A + advect * stabilisation_form_2
bcs = [DirichletBC(V.function_space().sub(0), 0, "on_boundary"),
PressureFixBC(V.function_space().sub(1), 0, 1),
]
return (A, bcs)
# Definition of outer Schur complement for the order ((u,p),(B,E))
class SchurPCBE(AuxiliaryOperatorPC):
def form(self, pc, V, U):
[B, E] = split(U)
[C, Ff] = split(V)
state = self.get_appctx(pc)['state']
[u_n, p_n, B_n, E_n] = split(state)
A = (
+ inner(E, Ff) * dx
+ inner(cross(u_n, B), Ff) * dx
- 1/Rem * inner(B, curl(Ff)) * dx
+ inner(curl(E), C) * dx
+ 1/Rem * inner(div(B), div(C)) * dx
+ gamma2 * inner(div(B), div(C)) * dx
)
bcs = [DirichletBC(V.function_space().sub(0), 0, "on_boundary"),
DirichletBC(V.function_space().sub(1), 0, "on_boundary"),
]
return (A, bcs)
# We fix the pressure at one vertex on every level
class PressureFixBC(DirichletBC):
def __init__(self, V, val, subdomain, method="topological"):
super().__init__(V, val, subdomain, method)
sec = V.dm.getDefaultSection()
dm = V.mesh().topology_dm
coordsSection = dm.getCoordinateSection()
dim = dm.getCoordinateDim()
coordsVec = dm.getCoordinatesLocal()
(vStart, vEnd) = dm.getDepthStratum(0)
indices = []
for pt in range(vStart, vEnd):
x = dm.getVecClosure(coordsSection, coordsVec, pt).reshape(-1, dim).mean(axis=0)
if x.dot(x) == 0.0: # fix [0, 0] in original mesh coordinates (bottom left corner)
if dm.getLabelValue("pyop2_ghost", pt) == -1:
indices = [pt]
break
nodes = []
for i in indices:
if sec.getDof(i) > 0:
nodes.append(sec.getOffset(i))
if V.mesh().comm.rank == 0:
nodes = [0]
else:
nodes = []
self.nodes = numpy.asarray(nodes, dtype=IntType)
if len(self.nodes) > 0:
print("Fixing nodes %s" % self.nodes)
import sys
sys.stdout.flush()
# Definition of different solver parameters
# Increase buffer size for MUMPS solver
ICNTL_14 = 5000
tele_reduc_fac = int(MPI.COMM_WORLD.size/4)
if tele_reduc_fac < 1:
tele_reduc_fac = 1
# LU solver
lu = {
"snes_type": "newtonls",
"snes_monitor": None,
"snes_atol": 1.0e-8,
"snes_rtol": 1.0e-13,
"snes_linesearch_type": "l2",
"snes_linesearch_maxstep": 1,
"mat_type": "aij",
"ksp_type": "preonly",
"pc_type": "lu",
"pc_factor_mat_solver_type": "mumps",
"mat_mumps_icntl_14": ICNTL_14,
}
nsfsstar = {
"ksp_type": "fgmres",
"ksp_atol": 1.0e-7,
"ksp_rtol": 1.0e-7,
"ksp_max_it": 2,
"pc_type": "fieldsplit",
"pc_fieldsplit_type": "schur",
"pc_fieldsplit_schur_factorization_type": "upper",
"pc_fieldsplit_schur_precondition": "user",
"fieldsplit_0": {
"ksp_type": "gmres",
"ksp_max_it": 1,
"ksp_norm_type": "unpreconditioned",
"pc_type": "mg",
"pc_mg_cycle_type": "v",
"pc_mg_type": "full",
"mg_levels_ksp_type": "fgmres",
"mg_levels_ksp_convergence_test": "skip",
"mg_levels_ksp_max_it": 6,
"mg_levels_ksp_norm_type": "unpreconditioned",
"mg_levels_pc_type": "python",
"mg_levels_pc_python_type": "firedrake.ASMStarPC",
"mg_levels_pc_star_backend": "tinyasm",
# "mg_levels_pc_star_construct_dim": 0,
# "mg_levels_pc_star_sub_sub_ksp_type": "preonly",
# "mg_levels_pc_star_sub_sub_pc_type": "lu",
# "mg_levels_pc_star_sub_sub_pc_factor_mat_solver_type": "umfpack",
"mg_coarse_ksp_type": "richardson",
"mg_coarse_ksp_max_it": 1,
"mg_coarse_ksp_norm_type": "unpreconditioned",
"mg_coarse_pc_type": "python",
"mg_coarse_pc_python_type": "firedrake.AssembledPC",
"mg_coarse_assembled": {
"mat_type": "aij",
"pc_type": "telescope",
"pc_telescope_reduction_factor": tele_reduc_fac,
"pc_telescope_subcomm_type": "contiguous",
"telescope_pc_type": "lu",
"telescope_pc_factor_mat_solver_type": "mumps",
"telescope_pc_factor_mat_mumps_icntl_14": ICNTL_14,
}
},
"fieldsplit_1": {
"ksp_type": "preonly",
"pc_type": "python",
"pc_python_type": "alfi.solver.DGMassInv"
},
}
# Monolithic macrostar solver for (E,B)-block
nsfsmacrostar = {
"ksp_type": "fgmres",
"ksp_atol": 1.0e-7,
"ksp_rtol": 1.0e-7,
"ksp_max_it": 2,
"pc_type": "fieldsplit",
"pc_fieldsplit_type": "schur",
"pc_fieldsplit_schur_factorization_type": "full",
"pc_fieldsplit_schur_precondition": "user",
"fieldsplit_0": {
"ksp_type": "gmres",
"ksp_max_it": 1,
"ksp_norm_type": "unpreconditioned",
"pc_type": "mg",
"pc_mg_cycle_type": "v",
"pc_mg_type": "full",
"mg_levels_ksp_type": "fgmres",
"mg_levels_ksp_convergence_test": "skip",
"mg_levels_ksp_max_it": 6,
"mg_levels_ksp_norm_type": "unpreconditioned",
"mg_levels_pc_type": "python",
"mg_levels_pc_python_type": "firedrake.PatchPC",
"mg_levels_patch_pc_patch_save_operators": True,
"mg_levels_patch_pc_patch_partition_of_unity": False,
"mg_levels_patch_pc_patch_sub_mat_type": "seqaij",
"mg_levels_patch_pc_patch_construct_dim": 0,
"mg_levels_patch_pc_patch_construct_type": "python",
"mg_levels_patch_pc_patch_construct_python_type": "alfi.MacroStar",
"mg_levels_patch_sub_ksp_type": "preonly",
"mg_levels_patch_sub_pc_type": "lu",
"mg_levels_patch_sub_pc_factor_mat_solver_type": "umfpack",
"mg_coarse_ksp_type": "richardson",
"mg_coarse_ksp_max_it": 1,
"mg_coarse_ksp_norm_type": "unpreconditioned",
"mg_coarse_pc_type": "python",
"mg_coarse_pc_python_type": "firedrake.AssembledPC",
"mg_coarse_assembled": {
"mat_type": "aij",
"pc_type": "telescope",
"pc_telescope_reduction_factor": 1,
"pc_telescope_subcomm_type": "contiguous",
"telescope_pc_type": "lu",
"telescope_pc_factor_mat_solver_type": "mumps",
"telescope_pc_factor_mat_mumps_icntl_14": ICNTL_14,
}
},
"fieldsplit_1": {
"ksp_type": "preonly",
"pc_type": "python",
"pc_python_type": "alfi.solver.DGMassInv"
},
}
# Fieldsplit solver for outer Schur complement with monolithic star solver
outerschurstar = {
"ksp_type": "fgmres",
"ksp_atol": 1.0e-7,
"ksp_rtol": 1.0e-7,
"pc_type": "python",
"pc_python_type": "__main__.SchurPCBE",
"aux_mg_transfer_manager": __name__ + ".transfermanager",
"ksp_max_it": 2,
"aux": {
"ksp_type": "gmres",
"ksp_max_it": 1,
"ksp_norm_type": "unpreconditioned",
"pc_type": "mg",
"pc_mg_cycle_type": "v",
"pc_mg_type": "full",
"mg_levels_ksp_type": "fgmres",
"mg_levels_ksp_convergence_test": "skip",
"mg_levels_ksp_max_it": 6,
"mg_levels_ksp_norm_type": "unpreconditioned",
"mg_levels_pc_type": "python",
"mg_levels_pc_python_type": "firedrake.ASMStarPC",
"mg_levels_pc_star_backend": "tinyasm",
# "mg_levels_pc_python_type": "firedrake.PatchPC",
# "mg_levels_patch_pc_patch_save_operators": True,
# "mg_levels_patch_pc_patch_partition_of_unity": False,
# "mg_levels_patch_pc_patch_sub_mat_type": "seqaij",
# "mg_levels_patch_pc_patch_construct_dim": 0,
# "mg_levels_patch_pc_patch_construct_type": "star",
# "mg_levels_patch_sub_ksp_type": "preonly",
# "mg_levels_patch_sub_pc_type": "lu",
# "mg_levels_patch_sub_pc_factor_mat_solver_type": "umfpack",
"mg_coarse_ksp_type": "richardson",
"mg_coarse_ksp_max_it": 1,
"mg_coarse_ksp_norm_type": "unpreconditioned",
"mg_coarse_pc_type": "python",
"mg_coarse_pc_python_type": "firedrake.AssembledPC",
"mg_coarse_assembled": {
"mat_type": "aij",
"pc_type": "telescope",
"pc_telescope_reduction_factor": tele_reduc_fac,
"pc_telescope_subcomm_type": "contiguous",
"telescope_pc_type": "lu",
"telescope_pc_factor_mat_solver_type": "mumps",
"telescope_pc_factor_mat_mumps_icntl_14": ICNTL_14,
}
},
}
# Fieldsplit solver for outer Schur complement with monolithic macrostar solver
outerschurmacrostar = {
"ksp_type": "fgmres",
"ksp_atol": 1.0e-7,
"ksp_rtol": 1.0e-7,
"pc_type": "python",
"pc_python_type": "__main__.SchurPCBE",
"aux_mg_transfer_manager": __name__ + ".transfermanager",
"ksp_max_it": 2,
"aux": {
"ksp_type": "gmres",
"ksp_max_it": 1,
"ksp_norm_type": "unpreconditioned",
"pc_type": "mg",
"pc_mg_cycle_type": "v",
"pc_mg_type": "full",
"mg_levels_ksp_type": "fgmres",
"mg_levels_ksp_convergence_test": "skip",
"mg_levels_ksp_max_it": 6,
"mg_levels_ksp_norm_type": "unpreconditioned",
"mg_levels_pc_type": "python",
"mg_levels_pc_python_type": "firedrake.PatchPC",
"mg_levels_patch_pc_patch_save_operators": True,
"mg_levels_patch_pc_patch_partition_of_unity": False,
"mg_levels_patch_pc_patch_sub_mat_type": "seqaij",
"mg_levels_patch_pc_patch_construct_dim": 0,
"mg_levels_patch_pc_patch_construct_type": "python",
"mg_levels_patch_pc_patch_construct_python_type": "alfi.MacroStar",
"mg_levels_patch_sub_ksp_type": "preonly",
"mg_levels_patch_sub_pc_type": "lu",
"mg_levels_patch_sub_pc_factor_mat_solver_type": "umfpack",
"mg_coarse_ksp_type": "richardson",
"mg_coarse_ksp_max_it": 1,
"mg_coarse_ksp_norm_type": "unpreconditioned",
"mg_coarse_pc_type": "python",
"mg_coarse_pc_python_type": "firedrake.AssembledPC",
"mg_coarse_assembled": {
"mat_type": "aij",
"pc_type": "telescope",
"pc_telescope_reduction_factor": 1,
"pc_telescope_subcomm_type": "contiguous",
"telescope_pc_type": "lu",
"telescope_pc_factor_mat_solver_type": "mumps",
"telescope_pc_factor_mat_mumps_icntl_14": ICNTL_14,
}
},
}
# LU solver for outer Schur complement
outerschurlu = {
"ksp_type": "gmres",
"ksp_max_it": 2,
"pc_type": "python",
"pc_python_type": "__main__.SchurPCBE",
"aux_pc_type": "lu",
"aux_pc_factor_mat_solver_type": "mumps",
"aux_mat_mumps_icntl_14": ICNTL_14,
}
# Main solver
fs2by2 = {
"snes_type": "newtonls",
"snes_max_it": 10,
"snes_linesearch_type": "basic",
"snes_linesearch_maxstep": 1.0,
"snes_rtol": 1.0e-10,
"snes_atol": 4.0e-6,
"snes_monitor": None,
"ksp_type": "fgmres",
"ksp_max_it": 35,
"ksp_atol": 3.0e-6,
"ksp_rtol": 1.0e-7,
"ksp_monitor_true_residual": None,
"ksp_converged_reason": None,
"mat_type": "aij",
"pc_type": "fieldsplit",
"pc_fieldsplit_type": "schur",
"pc_fieldsplit_schur_factorization_type": "upper",
"pc_fieldsplit_0_fields": "0,1",
"pc_fieldsplit_1_fields": "2,3",
}
# Main solver with LU solver for (u,p)-block
fs2by2nslu = {
"snes_type": "newtonls",
"snes_max_it": 10,
"snes_linesearch_type": "l2",
"snes_linesearch_maxstep": 1.0,
"snes_rtol": 1.0e-10,
"snes_atol": 4.0e-6,
"snes_monitor": None,
"ksp_type": "fgmres",
"ksp_max_it": 35,
"ksp_atol": 3.0e-6,
"ksp_rtol": 1.0e-7,
"ksp_monitor_true_residual": None,
"ksp_converged_reason": None,
"mat_type": "aij",
"pc_type": "fieldsplit",
"pc_fieldsplit_type": "schur",
"pc_fieldsplit_schur_factorization_type": "upper",
"pc_fieldsplit_0_fields": "0,1",
"pc_fieldsplit_1_fields": "2,3",
"fieldsplit_0_ksp_type": "preonly",
"fieldsplit_0_pc_type": "lu",
"fieldsplit_0_pc_factor_mat_solver_type": "mumps",
"fieldsplit_0_mat_mumps_icntl_14": ICNTL_14,
}
# Main solver with LU solver for Schur complement
fs2by2slu = {
"snes_type": "newtonls",
"snes_max_it": 10,
"snes_linesearch_type": "basic",
"snes_linesearch_maxstep": 1.0,
"snes_rtol": 1.0e-10,
"snes_atol": 4.0e-6,
"snes_monitor": None,
"ksp_type": "fgmres",
"ksp_max_it": 35,
"ksp_atol": 3.0e-6,
"ksp_rtol": 1.0e-7,
"ksp_monitor_true_residual": None,
"ksp_converged_reason": None,
"mat_type": "aij",
"pc_type": "fieldsplit",
"pc_fieldsplit_type": "schur",
"pc_fieldsplit_schur_factorization_type": "upper",
"pc_fieldsplit_0_fields": "0,1",
"pc_fieldsplit_1_fields": "2,3",
"fieldsplit_1": outerschurlu,
}
# Main solver with LU solver for (E,B)-block and Schur complement
# Can be used to determine how well the outer Schur complement is approximated by the different linerizations
fs2by2lu = {
"snes_type": "newtonls",
"snes_max_it": 10,
"snes_linesearch_type": "basic",
"snes_linesearch_maxstep": 1.0,
"snes_rtol": 1.0e-10,
"snes_atol": 4.0e-6,
"snes_monitor": None,
"ksp_type": "fgmres",
"ksp_max_it": 35,
"ksp_atol": 3.0e-6,
"ksp_rtol": 1.0e-7,
"ksp_monitor_true_residual": None,
"ksp_converged_reason": None,
"mat_type": "aij",
"pc_type": "fieldsplit",
"pc_fieldsplit_type": "schur",
"pc_fieldsplit_schur_factorization_type": "upper",
"pc_fieldsplit_0_fields": "0,1",
"pc_fieldsplit_1_fields": "2,3",
"fieldsplit_0_ksp_type": "preonly",
"fieldsplit_0_pc_type": "lu",
"fieldsplit_0_pc_factor_mat_solver_type": "mumps",
"fieldsplit_0_mat_mumps_icntl_14": ICNTL_14,
"fieldsplit_1": outerschurlu,
}
solvers = {"lu": lu, "fs2by2": fs2by2, "fs2by2nslu": fs2by2nslu, "fs2by2slu": fs2by2slu, "fs2by2lu": fs2by2lu}
# Definition of problem parameters
parser = argparse.ArgumentParser(add_help=False)
parser.add_argument("--baseN", type=int, default=10)
parser.add_argument("--k", type=int, default=3)
parser.add_argument("--nref", type=int, default=1)
parser.add_argument("--Re", nargs='+', type=float, default=[1])
parser.add_argument("--Rem", nargs='+', type=float, default=[1])
parser.add_argument("--gamma", type=float, default=10000)
parser.add_argument("--gamma2", type=float, default=0)
parser.add_argument("--advect", type=float, default=1)
parser.add_argument("--S", nargs='+', type=float, default=[1])
parser.add_argument("--hierarchy", choices=["bary", "uniform"], default="bary")
parser.add_argument("--solver-type", choices=list(solvers.keys()), default="lu")
parser.add_argument("--testproblem", choices=["ldc", "3dgenerator"], default="Wathen")
parser.add_argument("--discr", choices=["rt", "bdm", "cg"], required=True)
parser.add_argument("--linearisation", choices=["picard", "mdp", "newton"], required=True)
parser.add_argument("--stab", default=False, action="store_true")
parser.add_argument("--checkpoint", default=False, action="store_true")
parser.add_argument("--output", default=False, action="store_true")
args, _ = parser.parse_known_args()
baseN = args.baseN
k = args.k
nref = args.nref
Re = Constant(args.Re[0])
Rem = Constant(args.Rem[0])
gamma = Constant(args.gamma)
gamma2 = Constant(args.gamma2)
S = Constant(args.S[0])
hierarchy = args.hierarchy
solver_type = args.solver_type
testproblem = args.testproblem
gamma2 = Constant(args.gamma2)
advect = Constant(args.advect)
linearisation = args.linearisation
discr = args.discr
stab = args.stab
checkpoint = args.checkpoint
output = args.output
# Stabilisation weight for BurmanStabilisation
stab_weight = Constant(3e-3)
if len(args.Re) != 1 and len(args.Rem) != 1 and len(args.S) != 1:
raise ValueError("Re, Rem and S cannot all contain more than one element at the same time")
if discr == "cg" and hierarchy != "bary":
raise ValueError("SV is only stable on barycentric refined grids")
if k < 3 and hierarchy == "bary":
raise ValueError("Scott Vogelius is not stable for k<3")
# Define the base Mesh for the different problems
if testproblem == "3dgenerator":
base = BoxMesh(5*baseN, baseN, baseN, 5, 1, 1, distribution_parameters=distribution_parameters)
else:
base = UnitCubeMesh(baseN, baseN, baseN, distribution_parameters=distribution_parameters)
# Callbacks called before and after mesh refinement
def before(dm, i):
for p in range(*dm.getHeightStratum(1)):
dm.setLabelValue("prolongation", p, i+1)
def after(dm, i):
for p in range(*dm.getHeightStratum(1)):
dm.setLabelValue("prolongation", p, i+2)
# Create Mesh hierarchy
if hierarchy == "bary":
mh = alfi.BaryMeshHierarchy(base, nref, callbacks=(before, after))
elif hierarchy == "uniformbary":
bmesh = Mesh(bary(base._plex), distribution_parameters={"partition": False})
mh = MeshHierarchy(bmesh, nref, reorder=True, callbacks=(before, after),
distribution_parameters=distribution_parameters)
elif hierarchy == "uniform":
mh = MeshHierarchy(base, nref, reorder=True, callbacks=(before, after),
distribution_parameters=distribution_parameters)
else:
raise NotImplementedError("Only know bary, uniformbary and uniform for the hierarchy.")
# Change mesh from [0,1]^3 to [-0.5,0.5]^3
#if testproblem != "3dgenerator":
# for m in mh:
# m.coordinates.dat.data[:, 0] -= 0.5
# m.coordinates.dat.data[:, 1] -= 0.5
# m.coordinates.dat.data[:, 2] -= 0.5
mesh = mh[-1]
area = assemble(Constant(1, domain=mh[0])*dx)
def message(msg):
if mesh.comm.rank == 0:
warning(msg)
# Define mixed function spaces
if discr == "rt":
Vel = FiniteElement("N1div", mesh.ufl_cell(), k, variant="integral(3)")
V = FunctionSpace(mesh, Vel)
elif discr == "bdm":
Vel = FiniteElement("N2div", mesh.ufl_cell(), k, variant="integral(3)")
V = FunctionSpace(mesh, Vel)
elif discr == "cg":
V = VectorFunctionSpace(mesh, "CG", k)
Q = FunctionSpace(mesh, "DG", k-1) # p
Rel = FiniteElement("N1curl", mesh.ufl_cell(), k, variant="integral(3)")
R = FunctionSpace(mesh, Rel) # E
Wel = FiniteElement("N1div", mesh.ufl_cell(), k, variant="integral(3)")
W = FunctionSpace(mesh, Wel)
Z = MixedFunctionSpace([V, Q, W, R])
z = Function(Z)
(u, p, B, E) = split(z)
(v, q, C, Ff) = split(TestFunction(Z))
# used for BurmanStabilisation
z_last_u = Function(V)
# For continuation we might want to start from checkpoint
if checkpoint:
try:
if len(args.S) == 1 or len(args.Rem) == 1:
chk = DumbCheckpoint("dump/"+str(float(S*Rem))+str(linearisation)+str(testproblem), mode=FILE_READ)
else:
chk = DumbCheckpoint("dump/"+str(float(Rem))+str(linearisation)+str(testproblem), mode=FILE_READ)
chk.load(z)
except Exception as e:
message(e)
(u_, p_, B_, E_) = z.split()
z_last_u.assign(u_)
(x, y, zz) = SpatialCoordinate(Z.mesh())
n = FacetNormal(mesh)
eps = lambda x: sym(grad(x))
# Base weak form of problem
F = (
2/Re * inner(eps(u), eps(v))*dx
# + advect * inner(dot(grad(u), u), v) * dx
+ gamma * inner(div(u), div(v)) * dx
+ S * inner(cross(B, E), v) * dx
+ S * inner(cross(B, cross(u, B)), v) * dx
- inner(p, div(v)) * dx
- inner(div(u), q) * dx
+ inner(E, Ff) * dx
+ inner(cross(u, B), Ff) * dx
- 1/Rem * inner(B, curl(Ff)) * dx
+ inner(curl(E), C) * dx
+ 1/Rem * inner(div(B), div(C)) * dx
+ gamma2 * inner(div(B), div(C)) * dx
)
# Compute RHS for Method of Manufactured Solution (MMS)
def compute_rhs(u_ex, B_ex, p_ex, E_ex):
E_ex_ = interpolate(E_ex, R)
f1 = (-2/Re * div(eps(u_ex)) + advect * dot(grad(u_ex), u_ex) - gamma * grad(div(u_ex))
+ grad(p_ex) + S * cross(B_ex, (E_ex + cross(u_ex, B_ex))))
f2 = + curl(E_ex_) - 1/Rem * grad(div(B_ex)) - gamma2 * grad(div(B_ex))
f3 = -1/Rem * curl(B_ex) + E_ex + cross(u_ex, B_ex)
return (f1, f2, f3)
if testproblem == "3dgenerator":
# Hartmann number
Ha = sqrt(S*Rem*Re)
# Problem parameter
B0 = Constant(1.0)
x_on = Constant(2.0)
x_off = Constant(2.5)
delta = Constant(0.1)
B_z_gen = (B0/2)*(tanh((x-x_on)/delta) - tanh((x-x_off)/delta))
# *_ex are just defined for boundary conditions. We don't know exact solution
u_ex = Constant((1, 0, 0), domain=mesh)
B_ex = as_vector([Constant(0, domain=mesh), Constant(0, domain=mesh), B_z_gen])
E_ex = 1/Rem * curl(B_ex)
p_ex = Constant(0, domain=mesh)
# Inflow BCs for u at x=0, outflow at x=1, no-slip on 4 sides
# On what ids of the boundary do we want to apply the boundary conditions
# This is needed for DG-Form of H(div)-L2 formulation
bcs_ids_apply = (1)
bcs_ids_dont_apply = (3, 4, 5, 6)
bcs = [DirichletBC(Z.sub(0), u_ex, bcs_ids_apply),
DirichletBC(Z.sub(0), Constant((0., 0., 0.)), bcs_ids_dont_apply),
DirichletBC(Z.sub(2), B_ex, "on_boundary"),
DirichletBC(Z.sub(3), 0, "on_boundary"),
PressureFixBC(Z.sub(1), 0, 1)]
rhs = None # because rhs is zero for this problem
# Do we know what the exact solution of the problem is?
solution_known = False
# Do the boundary conditions depend on parameter thats change during continuation
bc_varying = True
elif testproblem == "ldc":
# example taken from https://doi.org/10.1016/j.jcp.2016.04.019
u_ex = Constant((1, 0, 0), domain=mesh)
B_ex = Constant((0, 1, 0), domain=mesh)
B_ex = project(B_ex, W)
bcs_ids_apply = 4
bcs_ids_dont_apply = (1, 2, 3, 5, 6)
bcs = [DirichletBC(Z.sub(0), u_ex, bcs_ids_apply), # 4 == upper boundary (y==1)
DirichletBC(Z.sub(0), 0, bcs_ids_dont_apply),
DirichletBC(Z.sub(2), B_ex, "on_boundary"),
DirichletBC(Z.sub(3), 0, "on_boundary"),
PressureFixBC(Z.sub(1), 0, 1)]
rhs = None
solution_known = False
bc_varying = False
if solution_known:
u_ex_ = interpolate(u_ex, V)
B_ex_ = interpolate(B_ex, W)
p_ex_ = interpolate(p_ex, Q)
E_ex_ = interpolate(E_ex, R)
# Add Burman Stabilisation
if stab:
initial = interpolate(as_vector([sin(y), x, x]), V)
z_last_u.assign(initial)
stabilisation = BurmanStabilisation(Z.sub(0), state=z_last_u, h=FacetArea(mesh), weight=stab_weight)
stabilisation_form_u = stabilisation.form(u, v)
F += (advect * stabilisation_form_u)
# For H(div)-L2 discretization we have to add a DG-formulation of the advection and diffusion term
if discr in ["rt", "bdm"]:
h = CellVolume(mesh)/FacetArea(mesh)
sigma = Constant(1) * Z.sub(0).ufl_element().degree()**2 # penalty term
theta = Constant(1) # theta=1 means upwind discretization of nonlinear advection term
uflux_int = 0.5*(dot(u, n) + theta*abs(dot(u, n)))*u
uflux_ext_1 = 0.5*(inner(u, n) + theta*abs(inner(u, n)))*u
uflux_ext_2 = 0.5*(inner(u, n) - theta*abs(inner(u, n)))*u_ex
F_DG = (
- 1/Re * inner(avg(2*sym(grad(u))), 2*avg(outer(v, n))) * dS
- 1/Re * inner(avg(2*sym(grad(v))), 2*avg(outer(u, n))) * dS
+ 1/Re * sigma/avg(h) * inner(2*avg(outer(u, n)), 2*avg(outer(v, n))) * dS
- inner(outer(v, n), 2/Re*sym(grad(u))) * ds
- inner(outer(u-u_ex, n), 2/Re*sym(grad(v))) * ds(bcs_ids_apply)
+ 1/Re*(sigma/h)*inner(v, u-u_ex) * ds(bcs_ids_apply)
- advect * dot(u, div(outer(v, u))) * dx
+ advect * dot(v('+')-v('-'), uflux_int('+')-uflux_int('-')) * dS
+ advect * dot(v, uflux_ext_1) * ds
+ advect * dot(v, uflux_ext_2) * ds(bcs_ids_apply)
)
if bcs_ids_dont_apply is not None:
F_DG += (
- inner(outer(u, n), 2/Re*sym(grad(v))) * ds(bcs_ids_dont_apply)
+ 1/Re*(sigma/h)*inner(v, u) * ds(bcs_ids_dont_apply)
)
F += F_DG
elif discr == "cg":
F += advect * inner(dot(grad(u), u), v) * dx
if rhs is not None:
F -= inner(f1, v) * dx + inner(f3, Ff) * dx + inner(f2, C) * dx
# Definition of the three different linearizations
w = TrialFunction(Z)
[w_u, w_p, w_B, w_E] = split(w)
J_newton = ufl.algorithms.expand_derivatives(derivative(F, z, w))
if linearisation == "newton":
J = J_newton
elif linearisation == "mdp":
J_mdp = (
J_newton
- inner(cross(w_u, B), Ff) * dx # G
)
J = J_mdp
elif linearisation == "picard":
J_picard = (
J_newton
- S * inner(cross(w_B, E), v) * dx # J_tilde
- S * inner(cross(w_B, cross(u, B)), v) * dx # D_1_tilde
- S * inner(cross(B, cross(u, w_B)), v) * dx # D_2_tilde
- inner(cross(u, w_B), Ff) * dx # G_tilde
)
J = J_picard
else:
raise ValueError("only know newton, mdp and picard as linearisation method")
problem = NonlinearVariationalProblem(F, z, bcs, J=J)
appctx = {"Re": Re, "gamma": gamma, "nu": 1/Re, "Rem": Rem, "gamma2": gamma2}
params = solvers[args.solver_type]
# Depending on the Mesh Hierarchy we have to use star or macrostar solver
if args.solver_type in ["fs2by2", "fs2by2slu"]:
params["fieldsplit_0"] = nsfsstar if hierarchy == "uniform" else nsfsmacrostar
if args.solver_type in ["fs2by2", "fs2by2mlu"]:
params["fieldsplit_1"] = outerschurstar if hierarchy == "uniform" else outerschurmacrostar
if args.hierarchy == "bary":
params["snes_linesearch_type"] = "l2"
import pprint
pprint.pprint(params)
# Definition of solver and transfer operators
solver = NonlinearVariationalSolver(problem, solver_parameters=params, options_prefix="", appctx=appctx)
qtransfer = NullTransfer()
Etransfer = NullTransfer()
vtransfer = SVSchoeberlTransfer((1/Re, gamma), 2, hierarchy)
dgtransfer = DGInjection()
transfers = {
Q.ufl_element(): (prolong, restrict, qtransfer.inject),
VectorElement("DG", mesh.ufl_cell(), args.k): (dgtransfer.prolong, restrict, dgtransfer.inject),
VectorElement("DG", mesh.ufl_cell(), args.k-1): (dgtransfer.prolong, restrict, dgtransfer.inject),
}
# On barycentric refined grids we need special prolongation operators
if hierarchy == "bary":
transfers[V.ufl_element()] = (vtransfer.prolong, vtransfer.restrict, inject)
transfermanager = TransferManager(native_transfers=transfers)
solver.set_transfer_manager(transfermanager)
results = {}
res = args.Re
rems = args.Rem
Ss = args.S
def run(re, rem, s):
(u, p, B, E) = z.split()
Re.assign(re)
Rem.assign(rem)
S.assign(s)
# Indices for output depending on Re-S, Re-Rem or S-Rem table
if len(args.S) == 1 or len(args.Rem) == 1:
ind1 = re
ind2 = s*rem
else:
ind1 = re*s
ind2 = rem
if bc_varying:
global p_ex
u_ex_ = interpolate(u_ex, V)
B_ex_ = interpolate(B_ex, W)
E_ex_ = interpolate(E_ex, R)
pintegral = assemble(p_ex*dx)
p_ex = p_ex - Constant(pintegral/area)
# bcs[0] is u, bcs[1] is B, bcs[2] is E
bcs[0].function_arg = u_ex_
if checkpoint:
try:
chk = DumbCheckpoint("dump/"+str(float(ind2))+str(linearisation)+str(testproblem), mode=FILE_READ)
chk.load(z)
except Exception as e:
message(e)
if mesh.comm.rank == 0:
print(GREEN % ("Solving for #dofs = %s, Re = %s, Rem = %s, gamma = %s, S = %s, baseN = %s, nref = %s, "
"linearisation = %s, testproblem = %s, discr = %s, k = %s"
% (Z.dim(), float(re), float(rem), float(gamma), float(S), int(baseN), int(nref),
linearisation, testproblem, discr, int(float(k)))), flush=True)
# Update z_last_u in Burman Stabilisation
if stab:
stabilisation.update(z.split()[0])
z_last_u.assign(u)
# For the continuation in Re we set up the solver again
if float(Re) == 100:
if mesh.comm.rank == 0:
print("Setting up new solver", flush=True)
global solver
problem = NonlinearVariationalProblem(F, z, bcs, J=J)
solver = NonlinearVariationalSolver(problem, solver_parameters=params, options_prefix="", appctx=appctx)
solver.set_transfer_manager(transfermanager)
# Solve the problem and measure time
start = datetime.now()
solver.solve()
end = datetime.now()
# Iteration numbers
linear_its = solver.snes.getLinearSolveIterations()
nonlinear_its = solver.snes.getIterationNumber()
time = (end-start).total_seconds() / 60
if mesh.comm.rank == 0:
if nonlinear_its == 0:
nonlinear_its = 1
print(GREEN % ("Time taken: %.2f min in %d nonlinear iterations, %d linear iterations (%.2f Krylov iters per Newton step)"
% (time, nonlinear_its, linear_its, linear_its/float(nonlinear_its))), flush=True)
print("%.2f @ %d @ %d @ %.2f" % (time, nonlinear_its, linear_its, linear_its/float(nonlinear_its)), flush=True)
(u, p, B, E) = z.split()
# Make sure that average of p is 0
pintegral = assemble(p*dx)
p.assign(p - Constant(pintegral/area))
B.rename("MagneticField")
u.rename("VelocityField")
p.rename("Pressure")
E.rename("ElectricFieldf")
# Compute divergence of u and B
norm_div_u = sqrt(assemble(inner(div(u), div(u))*dx))
norm_div_B = sqrt(assemble(inner(div(B), div(B))*dx))
if mesh.comm.rank == 0:
print("||div(u)||_L^2 = %s" % norm_div_u, flush=True)
print("||div(B)||_L^2 = %s" % norm_div_B, flush=True)
if solution_known:
B_ex_ = interpolate(B_ex, B.function_space())
u_ex_ = interpolate(u_ex, u.function_space())
p_ex_ = interpolate(p_ex, p.function_space())
E_ex_ = interpolate(E_ex, E.function_space())
B_ex_.rename("ExactSolutionB")
u_ex_.rename("ExactSolutionu")
p_ex_.rename("ExactSolutionp")
E_ex_.rename("ExactSolutionE")
# Compute error for MMS
error_u = errornorm(u_ex, u, 'L2')
error_B = errornorm(B_ex, B, 'L2')
error_E = errornorm(E_ex, E, 'L2')
error_p = errornorm(p_ex, p, 'L2')
if mesh.comm.rank == 0:
print("Error ||u_ex - u||_L^2 = %s" % error_u, flush=True)
print("Error ||p_ex - p||_L^2 = %s" % error_p, flush=True)
print("Error ||B_ex - B||_L^2 = %s" % error_B, flush=True)
print("Error ||E_ex - E||_L^2 = %s" % error_E, flush=True)
# Write errors to file
f = open("error.txt", 'a+')
f.write("%s,%s,%s,%s\n" % (error_u, error_p, error_B, error_E))
f.close()
# Save plots of solution
if output:
File("output/mhd.pvd").write(u, u_ex_, p, p_ex_, B, B_ex_, E, E_ex_)
else:
info_dict = {
"Re": re,
"Rem": rem,
"S": s,
"krylov/nonlin": linear_its/nonlinear_its,
"nonlinear_iter": nonlinear_its,
}
if output:
File("output/mhd.pvd").write(u, p, B, E, time=float(rem))
message(BLUE % info_dict)
# Write iteration numbers to file
if mesh.comm.rank == 0:
dir = 'results/results'+str(linearisation)+str(testproblem)+'/'
if not os.path.exists(dir):
os.mkdir(dir)
f = open(dir+str(float(ind1))+str(float(ind2))+'.txt', 'w+')
f.write("({0:2.0f}){1:4.1f}".format(float(info_dict["nonlinear_iter"]), float(info_dict["krylov/nonlin"])))
f.close()
z_last_u.assign(u)
# Create Checkpoints of solution
chk = DumbCheckpoint("dump/"+str(float(ind2))+str(linearisation)+str(testproblem), mode=FILE_CREATE)
chk.store(z)
# Loop over parameters