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