1 static char help[] = "This example illustrates the use of PCBDDC/FETI-DP with 2D/3D DMDA.\n\
2 It solves the constant coefficient Poisson problem or the Elasticity problem \n\
3 on a uniform grid of [0,cells_x] x [0,cells_y] x [0,cells_z]\n\n";
4
5 /* Contributed by Wim Vanroose <wim@vanroo.se> */
6
7 #include <petscksp.h>
8 #include <petscpc.h>
9 #include <petscdm.h>
10 #include <petscdmda.h>
11 #include <petscdmplex.h>
12
13 static PetscScalar poiss_1D_emat[] = {1.0000000000000000e+00, -1.0000000000000000e+00, -1.0000000000000000e+00, 1.0000000000000000e+00};
14 static PetscScalar poiss_2D_emat[] = {6.6666666666666674e-01, -1.6666666666666666e-01, -1.6666666666666666e-01, -3.3333333333333337e-01, -1.6666666666666666e-01, 6.6666666666666674e-01, -3.3333333333333337e-01, -1.6666666666666666e-01,
15 -1.6666666666666666e-01, -3.3333333333333337e-01, 6.6666666666666674e-01, -1.6666666666666666e-01, -3.3333333333333337e-01, -1.6666666666666666e-01, -1.6666666666666666e-01, 6.6666666666666674e-01};
16 static PetscScalar poiss_3D_emat[] = {3.3333333333333348e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, -8.3333333333333343e-02, 0.0000000000000000e+00, -8.3333333333333343e-02, -8.3333333333333343e-02, -8.3333333333333356e-02,
17 0.0000000000000000e+00, 3.3333333333333337e-01, -8.3333333333333343e-02, 0.0000000000000000e+00, -8.3333333333333343e-02, 0.0000000000000000e+00, -8.3333333333333356e-02, -8.3333333333333343e-02,
18 0.0000000000000000e+00, -8.3333333333333343e-02, 3.3333333333333337e-01, 0.0000000000000000e+00, -8.3333333333333343e-02, -8.3333333333333356e-02, 0.0000000000000000e+00, -8.3333333333333343e-02,
19 -8.3333333333333343e-02, 0.0000000000000000e+00, 0.0000000000000000e+00, 3.3333333333333348e-01, -8.3333333333333356e-02, -8.3333333333333343e-02, -8.3333333333333343e-02, 0.0000000000000000e+00,
20 0.0000000000000000e+00, -8.3333333333333343e-02, -8.3333333333333343e-02, -8.3333333333333356e-02, 3.3333333333333337e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, -8.3333333333333343e-02,
21 -8.3333333333333343e-02, 0.0000000000000000e+00, -8.3333333333333356e-02, -8.3333333333333343e-02, 0.0000000000000000e+00, 3.3333333333333337e-01, -8.3333333333333343e-02, 0.0000000000000000e+00,
22 -8.3333333333333343e-02, -8.3333333333333356e-02, 0.0000000000000000e+00, -8.3333333333333343e-02, 0.0000000000000000e+00, -8.3333333333333343e-02, 3.3333333333333337e-01, 0.0000000000000000e+00,
23 -8.3333333333333356e-02, -8.3333333333333343e-02, -8.3333333333333343e-02, 0.0000000000000000e+00, -8.3333333333333343e-02, 0.0000000000000000e+00, 0.0000000000000000e+00, 3.3333333333333337e-01};
24 static PetscScalar elast_1D_emat[] = {3.0000000000000000e+00, -3.0000000000000000e+00, -3.0000000000000000e+00, 3.0000000000000000e+00};
25 static PetscScalar elast_2D_emat[] = {1.3333333333333335e+00, 5.0000000000000000e-01, -8.3333333333333337e-01, 0.0000000000000000e+00, 1.6666666666666671e-01, 0.0000000000000000e+00, -6.6666666666666674e-01, -5.0000000000000000e-01,
26 5.0000000000000000e-01, 1.3333333333333335e+00, 0.0000000000000000e+00, 1.6666666666666671e-01, 0.0000000000000000e+00, -8.3333333333333337e-01, -5.0000000000000000e-01, -6.6666666666666674e-01,
27 -8.3333333333333337e-01, 0.0000000000000000e+00, 1.3333333333333335e+00, -5.0000000000000000e-01, -6.6666666666666674e-01, 5.0000000000000000e-01, 1.6666666666666674e-01, 0.0000000000000000e+00,
28 0.0000000000000000e+00, 1.6666666666666671e-01, -5.0000000000000000e-01, 1.3333333333333335e+00, 5.0000000000000000e-01, -6.6666666666666674e-01, 0.0000000000000000e+00, -8.3333333333333337e-01,
29 1.6666666666666671e-01, 0.0000000000000000e+00, -6.6666666666666674e-01, 5.0000000000000000e-01, 1.3333333333333335e+00, -5.0000000000000000e-01, -8.3333333333333337e-01, 0.0000000000000000e+00,
30 0.0000000000000000e+00, -8.3333333333333337e-01, 5.0000000000000000e-01, -6.6666666666666674e-01, -5.0000000000000000e-01, 1.3333333333333335e+00, 0.0000000000000000e+00, 1.6666666666666674e-01,
31 -6.6666666666666674e-01, -5.0000000000000000e-01, 1.6666666666666674e-01, 0.0000000000000000e+00, -8.3333333333333337e-01, 0.0000000000000000e+00, 1.3333333333333335e+00, 5.0000000000000000e-01,
32 -5.0000000000000000e-01, -6.6666666666666674e-01, 0.0000000000000000e+00, -8.3333333333333337e-01, 0.0000000000000000e+00, 1.6666666666666674e-01, 5.0000000000000000e-01, 1.3333333333333335e+00};
33 static PetscScalar elast_3D_emat[] =
34 {5.5555555555555558e-01, 1.6666666666666666e-01, 1.6666666666666666e-01, -2.2222222222222232e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, 1.1111111111111113e-01, 0.0000000000000000e+00, 8.3333333333333356e-02,
35 -1.9444444444444442e-01, -1.6666666666666669e-01, 0.0000000000000000e+00, 1.1111111111111112e-01, 8.3333333333333356e-02, 0.0000000000000000e+00, -1.9444444444444445e-01, 0.0000000000000000e+00, -1.6666666666666669e-01,
36 -2.7777777777777769e-02, 0.0000000000000000e+00, 0.0000000000000000e+00, -1.3888888888888887e-01, -8.3333333333333356e-02, -8.3333333333333356e-02, 1.6666666666666666e-01, 5.5555555555555558e-01, 1.6666666666666666e-01,
37 0.0000000000000000e+00, 1.1111111111111113e-01, 8.3333333333333356e-02, 0.0000000000000000e+00, -2.2222222222222232e-01, 0.0000000000000000e+00, -1.6666666666666669e-01, -1.9444444444444442e-01, 0.0000000000000000e+00,
38 8.3333333333333356e-02, 1.1111111111111112e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, -2.7777777777777769e-02, 0.0000000000000000e+00, 0.0000000000000000e+00, -1.9444444444444445e-01, -1.6666666666666669e-01,
39 -8.3333333333333356e-02, -1.3888888888888887e-01, -8.3333333333333356e-02, 1.6666666666666666e-01, 1.6666666666666666e-01, 5.5555555555555558e-01, 0.0000000000000000e+00, 8.3333333333333356e-02, 1.1111111111111112e-01,
40 8.3333333333333356e-02, 0.0000000000000000e+00, 1.1111111111111112e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, -2.7777777777777769e-02, 0.0000000000000000e+00, 0.0000000000000000e+00, -2.2222222222222229e-01,
41 -1.6666666666666669e-01, 0.0000000000000000e+00, -1.9444444444444445e-01, 0.0000000000000000e+00, -1.6666666666666669e-01, -1.9444444444444445e-01, -8.3333333333333356e-02, -8.3333333333333356e-02, -1.3888888888888887e-01,
42 -2.2222222222222232e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, 5.5555555555555558e-01, -1.6666666666666666e-01, -1.6666666666666666e-01, -1.9444444444444442e-01, 1.6666666666666669e-01, 0.0000000000000000e+00,
43 1.1111111111111113e-01, 0.0000000000000000e+00, -8.3333333333333356e-02, -1.9444444444444445e-01, 0.0000000000000000e+00, 1.6666666666666669e-01, 1.1111111111111113e-01, -8.3333333333333356e-02, 0.0000000000000000e+00,
44 -1.3888888888888887e-01, 8.3333333333333356e-02, 8.3333333333333356e-02, -2.7777777777777769e-02, 0.0000000000000000e+00, 0.0000000000000000e+00, 0.0000000000000000e+00, 1.1111111111111113e-01, 8.3333333333333356e-02,
45 -1.6666666666666666e-01, 5.5555555555555558e-01, 1.6666666666666669e-01, 1.6666666666666669e-01, -1.9444444444444442e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, -2.2222222222222229e-01, 0.0000000000000000e+00,
46 0.0000000000000000e+00, -2.7777777777777769e-02, 0.0000000000000000e+00, -8.3333333333333356e-02, 1.1111111111111112e-01, 0.0000000000000000e+00, 8.3333333333333356e-02, -1.3888888888888887e-01, -8.3333333333333356e-02,
47 0.0000000000000000e+00, -1.9444444444444448e-01, -1.6666666666666666e-01, 0.0000000000000000e+00, 8.3333333333333356e-02, 1.1111111111111112e-01, -1.6666666666666666e-01, 1.6666666666666669e-01, 5.5555555555555558e-01,
48 0.0000000000000000e+00, 0.0000000000000000e+00, -2.7777777777777769e-02, -8.3333333333333356e-02, 0.0000000000000000e+00, 1.1111111111111112e-01, 1.6666666666666669e-01, 0.0000000000000000e+00, -1.9444444444444445e-01,
49 0.0000000000000000e+00, 0.0000000000000000e+00, -2.2222222222222227e-01, 8.3333333333333356e-02, -8.3333333333333356e-02, -1.3888888888888887e-01, 0.0000000000000000e+00, -1.6666666666666666e-01, -1.9444444444444448e-01,
50 1.1111111111111113e-01, 0.0000000000000000e+00, 8.3333333333333356e-02, -1.9444444444444442e-01, 1.6666666666666669e-01, 0.0000000000000000e+00, 5.5555555555555569e-01, -1.6666666666666666e-01, 1.6666666666666669e-01,
51 -2.2222222222222229e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, -2.7777777777777769e-02, 0.0000000000000000e+00, 0.0000000000000000e+00, -1.3888888888888887e-01, 8.3333333333333356e-02, -8.3333333333333356e-02,
52 1.1111111111111112e-01, -8.3333333333333343e-02, 0.0000000000000000e+00, -1.9444444444444448e-01, 0.0000000000000000e+00, -1.6666666666666669e-01, 0.0000000000000000e+00, -2.2222222222222232e-01, 0.0000000000000000e+00,
53 1.6666666666666669e-01, -1.9444444444444442e-01, 0.0000000000000000e+00, -1.6666666666666666e-01, 5.5555555555555558e-01, -1.6666666666666669e-01, 0.0000000000000000e+00, 1.1111111111111113e-01, -8.3333333333333343e-02,
54 0.0000000000000000e+00, -1.9444444444444445e-01, 1.6666666666666669e-01, 8.3333333333333356e-02, -1.3888888888888887e-01, 8.3333333333333356e-02, -8.3333333333333343e-02, 1.1111111111111113e-01, 0.0000000000000000e+00,
55 0.0000000000000000e+00, -2.7777777777777769e-02, 0.0000000000000000e+00, 8.3333333333333356e-02, 0.0000000000000000e+00, 1.1111111111111112e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, -2.7777777777777769e-02,
56 1.6666666666666669e-01, -1.6666666666666669e-01, 5.5555555555555558e-01, 0.0000000000000000e+00, -8.3333333333333343e-02, 1.1111111111111112e-01, 0.0000000000000000e+00, 1.6666666666666669e-01, -1.9444444444444445e-01,
57 -8.3333333333333356e-02, 8.3333333333333356e-02, -1.3888888888888887e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, -2.2222222222222227e-01, -1.6666666666666669e-01, 0.0000000000000000e+00, -1.9444444444444448e-01,
58 -1.9444444444444442e-01, -1.6666666666666669e-01, 0.0000000000000000e+00, 1.1111111111111113e-01, 0.0000000000000000e+00, -8.3333333333333356e-02, -2.2222222222222229e-01, 0.0000000000000000e+00, 0.0000000000000000e+00,
59 5.5555555555555558e-01, 1.6666666666666669e-01, -1.6666666666666666e-01, -1.3888888888888887e-01, -8.3333333333333356e-02, 8.3333333333333356e-02, -2.7777777777777769e-02, 0.0000000000000000e+00, 0.0000000000000000e+00,
60 -1.9444444444444448e-01, 0.0000000000000000e+00, 1.6666666666666669e-01, 1.1111111111111112e-01, 8.3333333333333343e-02, 0.0000000000000000e+00, -1.6666666666666669e-01, -1.9444444444444442e-01, 0.0000000000000000e+00,
61 0.0000000000000000e+00, -2.2222222222222229e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, 1.1111111111111113e-01, -8.3333333333333343e-02, 1.6666666666666669e-01, 5.5555555555555558e-01, -1.6666666666666669e-01,
62 -8.3333333333333356e-02, -1.3888888888888887e-01, 8.3333333333333356e-02, 0.0000000000000000e+00, -1.9444444444444448e-01, 1.6666666666666669e-01, 0.0000000000000000e+00, -2.7777777777777769e-02, 0.0000000000000000e+00,
63 8.3333333333333343e-02, 1.1111111111111112e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, 0.0000000000000000e+00, -2.7777777777777769e-02, -8.3333333333333356e-02, 0.0000000000000000e+00, 1.1111111111111112e-01,
64 0.0000000000000000e+00, -8.3333333333333343e-02, 1.1111111111111112e-01, -1.6666666666666666e-01, -1.6666666666666669e-01, 5.5555555555555558e-01, 8.3333333333333356e-02, 8.3333333333333356e-02, -1.3888888888888887e-01,
65 0.0000000000000000e+00, 1.6666666666666669e-01, -1.9444444444444448e-01, 1.6666666666666669e-01, 0.0000000000000000e+00, -1.9444444444444448e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, -2.2222222222222227e-01,
66 1.1111111111111112e-01, 8.3333333333333356e-02, 0.0000000000000000e+00, -1.9444444444444445e-01, 0.0000000000000000e+00, 1.6666666666666669e-01, -2.7777777777777769e-02, 0.0000000000000000e+00, 0.0000000000000000e+00,
67 -1.3888888888888887e-01, -8.3333333333333356e-02, 8.3333333333333356e-02, 5.5555555555555569e-01, 1.6666666666666669e-01, -1.6666666666666669e-01, -2.2222222222222227e-01, 0.0000000000000000e+00, 0.0000000000000000e+00,
68 1.1111111111111112e-01, 0.0000000000000000e+00, -8.3333333333333343e-02, -1.9444444444444448e-01, -1.6666666666666669e-01, 0.0000000000000000e+00, 8.3333333333333356e-02, 1.1111111111111112e-01, 0.0000000000000000e+00,
69 0.0000000000000000e+00, -2.7777777777777769e-02, 0.0000000000000000e+00, 0.0000000000000000e+00, -1.9444444444444445e-01, 1.6666666666666669e-01, -8.3333333333333356e-02, -1.3888888888888887e-01, 8.3333333333333356e-02,
70 1.6666666666666669e-01, 5.5555555555555558e-01, -1.6666666666666669e-01, 0.0000000000000000e+00, 1.1111111111111112e-01, -8.3333333333333343e-02, 0.0000000000000000e+00, -2.2222222222222227e-01, 0.0000000000000000e+00,
71 -1.6666666666666669e-01, -1.9444444444444448e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, 0.0000000000000000e+00, -2.2222222222222229e-01, 1.6666666666666669e-01, 0.0000000000000000e+00, -1.9444444444444445e-01,
72 0.0000000000000000e+00, 1.6666666666666669e-01, -1.9444444444444445e-01, 8.3333333333333356e-02, 8.3333333333333356e-02, -1.3888888888888887e-01, -1.6666666666666669e-01, -1.6666666666666669e-01, 5.5555555555555558e-01,
73 0.0000000000000000e+00, -8.3333333333333343e-02, 1.1111111111111113e-01, -8.3333333333333343e-02, 0.0000000000000000e+00, 1.1111111111111113e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, -2.7777777777777769e-02,
74 -1.9444444444444445e-01, 0.0000000000000000e+00, -1.6666666666666669e-01, 1.1111111111111113e-01, -8.3333333333333356e-02, 0.0000000000000000e+00, -1.3888888888888887e-01, 8.3333333333333356e-02, -8.3333333333333356e-02,
75 -2.7777777777777769e-02, 0.0000000000000000e+00, 0.0000000000000000e+00, -2.2222222222222227e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, 5.5555555555555558e-01, -1.6666666666666669e-01, 1.6666666666666669e-01,
76 -1.9444444444444448e-01, 1.6666666666666669e-01, 0.0000000000000000e+00, 1.1111111111111112e-01, 0.0000000000000000e+00, 8.3333333333333343e-02, 0.0000000000000000e+00, -2.7777777777777769e-02, 0.0000000000000000e+00,
77 -8.3333333333333356e-02, 1.1111111111111112e-01, 0.0000000000000000e+00, 8.3333333333333356e-02, -1.3888888888888887e-01, 8.3333333333333356e-02, 0.0000000000000000e+00, -1.9444444444444448e-01, 1.6666666666666669e-01,
78 0.0000000000000000e+00, 1.1111111111111112e-01, -8.3333333333333343e-02, -1.6666666666666669e-01, 5.5555555555555558e-01, -1.6666666666666666e-01, 1.6666666666666669e-01, -1.9444444444444448e-01, 0.0000000000000000e+00,
79 0.0000000000000000e+00, -2.2222222222222227e-01, 0.0000000000000000e+00, -1.6666666666666669e-01, 0.0000000000000000e+00, -1.9444444444444445e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, -2.2222222222222227e-01,
80 -8.3333333333333356e-02, 8.3333333333333356e-02, -1.3888888888888887e-01, 0.0000000000000000e+00, 1.6666666666666669e-01, -1.9444444444444448e-01, 0.0000000000000000e+00, -8.3333333333333343e-02, 1.1111111111111113e-01,
81 1.6666666666666669e-01, -1.6666666666666666e-01, 5.5555555555555558e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, -2.7777777777777769e-02, 8.3333333333333343e-02, 0.0000000000000000e+00, 1.1111111111111113e-01,
82 -2.7777777777777769e-02, 0.0000000000000000e+00, 0.0000000000000000e+00, -1.3888888888888887e-01, 8.3333333333333356e-02, 8.3333333333333356e-02, 1.1111111111111112e-01, -8.3333333333333343e-02, 0.0000000000000000e+00,
83 -1.9444444444444448e-01, 0.0000000000000000e+00, 1.6666666666666669e-01, 1.1111111111111112e-01, 0.0000000000000000e+00, -8.3333333333333343e-02, -1.9444444444444448e-01, 1.6666666666666669e-01, 0.0000000000000000e+00,
84 5.5555555555555558e-01, -1.6666666666666669e-01, -1.6666666666666669e-01, -2.2222222222222227e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, 0.0000000000000000e+00, -1.9444444444444445e-01, -1.6666666666666669e-01,
85 8.3333333333333356e-02, -1.3888888888888887e-01, -8.3333333333333356e-02, -8.3333333333333343e-02, 1.1111111111111113e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, -2.7777777777777769e-02, 0.0000000000000000e+00,
86 0.0000000000000000e+00, -2.2222222222222227e-01, 0.0000000000000000e+00, 1.6666666666666669e-01, -1.9444444444444448e-01, 0.0000000000000000e+00, -1.6666666666666669e-01, 5.5555555555555558e-01, 1.6666666666666669e-01,
87 0.0000000000000000e+00, 1.1111111111111112e-01, 8.3333333333333343e-02, 0.0000000000000000e+00, -1.6666666666666669e-01, -1.9444444444444445e-01, 8.3333333333333356e-02, -8.3333333333333356e-02, -1.3888888888888887e-01,
88 0.0000000000000000e+00, 0.0000000000000000e+00, -2.2222222222222227e-01, 1.6666666666666669e-01, 0.0000000000000000e+00, -1.9444444444444448e-01, -8.3333333333333343e-02, 0.0000000000000000e+00, 1.1111111111111113e-01,
89 0.0000000000000000e+00, 0.0000000000000000e+00, -2.7777777777777769e-02, -1.6666666666666669e-01, 1.6666666666666669e-01, 5.5555555555555558e-01, 0.0000000000000000e+00, 8.3333333333333343e-02, 1.1111111111111113e-01,
90 -1.3888888888888887e-01, -8.3333333333333356e-02, -8.3333333333333356e-02, -2.7777777777777769e-02, 0.0000000000000000e+00, 0.0000000000000000e+00, -1.9444444444444448e-01, 0.0000000000000000e+00, -1.6666666666666669e-01,
91 1.1111111111111112e-01, 8.3333333333333343e-02, 0.0000000000000000e+00, -1.9444444444444448e-01, -1.6666666666666669e-01, 0.0000000000000000e+00, 1.1111111111111112e-01, 0.0000000000000000e+00, 8.3333333333333343e-02,
92 -2.2222222222222227e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, 5.5555555555555558e-01, 1.6666666666666669e-01, 1.6666666666666669e-01, -8.3333333333333356e-02, -1.3888888888888887e-01, -8.3333333333333356e-02,
93 0.0000000000000000e+00, -1.9444444444444448e-01, -1.6666666666666666e-01, 0.0000000000000000e+00, -2.7777777777777769e-02, 0.0000000000000000e+00, 8.3333333333333343e-02, 1.1111111111111112e-01, 0.0000000000000000e+00,
94 -1.6666666666666669e-01, -1.9444444444444448e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, -2.2222222222222227e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, 1.1111111111111112e-01, 8.3333333333333343e-02,
95 1.6666666666666669e-01, 5.5555555555555558e-01, 1.6666666666666669e-01, -8.3333333333333356e-02, -8.3333333333333356e-02, -1.3888888888888887e-01, 0.0000000000000000e+00, -1.6666666666666666e-01, -1.9444444444444448e-01,
96 -1.6666666666666669e-01, 0.0000000000000000e+00, -1.9444444444444448e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, -2.2222222222222227e-01, 0.0000000000000000e+00, 0.0000000000000000e+00, -2.7777777777777769e-02,
97 8.3333333333333343e-02, 0.0000000000000000e+00, 1.1111111111111113e-01, 0.0000000000000000e+00, 8.3333333333333343e-02, 1.1111111111111113e-01, 1.6666666666666669e-01, 1.6666666666666669e-01, 5.5555555555555558e-01};
98
99 typedef enum {
100 PDE_POISSON,
101 PDE_ELASTICITY
102 } PDEType;
103
104 typedef struct {
105 PDEType pde;
106 PetscInt dim;
107 PetscInt dof;
108 PetscInt cells[3];
109 PetscBool useglobal;
110 PetscBool multi_element;
111 PetscBool dirbc;
112 PetscBool per[3];
113 PetscBool test;
114 PetscScalar *elemMat;
115 PetscBool use_composite_pc;
116 PetscBool random_initial_guess;
117 PetscBool random_real;
118 } AppCtx;
119
ProcessOptions(MPI_Comm comm,AppCtx * options)120 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
121 {
122 const char *pdeTypes[2] = {"Poisson", "Elasticity"};
123 PetscInt n, pde;
124
125 PetscFunctionBeginUser;
126 options->pde = PDE_POISSON;
127 options->elemMat = NULL;
128 options->dim = 1;
129 options->cells[0] = 8;
130 options->cells[1] = 6;
131 options->cells[2] = 4;
132 options->useglobal = PETSC_FALSE;
133 options->multi_element = PETSC_FALSE;
134 options->dirbc = PETSC_TRUE;
135 options->test = PETSC_FALSE;
136 options->per[0] = PETSC_FALSE;
137 options->per[1] = PETSC_FALSE;
138 options->per[2] = PETSC_FALSE;
139 options->use_composite_pc = PETSC_FALSE;
140 options->random_initial_guess = PETSC_FALSE;
141 options->random_real = PETSC_FALSE;
142
143 PetscOptionsBegin(comm, NULL, "Problem Options", NULL);
144 pde = options->pde;
145 PetscCall(PetscOptionsEList("-pde_type", "The PDE type", __FILE__, pdeTypes, 2, pdeTypes[options->pde], &pde, NULL));
146 options->pde = (PDEType)pde;
147 PetscCall(PetscOptionsInt("-dim", "The topological mesh dimension", __FILE__, options->dim, &options->dim, NULL));
148 PetscCall(PetscOptionsIntArray("-cells", "The mesh division", __FILE__, options->cells, (n = 3, &n), NULL));
149 PetscCall(PetscOptionsBoolArray("-periodicity", "The mesh periodicity", __FILE__, options->per, (n = 3, &n), NULL));
150 PetscCall(PetscOptionsBool("-use_global", "Test MatSetValues", __FILE__, options->useglobal, &options->useglobal, NULL));
151 PetscCall(PetscOptionsBool("-multi_element", "Use multi-element BDDC", __FILE__, options->multi_element, &options->multi_element, NULL));
152 PetscCall(PetscOptionsBool("-dirichlet", "Use dirichlet BC", __FILE__, options->dirbc, &options->dirbc, NULL));
153 PetscCall(PetscOptionsBool("-use_composite_pc", "Multiplicative composite with BDDC + Richardson/Jacobi", __FILE__, options->use_composite_pc, &options->use_composite_pc, NULL));
154 PetscCall(PetscOptionsBool("-random_initial_guess", "Solve A x = 0 with random initial guess, instead of A x = b with random b", __FILE__, options->random_initial_guess, &options->random_initial_guess, NULL));
155 PetscCall(PetscOptionsBool("-random_real", "Use real-valued b (or x, if -random_initial_guess) instead of default scalar type", __FILE__, options->random_real, &options->random_real, NULL));
156 PetscOptionsEnd();
157
158 for (n = options->dim; n < 3; n++) options->cells[n] = 0;
159 if (options->per[0]) options->dirbc = PETSC_FALSE;
160
161 /* element matrices */
162 switch (options->pde) {
163 case PDE_ELASTICITY:
164 options->dof = options->dim;
165 switch (options->dim) {
166 case 1:
167 options->elemMat = elast_1D_emat;
168 break;
169 case 2:
170 options->elemMat = elast_2D_emat;
171 break;
172 case 3:
173 options->elemMat = elast_3D_emat;
174 break;
175 default:
176 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Unsupported dimension %" PetscInt_FMT, options->dim);
177 }
178 break;
179 case PDE_POISSON:
180 options->dof = 1;
181 switch (options->dim) {
182 case 1:
183 options->elemMat = poiss_1D_emat;
184 break;
185 case 2:
186 options->elemMat = poiss_2D_emat;
187 break;
188 case 3:
189 options->elemMat = poiss_3D_emat;
190 break;
191 default:
192 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Unsupported dimension %" PetscInt_FMT, options->dim);
193 }
194 break;
195 default:
196 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Unsupported PDE %d", options->pde);
197 }
198 PetscFunctionReturn(PETSC_SUCCESS);
199 }
200
main(int argc,char ** args)201 int main(int argc, char **args)
202 {
203 AppCtx user;
204 KSP ksp;
205 PC pc;
206 Mat A;
207 DM da;
208 Vec x, b, xcoor, xcoorl;
209 IS zero;
210 ISLocalToGlobalMapping map;
211 MatNullSpace nullsp = NULL;
212 PetscInt i;
213 PetscInt nel, nen; /* Number of elements & element nodes */
214 const PetscInt *e_loc; /* Local indices of element nodes (in local element order) */
215 PetscInt *e_glo = NULL; /* Global indices of element nodes (in local element order) */
216 PetscInt nodes[3];
217 PetscBool ismatis, flg;
218 PetscLogStage stages[2];
219
220 PetscFunctionBeginUser;
221 PetscCall(PetscInitialize(&argc, &args, NULL, help));
222 PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
223 for (i = 0; i < 3; i++) nodes[i] = user.cells[i] + !user.per[i];
224 switch (user.dim) {
225 case 3:
226 PetscCall(DMDACreate3d(PETSC_COMM_WORLD, user.per[0] ? DM_BOUNDARY_PERIODIC : DM_BOUNDARY_NONE, user.per[1] ? DM_BOUNDARY_PERIODIC : DM_BOUNDARY_NONE, user.per[2] ? DM_BOUNDARY_PERIODIC : DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, nodes[0], nodes[1], nodes[2], PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE,
227 user.dof, 1, NULL, NULL, NULL, &da));
228 break;
229 case 2:
230 PetscCall(DMDACreate2d(PETSC_COMM_WORLD, user.per[0] ? DM_BOUNDARY_PERIODIC : DM_BOUNDARY_NONE, user.per[1] ? DM_BOUNDARY_PERIODIC : DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, nodes[0], nodes[1], PETSC_DECIDE, PETSC_DECIDE, user.dof, 1, NULL, NULL, &da));
231 break;
232 case 1:
233 PetscCall(DMDACreate1d(PETSC_COMM_WORLD, user.per[0] ? DM_BOUNDARY_PERIODIC : DM_BOUNDARY_NONE, nodes[0], user.dof, 1, NULL, &da));
234 break;
235 default:
236 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Unsupported dimension %" PetscInt_FMT, user.dim);
237 }
238
239 PetscCall(PetscLogStageRegister("KSPSetUp", &stages[0]));
240 PetscCall(PetscLogStageRegister("KSPSolve", &stages[1]));
241
242 PetscCall(DMSetMatType(da, MATIS));
243 PetscCall(DMSetFromOptions(da));
244 PetscCall(DMDASetElementType(da, DMDA_ELEMENT_Q1));
245 PetscCall(DMSetUp(da));
246 {
247 PetscInt M, N, P;
248 PetscCall(DMDAGetInfo(da, 0, &M, &N, &P, 0, 0, 0, 0, 0, 0, 0, 0, 0));
249 switch (user.dim) {
250 case 3:
251 user.cells[2] = P - !user.per[2]; /* fall through */
252 case 2:
253 user.cells[1] = N - !user.per[1]; /* fall through */
254 case 1:
255 user.cells[0] = M - !user.per[0];
256 break;
257 default:
258 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Unsupported dimension %" PetscInt_FMT, user.dim);
259 }
260 }
261 PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0 * user.cells[0], 0.0, 1.0 * user.cells[1], 0.0, 1.0 * user.cells[2]));
262 PetscCall(DMGetCoordinates(da, &xcoor));
263
264 PetscCall(DMCreateMatrix(da, &A));
265 PetscCall(MatSetFromOptions(A));
266 PetscCall(DMGetLocalToGlobalMapping(da, &map));
267 PetscCall(DMDAGetElements(da, &nel, &nen, &e_loc));
268 if (user.useglobal) {
269 PetscCall(PetscMalloc1(nel * nen, &e_glo));
270 PetscCall(ISLocalToGlobalMappingApplyBlock(map, nen * nel, e_loc, e_glo));
271 }
272
273 if (user.multi_element) {
274 ISLocalToGlobalMapping mapn;
275 PetscInt *el_glo = NULL, m, n, M, N, *el_sizes;
276 Mat lA;
277
278 PetscCall(PetscMalloc1(nel * nen, &el_glo));
279 PetscCall(ISLocalToGlobalMappingApplyBlock(map, nen * nel, e_loc, el_glo));
280 PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)map), user.dof, nen * nel, el_glo, PETSC_OWN_POINTER, &mapn));
281 PetscCall(MatGetLocalSize(A, &m, &n));
282 PetscCall(MatGetSize(A, &M, &N));
283 PetscCall(MatDestroy(&A));
284 PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
285 PetscCall(MatSetSizes(A, m, n, M, N));
286 PetscCall(MatSetBlockSize(A, user.dof));
287 PetscCall(MatSetType(A, MATIS));
288 PetscCall(MatISSetAllowRepeated(A, PETSC_TRUE));
289 PetscCall(MatSetLocalToGlobalMapping(A, mapn, mapn));
290 PetscCall(MatISSetPreallocation(A, user.dof * nen, NULL, user.dof * nen, NULL));
291 PetscCall(ISLocalToGlobalMappingViewFromOptions(mapn, NULL, "-multi_view"));
292 PetscCall(ISLocalToGlobalMappingDestroy(&mapn));
293
294 /* The information set with MatSetVariableBlockSizes on the local mat
295 can be used to detect the local elements instead of having to analyze
296 the sparsity pattern of the local matrix */
297 PetscCall(MatISGetLocalMat(A, &lA));
298 PetscCall(PetscMalloc1(nel, &el_sizes));
299 for (i = 0; i < nel; i++) el_sizes[i] = user.dof * nen;
300 PetscCall(MatSetVariableBlockSizes(lA, nel, el_sizes));
301 PetscCall(PetscFree(el_sizes));
302 }
303
304 /* we reorder the indices since the element matrices are given in lexicographic order,
305 whereas the elements indices returned by DMDAGetElements follow the usual FEM ordering
306 i.e., element matrices DMDA ordering
307 2---3 3---2
308 / / / /
309 0---1 0---1
310 */
311 for (i = 0; i < nel; ++i) {
312 PetscInt ord[8] = {0, 1, 3, 2, 4, 5, 7, 6};
313 PetscInt j, idxs[8];
314
315 PetscCheck(nen <= 8, PETSC_COMM_WORLD, PETSC_ERR_SUP, "Not coded");
316 if (!user.useglobal) {
317 if (user.multi_element) {
318 for (j = 0; j < nen; j++) idxs[j] = i * nen + ord[j];
319 } else {
320 for (j = 0; j < nen; j++) idxs[j] = e_loc[i * nen + ord[j]];
321 }
322 PetscCall(MatSetValuesBlockedLocal(A, nen, idxs, nen, idxs, user.elemMat, ADD_VALUES));
323 } else {
324 for (j = 0; j < nen; j++) idxs[j] = e_glo[i * nen + ord[j]];
325 PetscCall(MatSetValuesBlocked(A, nen, idxs, nen, idxs, user.elemMat, ADD_VALUES));
326 }
327 }
328 PetscCall(DMDARestoreElements(da, &nel, &nen, &e_loc));
329 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
330 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
331 PetscCall(MatViewFromOptions(A, NULL, "-A_mat_view"));
332 PetscCall(MatSetOption(A, MAT_SPD, PETSC_TRUE));
333 PetscCall(MatSetOption(A, MAT_SPD_ETERNAL, PETSC_TRUE));
334
335 /* Boundary conditions */
336 zero = NULL;
337 if (user.dirbc) { /* fix one side of DMDA */
338 Vec nat, glob;
339 PetscScalar *vals;
340 PetscInt n, *idx, j, st;
341
342 n = PetscGlobalRank ? 0 : (user.cells[1] + 1) * (user.cells[2] + 1);
343 PetscCall(ISCreateStride(PETSC_COMM_WORLD, n, 0, user.cells[0] + 1, &zero));
344 if (user.dof > 1) { /* zero all components */
345 const PetscInt *idx;
346 IS bzero;
347
348 PetscCall(ISGetIndices(zero, &idx));
349 PetscCall(ISCreateBlock(PETSC_COMM_WORLD, user.dof, n, idx, PETSC_COPY_VALUES, &bzero));
350 PetscCall(ISRestoreIndices(zero, &idx));
351 PetscCall(ISDestroy(&zero));
352 zero = bzero;
353 }
354 /* map indices from natural to global */
355 PetscCall(DMDACreateNaturalVector(da, &nat));
356 PetscCall(ISGetLocalSize(zero, &n));
357 PetscCall(PetscMalloc1(n, &vals));
358 for (i = 0; i < n; i++) vals[i] = 1.0;
359 PetscCall(ISGetIndices(zero, (const PetscInt **)&idx));
360 PetscCall(VecSetValues(nat, n, idx, vals, INSERT_VALUES));
361 PetscCall(ISRestoreIndices(zero, (const PetscInt **)&idx));
362 PetscCall(PetscFree(vals));
363 PetscCall(VecAssemblyBegin(nat));
364 PetscCall(VecAssemblyEnd(nat));
365 PetscCall(DMCreateGlobalVector(da, &glob));
366 PetscCall(DMDANaturalToGlobalBegin(da, nat, INSERT_VALUES, glob));
367 PetscCall(DMDANaturalToGlobalEnd(da, nat, INSERT_VALUES, glob));
368 PetscCall(VecDestroy(&nat));
369 PetscCall(ISDestroy(&zero));
370 PetscCall(VecGetLocalSize(glob, &n));
371 PetscCall(PetscMalloc1(n, &idx));
372 PetscCall(VecGetOwnershipRange(glob, &st, NULL));
373 PetscCall(VecGetArray(glob, &vals));
374 for (i = 0, j = 0; i < n; i++)
375 if (PetscRealPart(vals[i]) == 1.0) idx[j++] = i + st;
376 PetscCall(VecRestoreArray(glob, &vals));
377 PetscCall(VecDestroy(&glob));
378 PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, j, idx, PETSC_OWN_POINTER, &zero));
379 PetscCall(MatZeroRowsColumnsIS(A, zero, 1.0, NULL, NULL));
380 } else {
381 switch (user.pde) {
382 case PDE_POISSON:
383 PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, NULL, &nullsp));
384 break;
385 case PDE_ELASTICITY:
386 PetscCall(MatNullSpaceCreateRigidBody(xcoor, &nullsp));
387 break;
388 }
389 /* with periodic BC and Elasticity, just the displacements are in the nullspace
390 this is no harm since we eliminate all the components of the rhs */
391 PetscCall(MatSetNullSpace(A, nullsp));
392 }
393
394 PetscCall(PetscOptionsHasName(NULL, NULL, "-assembled_view", &flg));
395 if (flg) {
396 Mat AA;
397
398 PetscCall(MatConvert(A, MATAIJ, MAT_INITIAL_MATRIX, &AA));
399 PetscCall(MatViewFromOptions(AA, NULL, "-assembled_view"));
400 PetscCall(MatDestroy(&AA));
401 }
402
403 /* Attach near null space for elasticity */
404 if (user.pde == PDE_ELASTICITY) {
405 MatNullSpace nearnullsp;
406
407 PetscCall(MatNullSpaceCreateRigidBody(xcoor, &nearnullsp));
408 PetscCall(MatSetNearNullSpace(A, nearnullsp));
409 PetscCall(MatNullSpaceDestroy(&nearnullsp));
410 }
411
412 /* we may want to use MG for the local solvers: attach local nearnullspace to the local matrices */
413 PetscCall(DMGetCoordinatesLocal(da, &xcoorl));
414 PetscCall(PetscObjectTypeCompare((PetscObject)A, MATIS, &ismatis));
415 if (ismatis) {
416 MatNullSpace lnullsp = NULL;
417 Mat lA;
418
419 PetscCall(MatISGetLocalMat(A, &lA));
420 if (user.pde == PDE_ELASTICITY) {
421 Vec lc;
422 ISLocalToGlobalMapping l2l;
423 IS is;
424 const PetscScalar *a;
425 const PetscInt *idxs;
426 PetscInt n, bs;
427
428 /* when using a DMDA, the local matrices have an additional local-to-local map
429 that maps from the DA local ordering to the ordering induced by the elements */
430 PetscCall(MatGetLocalToGlobalMapping(lA, &l2l, NULL));
431 if (l2l) {
432 PetscCall(MatCreateVecs(lA, &lc, NULL));
433 PetscCall(VecSetLocalToGlobalMapping(lc, l2l));
434
435 PetscCall(VecSetOption(lc, VEC_IGNORE_NEGATIVE_INDICES, PETSC_TRUE));
436 PetscCall(VecGetLocalSize(xcoorl, &n));
437 PetscCall(VecGetBlockSize(xcoorl, &bs));
438 PetscCall(ISCreateStride(PETSC_COMM_SELF, n / bs, 0, 1, &is));
439 PetscCall(ISGetIndices(is, &idxs));
440 PetscCall(VecGetArrayRead(xcoorl, &a));
441 PetscCall(VecSetValuesBlockedLocal(lc, n / bs, idxs, a, INSERT_VALUES));
442 PetscCall(VecAssemblyBegin(lc));
443 PetscCall(VecAssemblyEnd(lc));
444 PetscCall(VecRestoreArrayRead(xcoorl, &a));
445 PetscCall(ISRestoreIndices(is, &idxs));
446 PetscCall(ISDestroy(&is));
447 PetscCall(MatNullSpaceCreateRigidBody(lc, &lnullsp));
448 PetscCall(VecDestroy(&lc));
449 }
450 } else if (user.pde == PDE_POISSON) {
451 PetscCall(MatNullSpaceCreate(PETSC_COMM_SELF, PETSC_TRUE, 0, NULL, &lnullsp));
452 }
453 PetscCall(MatSetNearNullSpace(lA, lnullsp));
454 PetscCall(MatNullSpaceDestroy(&lnullsp));
455 PetscCall(MatISRestoreLocalMat(A, &lA));
456 }
457
458 PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
459 PetscCall(KSPSetOperators(ksp, A, A));
460 PetscCall(KSPSetType(ksp, KSPCG));
461 PetscCall(KSPGetPC(ksp, &pc));
462 if (ismatis) {
463 if (user.use_composite_pc) {
464 PC pcksp, pcjacobi;
465 KSP ksprich;
466 PetscCall(PCSetType(pc, PCCOMPOSITE));
467 PetscCall(PCCompositeSetType(pc, PC_COMPOSITE_MULTIPLICATIVE));
468 PetscCall(PCCompositeAddPCType(pc, PCBDDC));
469 PetscCall(PCCompositeAddPCType(pc, PCKSP));
470 PetscCall(PCCompositeGetPC(pc, 1, &pcksp));
471 PetscCall(PCKSPGetKSP(pcksp, &ksprich));
472 PetscCall(KSPSetType(ksprich, KSPRICHARDSON));
473 PetscCall(KSPSetTolerances(ksprich, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT, 1));
474 PetscCall(KSPSetNormType(ksprich, KSP_NORM_NONE));
475 PetscCall(KSPSetConvergenceTest(ksprich, KSPConvergedSkip, NULL, NULL));
476 PetscCall(KSPGetPC(ksprich, &pcjacobi));
477 PetscCall(PCSetType(pcjacobi, PCJACOBI));
478 } else {
479 PetscCall(PCSetType(pc, PCBDDC));
480 }
481 }
482 PetscCall(KSPSetFromOptions(ksp));
483 PetscCall(PetscLogStagePush(stages[0]));
484 PetscCall(KSPSetUp(ksp));
485 PetscCall(PetscLogStagePop());
486
487 PetscCall(MatCreateVecs(A, &x, &b));
488 if (user.random_initial_guess) {
489 /* Solving A x = 0 with random initial guess allows Arnoldi to run for more iterations, thereby yielding a more
490 * complete Hessenberg matrix and more accurate eigenvalues. */
491 PetscCall(VecZeroEntries(b));
492 PetscCall(VecSetRandom(x, NULL));
493 if (user.random_real) PetscCall(VecRealPart(x));
494 if (nullsp) PetscCall(MatNullSpaceRemove(nullsp, x));
495 PetscCall(KSPSetInitialGuessNonzero(ksp, PETSC_TRUE));
496 PetscCall(KSPSetComputeEigenvalues(ksp, PETSC_TRUE));
497 PetscCall(KSPGMRESSetRestart(ksp, 100));
498 } else {
499 PetscCall(VecSetRandom(b, NULL));
500 if (user.random_real) PetscCall(VecRealPart(x));
501 if (nullsp) PetscCall(MatNullSpaceRemove(nullsp, b));
502 }
503 PetscCall(PetscLogStagePush(stages[1]));
504 PetscCall(KSPSolve(ksp, b, x));
505 PetscCall(PetscLogStagePop());
506
507 /* cleanup */
508 PetscCall(VecDestroy(&x));
509 PetscCall(VecDestroy(&b));
510 PetscCall(ISDestroy(&zero));
511 PetscCall(PetscFree(e_glo));
512 PetscCall(MatNullSpaceDestroy(&nullsp));
513 PetscCall(KSPDestroy(&ksp));
514 PetscCall(MatDestroy(&A));
515 PetscCall(DMDestroy(&da));
516 PetscCall(PetscFinalize());
517 return 0;
518 }
519
520 /*TEST
521
522 test:
523 nsize: 8
524 filter: grep -v "variant HERMITIAN"
525 suffix: bddc_1
526 args: -pde_type Poisson -dim 3 -dirichlet 0 -ksp_view -pc_bddc_coarse_redundant_pc_type svd -ksp_error_if_not_converged
527 test:
528 nsize: 8
529 filter: grep -v "variant HERMITIAN"
530 suffix: bddc_2
531 args: -pde_type Poisson -dim 3 -dirichlet 0 -ksp_view -use_global -pc_bddc_coarse_redundant_pc_type svd -ksp_error_if_not_converged
532 test:
533 nsize: 8
534 filter: grep -v "variant HERMITIAN"
535 suffix: bddc_elast
536 args: -pde_type Elasticity -cells 7,9,8 -dim 3 -ksp_view -pc_bddc_coarse_redundant_pc_type svd -ksp_error_if_not_converged -pc_bddc_monolithic
537 test:
538 nsize: 8
539 filter: grep -v "variant HERMITIAN"
540 suffix: bddc_elast_3lev
541 args: -pde_type Elasticity -cells 7,9,8 -dim 3 -ksp_view -pc_bddc_levels 1 -pc_bddc_coarsening_ratio 1 -ksp_error_if_not_converged -pc_bddc_monolithic -pc_bddc_use_faces -pc_bddc_coarse_pc_bddc_corner_selection
542 testset:
543 nsize: 8
544 requires: hpddm slepc defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
545 # on some architectures, this test will converge in 19 or 21 iterations
546 filter: grep -v "variant HERMITIAN" | grep -v " tolerance" | sed -e "s/CONVERGED_RTOL iterations [1-2][91]\{0,1\}$/CONVERGED_RTOL iterations 20/g"
547 args: -pde_type Elasticity -cells 7,9,8 -dim 3 -ksp_view -pc_bddc_levels 1 -pc_bddc_coarsening_ratio 1 -ksp_error_if_not_converged -pc_bddc_monolithic -pc_bddc_use_faces -pc_bddc_coarse_pc_type hpddm -prefix_push pc_bddc_coarse_ -pc_hpddm_levels_1_sub_pc_type cholesky -pc_hpddm_levels_1_eps_nev 6 -pc_hpddm_levels_1_st_pc_factor_shift_type INBLOCKS -prefix_pop -ksp_type fgmres -ksp_max_it 50 -ksp_converged_reason
548 test:
549 args: -pc_bddc_coarse_pc_hpddm_coarse_mat_type baij -options_left no
550 suffix: bddc_elast_3lev_hpddm_baij
551 test:
552 requires: !complex
553 suffix: bddc_elast_3lev_hpddm
554 test:
555 nsize: 8
556 requires: !single
557 filter: grep -v "variant HERMITIAN"
558 suffix: bddc_elast_4lev
559 args: -pde_type Elasticity -cells 7,9,8 -dim 3 -ksp_view -pc_bddc_levels 2 -pc_bddc_coarsening_ratio 2 -ksp_error_if_not_converged -pc_bddc_monolithic -pc_bddc_use_faces -pc_bddc_coarse_pc_bddc_corner_selection -pc_bddc_coarse_l1_pc_bddc_corner_selection -mat_partitioning_type average -options_left 0
560 test:
561 nsize: 8
562 filter: grep -v "variant HERMITIAN"
563 suffix: bddc_elast_deluxe_layers
564 args: -pde_type Elasticity -cells 7,9,8 -dim 3 -ksp_view -pc_bddc_coarse_redundant_pc_type svd -ksp_error_if_not_converged -pc_bddc_monolithic -pc_bddc_use_deluxe_scaling -pc_bddc_schur_layers 1
565 test:
566 nsize: 8
567 filter: grep -v "variant HERMITIAN" | sed -e "s/iterations 1[0-9]/iterations 10/g"
568 suffix: bddc_elast_dir_approx
569 args: -pde_type Elasticity -cells 7,9,8 -dim 3 -ksp_view -pc_bddc_coarse_redundant_pc_type svd -ksp_error_if_not_converged -pc_bddc_monolithic -pc_bddc_dirichlet_pc_type gamg -pc_bddc_dirichlet_pc_gamg_esteig_ksp_max_it 10 -ksp_converged_reason -pc_bddc_dirichlet_approximate
570 test:
571 nsize: 8
572 filter: grep -v "variant HERMITIAN" | sed -e "s/iterations 1[0-9]/iterations 10/g"
573 suffix: bddc_elast_neu_approx
574 args: -pde_type Elasticity -cells 7,9,8 -dim 3 -ksp_view -pc_bddc_coarse_redundant_pc_type svd -ksp_error_if_not_converged -pc_bddc_monolithic -pc_bddc_neumann_pc_type gamg -pc_bddc_neumann_pc_gamg_esteig_ksp_max_it 10 -ksp_converged_reason -pc_bddc_neumann_approximate
575 test:
576 nsize: 8
577 filter: grep -v "variant HERMITIAN" | sed -e "s/iterations 1[0-9]/iterations 10/g"
578 suffix: bddc_elast_both_approx
579 args: -pde_type Elasticity -cells 7,9,8 -dim 3 -ksp_view -pc_bddc_coarse_redundant_pc_type svd -ksp_error_if_not_converged -pc_bddc_monolithic -pc_bddc_dirichlet_pc_type gamg -pc_bddc_dirichlet_pc_gamg_esteig_ksp_max_it 10 -pc_bddc_neumann_pc_type gamg -pc_bddc_neumann_pc_gamg_esteig_ksp_max_it 10 -ksp_converged_reason -pc_bddc_neumann_approximate -pc_bddc_dirichlet_approximate
580 test:
581 nsize: 8
582 filter: grep -v "variant HERMITIAN"
583 suffix: fetidp_1
584 args: -pde_type Poisson -dim 3 -dirichlet 0 -ksp_view -ksp_type fetidp -fetidp_ksp_type cg -fetidp_bddc_pc_bddc_coarse_redundant_pc_type svd -ksp_fetidp_fullyredundant -ksp_error_if_not_converged
585 test:
586 nsize: 8
587 filter: grep -v "variant HERMITIAN"
588 suffix: fetidp_2
589 args: -pde_type Poisson -dim 3 -dirichlet 0 -ksp_view -use_global -ksp_type fetidp -fetidp_ksp_type cg -fetidp_bddc_pc_bddc_coarse_redundant_pc_type svd -ksp_fetidp_fullyredundant -ksp_error_if_not_converged
590 test:
591 nsize: 8
592 filter: grep -v "variant HERMITIAN"
593 suffix: fetidp_elast
594 args: -pde_type Elasticity -cells 9,7,8 -dim 3 -ksp_view -ksp_type fetidp -fetidp_ksp_type cg -fetidp_bddc_pc_bddc_coarse_redundant_pc_type svd -ksp_fetidp_fullyredundant -ksp_error_if_not_converged -fetidp_bddc_pc_bddc_monolithic
595 testset:
596 nsize: 8
597 requires: hpddm slepc defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
598 args: -pde_type Elasticity -cells 12,12 -dim 2 -ksp_converged_reason -pc_type hpddm -pc_hpddm_coarse_correction balanced -pc_hpddm_levels_1_pc_type asm -pc_hpddm_levels_1_pc_asm_overlap 1 -pc_hpddm_levels_1_pc_asm_type basic -pc_hpddm_levels_1_sub_pc_type cholesky -pc_hpddm_levels_1_eps_nev 10 -pc_hpddm_levels_1_st_pc_factor_shift_type INBLOCKS
599 test:
600 args: -mat_is_localmat_type {{aij baij sbaij}shared output} -pc_hpddm_coarse_mat_type {{baij sbaij}shared output}
601 suffix: hpddm
602 output_file: output/ex71_hpddm.out
603 filter: sed -e "s/CONVERGED_RTOL iterations 15/CONVERGED_RTOL iterations 14/g"
604 test:
605 args: -mat_is_localmat_type sbaij -pc_hpddm_coarse_mat_type sbaij -pc_hpddm_levels_1_st_share_sub_ksp -pc_hpddm_levels_1_eps_type lapack -pc_hpddm_levels_1_eps_smallest_magnitude -pc_hpddm_levels_1_st_type shift
606 suffix: hpddm_lapack
607 output_file: output/ex71_hpddm.out
608 testset:
609 nsize: 9
610 args: -assembled_view -pc_bddc_coarse_redundant_pc_type svd -ksp_error_if_not_converged
611 test:
612 args: -dim 1 -cells 12 -pde_type Poisson
613 suffix: dmda_matis_poiss_1d_loc
614 output_file: output/ex71_dmda_matis_poiss_1d.out
615 test:
616 args: -dim 1 -cells 12 -pde_type Poisson -use_global
617 suffix: dmda_matis_poiss_1d_glob
618 output_file: output/ex71_dmda_matis_poiss_1d.out
619 test:
620 args: -dim 1 -cells 12 -pde_type Elasticity
621 suffix: dmda_matis_elast_1d_loc
622 output_file: output/ex71_dmda_matis_elast_1d.out
623 test:
624 args: -dim 1 -cells 12 -pde_type Elasticity -use_global
625 suffix: dmda_matis_elast_1d_glob
626 output_file: output/ex71_dmda_matis_elast_1d.out
627 test:
628 args: -dim 2 -cells 5,7 -pde_type Poisson
629 suffix: dmda_matis_poiss_2d_loc
630 output_file: output/ex71_dmda_matis_poiss_2d.out
631 test:
632 args: -dim 2 -cells 5,7 -pde_type Poisson -use_global
633 suffix: dmda_matis_poiss_2d_glob
634 output_file: output/ex71_dmda_matis_poiss_2d.out
635 test:
636 args: -dim 2 -cells 5,7 -pde_type Elasticity
637 suffix: dmda_matis_elast_2d_loc
638 output_file: output/ex71_dmda_matis_elast_2d.out
639 test:
640 args: -dim 2 -cells 5,7 -pde_type Elasticity -use_global
641 suffix: dmda_matis_elast_2d_glob
642 output_file: output/ex71_dmda_matis_elast_2d.out
643 test:
644 args: -dim 3 -cells 3,3,3 -pde_type Poisson
645 suffix: dmda_matis_poiss_3d_loc
646 output_file: output/ex71_dmda_matis_poiss_3d.out
647 test:
648 args: -dim 3 -cells 3,3,3 -pde_type Poisson -use_global
649 suffix: dmda_matis_poiss_3d_glob
650 output_file: output/ex71_dmda_matis_poiss_3d.out
651 test:
652 args: -dim 3 -cells 3,3,3 -pde_type Elasticity
653 suffix: dmda_matis_elast_3d_loc
654 output_file: output/ex71_dmda_matis_elast_3d.out
655 test:
656 args: -dim 3 -cells 3,3,3 -pde_type Elasticity -use_global
657 suffix: dmda_matis_elast_3d_glob
658 output_file: output/ex71_dmda_matis_elast_3d.out
659 test:
660 nsize: 8
661 filter: grep -v "variant HERMITIAN" | sed -e "s/CONVERGED_RTOL iterations 1[0-9]/CONVERGED_RTOL iterations 13/g"
662 suffix: bddc_elast_deluxe_layers_adapt
663 requires: mumps !complex
664 args: -pde_type Elasticity -cells 7,9,8 -dim 3 -ksp_converged_reason -pc_bddc_coarse_redundant_pc_type svd -ksp_error_if_not_converged -pc_bddc_monolithic -sub_schurs_mat_solver_type mumps -pc_bddc_use_deluxe_scaling -pc_bddc_adaptive_threshold 2.0 -pc_bddc_schur_layers {{1 10}separate_output} -pc_bddc_adaptive_userdefined {{0 1}separate output} -sub_schurs_schur_mat_type seqdense
665 # gitlab runners have a quite old MKL (2016) which interacts badly with AMD machines (not Intel-based ones!)
666 # this is the reason behind the filtering rule
667 test:
668 nsize: 8
669 suffix: bddc_elast_deluxe_layers_adapt_mkl_pardiso
670 filter: sed -e "s/CONVERGED_RTOL iterations [1-2][0-9]/CONVERGED_RTOL iterations 13/g" | sed -e "s/CONVERGED_RTOL iterations 6/CONVERGED_RTOL iterations 5/g"
671 requires: mkl_pardiso !complex
672 args: -pde_type Elasticity -cells 7,9,8 -dim 3 -ksp_converged_reason -pc_bddc_coarse_redundant_pc_type svd -ksp_error_if_not_converged -pc_bddc_monolithic -sub_schurs_mat_solver_type mkl_pardiso -sub_schurs_mat_mkl_pardiso_65 1 -pc_bddc_use_deluxe_scaling -pc_bddc_adaptive_threshold 2.0 -pc_bddc_schur_layers {{1 10}separate_output} -pc_bddc_adaptive_userdefined {{0 1}separate output} -sub_schurs_schur_mat_type seqdense
673 test:
674 nsize: 8
675 filter: grep -v "variant HERMITIAN"
676 suffix: bddc_cusparse
677 # no kokkos since it seems kokkos's resource demand is too much with 8 ranks and the test will fail on cuda related initialization.
678 requires: cuda !kokkos
679 args: -pde_type Poisson -cells 7,9,8 -dim 3 -ksp_view -pc_bddc_coarse_redundant_pc_type svd -ksp_error_if_not_converged -pc_bddc_dirichlet_pc_type cholesky -pc_bddc_dirichlet_pc_factor_mat_solver_type cusparse -pc_bddc_dirichlet_pc_factor_mat_ordering_type nd -pc_bddc_neumann_pc_type cholesky -pc_bddc_neumann_pc_factor_mat_solver_type cusparse -pc_bddc_neumann_pc_factor_mat_ordering_type nd -mat_is_localmat_type aijcusparse
680 test:
681 nsize: 8
682 filter: grep -v "variant HERMITIAN"
683 suffix: bddc_elast_deluxe_layers_adapt_cuda
684 requires: !complex mumps cuda defined(PETSC_HAVE_CUSOLVERDNDPOTRI)
685 args: -pde_type Elasticity -cells 7,9,8 -dim 3 -ksp_converged_reason -pc_bddc_coarse_redundant_pc_type svd -ksp_error_if_not_converged -pc_bddc_monolithic -sub_schurs_mat_solver_type mumps -pc_bddc_use_deluxe_scaling -pc_bddc_adaptive_threshold 2.0 -pc_bddc_schur_layers {{1 10}separate_output} -pc_bddc_adaptive_userdefined {{0 1}separate output} -mat_is_localmat_type seqaijcusparse -sub_schurs_schur_mat_type {{seqdensecuda seqdense}}
686 test:
687 nsize: 8
688 filter: grep -v "variant HERMITIAN" | grep -v "I-node routines" | sed -e "s/seqaijcusparse/seqaij/g"
689 suffix: bddc_elast_deluxe_layers_adapt_cuda_approx
690 requires: !complex mumps cuda defined(PETSC_HAVE_CUSOLVERDNDPOTRI)
691 args: -pde_type Elasticity -cells 7,9,8 -dim 3 -ksp_view -pc_bddc_coarse_redundant_pc_type svd -ksp_error_if_not_converged -pc_bddc_monolithic -sub_schurs_mat_solver_type mumps -pc_bddc_use_deluxe_scaling -pc_bddc_adaptive_threshold 2.0 -pc_bddc_schur_layers 1 -mat_is_localmat_type {{seqaij seqaijcusparse}separate output} -sub_schurs_schur_mat_type {{seqdensecuda seqdense}} -pc_bddc_dirichlet_pc_type gamg -pc_bddc_dirichlet_approximate -pc_bddc_neumann_pc_type gamg -pc_bddc_neumann_approximate -pc_bddc_dirichlet_pc_gamg_esteig_ksp_max_it 10 -pc_bddc_neumann_pc_gamg_esteig_ksp_max_it 10
692 test:
693 nsize: 8
694 suffix: bddc_elast_deluxe_layers_adapt_mkl_pardiso_cuda
695 requires: !complex mkl_pardiso cuda defined(PETSC_HAVE_CUSOLVERDNDPOTRI)
696 filter: sed -e "s/CONVERGED_RTOL iterations 6/CONVERGED_RTOL iterations 5/g"
697 args: -pde_type Elasticity -cells 7,9,8 -dim 3 -ksp_converged_reason -pc_bddc_coarse_redundant_pc_type svd -ksp_error_if_not_converged -pc_bddc_monolithic -sub_schurs_mat_solver_type mkl_pardiso -sub_schurs_mat_mkl_pardiso_65 1 -pc_bddc_use_deluxe_scaling -pc_bddc_adaptive_threshold 2.0 -pc_bddc_schur_layers {{1 10}separate_output} -pc_bddc_adaptive_userdefined {{0 1}separate output} -mat_is_localmat_type seqaijcusparse -sub_schurs_schur_mat_type {{seqdensecuda seqdense}}
698
699 testset:
700 nsize: 2
701 output_file: output/ex71_aij_dmda_preall.out
702 filter: sed -e "s/CONVERGED_RTOL iterations 7/CONVERGED_RTOL iterations 6/g"
703 args: -pde_type Poisson -dim 1 -cells 6 -pc_type none -ksp_converged_reason
704 test:
705 suffix: aijviennacl_dmda_preall
706 requires: viennacl
707 args: -dm_mat_type aijviennacl -dm_preallocate_only {{0 1}} -dirichlet {{0 1}}
708 # -dm_preallocate_only 0 is broken
709 test:
710 suffix: aijcusparse_dmda_preall
711 requires: cuda
712 args: -dm_mat_type aijcusparse -dm_preallocate_only -dirichlet {{0 1}}
713 test:
714 suffix: aij_dmda_preall
715 args: -dm_mat_type aij -dm_preallocate_only {{0 1}} -dirichlet {{0 1}}
716 testset:
717 nsize: 4
718 args: -dim 2 -cells 16,16 -periodicity 1,1 -random_initial_guess -random_real -sub_0_pc_bddc_switch_static -use_composite_pc -ksp_monitor -ksp_converged_reason -ksp_type gmres -ksp_view_singularvalues -ksp_view_eigenvalues -sub_0_pc_bddc_use_edges 0 -sub_0_pc_bddc_coarse_pc_type svd -sub_1_ksp_ksp_max_it 1 -sub_1_ksp_ksp_richardson_scale 2.3
719 test:
720 args: -sub_0_pc_bddc_interface_ext_type lump
721 suffix: composite_bddc_lumped
722 test:
723 requires: !single
724 args: -sub_0_pc_bddc_interface_ext_type dirichlet
725 suffix: composite_bddc_dirichlet
726
727 # GDSW tests
728 testset:
729 nsize: 8
730 filter: grep -v "variant HERMITIAN"
731 args: -cells 7,9,8 -dim 3 -ksp_view -ksp_error_if_not_converged -pc_type mg -pc_mg_levels 2 -pc_mg_galerkin -pc_mg_adapt_interp_coarse_space gdsw -mg_levels_pc_type asm -mg_levels_sub_pc_type icc -mg_coarse_redundant_pc_type cholesky
732 test:
733 suffix: gdsw_poisson
734 args: -pde_type Poisson
735 test:
736 requires: mumps !complex
737 suffix: gdsw_poisson_adaptive
738 args: -pde_type Poisson -mg_levels_gdsw_tolerance 0.01 -ksp_monitor_singular_value -mg_levels_gdsw_userdefined {{0 1}separate output} -mg_levels_gdsw_pseudo_pc_type qr
739 test:
740 suffix: gdsw_elast
741 args: -pde_type Elasticity
742 test:
743 requires: hpddm
744 suffix: gdsw_elast_hpddm
745 args: -pde_type Elasticity -mg_levels_gdsw_ksp_type hpddm -mg_levels_gdsw_ksp_hpddm_type cg
746 test:
747 requires: mumps !complex
748 suffix: gdsw_elast_adaptive
749 args: -pde_type Elasticity -mg_levels_gdsw_tolerance 0.01 -ksp_monitor_singular_value -mg_levels_gdsw_userdefined {{0 1}separate output}
750
751 # Multi-Element tests
752 test:
753 nsize: {{1 2 3}}
754 suffix: bddc_multi_element
755 args: -cells 3,3,3 -dim 3 -ksp_error_if_not_converged -multi_element -pde_type {{Poisson Elasticity}} -ksp_converged_reason
756
757 test:
758 suffix: bddc_multi_square
759 output_file: output/ex71_bddc_multi_element.out
760 args: -cells 2,2 -dim 2 -ksp_error_if_not_converged -multi_element -pc_bddc_local_mat_graph_square 4 -ksp_converged_reason
761
762 TEST*/
763