Actual source code: ex71.c
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";
5: /* Contributed by Wim Vanroose <wim@vanroo.se> */
7: #include <petscksp.h>
8: #include <petscpc.h>
9: #include <petscdm.h>
10: #include <petscdmda.h>
11: #include <petscdmplex.h>
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};
99: typedef enum {
100: PDE_POISSON,
101: PDE_ELASTICITY
102: } PDEType;
104: typedef struct {
105: PDEType pde;
106: PetscInt dim;
107: PetscInt dof;
108: PetscInt cells[3];
109: PetscBool useglobal;
110: PetscBool dirbc;
111: PetscBool per[3];
112: PetscBool test;
113: PetscScalar *elemMat;
114: PetscBool use_composite_pc;
115: PetscBool random_initial_guess;
116: PetscBool random_real;
117: } AppCtx;
119: static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
120: {
121: const char *pdeTypes[2] = {"Poisson", "Elasticity"};
122: PetscInt n, pde;
124: PetscFunctionBeginUser;
125: options->pde = PDE_POISSON;
126: options->elemMat = NULL;
127: options->dim = 1;
128: options->cells[0] = 8;
129: options->cells[1] = 6;
130: options->cells[2] = 4;
131: options->useglobal = PETSC_FALSE;
132: options->dirbc = PETSC_TRUE;
133: options->test = PETSC_FALSE;
134: options->per[0] = PETSC_FALSE;
135: options->per[1] = PETSC_FALSE;
136: options->per[2] = PETSC_FALSE;
137: options->use_composite_pc = PETSC_FALSE;
138: options->random_initial_guess = PETSC_FALSE;
139: options->random_real = PETSC_FALSE;
141: PetscOptionsBegin(comm, NULL, "Problem Options", NULL);
142: pde = options->pde;
143: PetscCall(PetscOptionsEList("-pde_type", "The PDE type", __FILE__, pdeTypes, 2, pdeTypes[options->pde], &pde, NULL));
144: options->pde = (PDEType)pde;
145: PetscCall(PetscOptionsInt("-dim", "The topological mesh dimension", __FILE__, options->dim, &options->dim, NULL));
146: PetscCall(PetscOptionsIntArray("-cells", "The mesh division", __FILE__, options->cells, (n = 3, &n), NULL));
147: PetscCall(PetscOptionsBoolArray("-periodicity", "The mesh periodicity", __FILE__, options->per, (n = 3, &n), NULL));
148: PetscCall(PetscOptionsBool("-use_global", "Test MatSetValues", __FILE__, options->useglobal, &options->useglobal, NULL));
149: PetscCall(PetscOptionsBool("-dirichlet", "Use dirichlet BC", __FILE__, options->dirbc, &options->dirbc, NULL));
150: PetscCall(PetscOptionsBool("-test_assembly", "Test MATIS assembly", __FILE__, options->test, &options->test, NULL));
151: PetscCall(PetscOptionsBool("-use_composite_pc", "Multiplicative composite with BDDC + Richardson/Jacobi", __FILE__, options->use_composite_pc, &options->use_composite_pc, NULL));
152: 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));
153: 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));
154: PetscOptionsEnd();
156: for (n = options->dim; n < 3; n++) options->cells[n] = 0;
157: if (options->per[0]) options->dirbc = PETSC_FALSE;
159: /* element matrices */
160: switch (options->pde) {
161: case PDE_ELASTICITY:
162: options->dof = options->dim;
163: switch (options->dim) {
164: case 1:
165: options->elemMat = elast_1D_emat;
166: break;
167: case 2:
168: options->elemMat = elast_2D_emat;
169: break;
170: case 3:
171: options->elemMat = elast_3D_emat;
172: break;
173: default:
174: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Unsupported dimension %" PetscInt_FMT, options->dim);
175: }
176: break;
177: case PDE_POISSON:
178: options->dof = 1;
179: switch (options->dim) {
180: case 1:
181: options->elemMat = poiss_1D_emat;
182: break;
183: case 2:
184: options->elemMat = poiss_2D_emat;
185: break;
186: case 3:
187: options->elemMat = poiss_3D_emat;
188: break;
189: default:
190: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Unsupported dimension %" PetscInt_FMT, options->dim);
191: }
192: break;
193: default:
194: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Unsupported PDE %d", options->pde);
195: }
196: PetscFunctionReturn(PETSC_SUCCESS);
197: }
199: int main(int argc, char **args)
200: {
201: AppCtx user;
202: KSP ksp;
203: PC pc;
204: Mat A;
205: DM da;
206: Vec x, b, xcoor, xcoorl;
207: IS zero;
208: ISLocalToGlobalMapping map;
209: MatNullSpace nullsp = NULL;
210: PetscInt i;
211: PetscInt nel, nen; /* Number of elements & element nodes */
212: const PetscInt *e_loc; /* Local indices of element nodes (in local element order) */
213: PetscInt *e_glo = NULL; /* Global indices of element nodes (in local element order) */
214: PetscInt nodes[3];
215: PetscBool ismatis;
216: #if defined(PETSC_USE_LOG)
217: PetscLogStage stages[2];
218: #endif
220: PetscFunctionBeginUser;
221: PetscCall(PetscInitialize(&argc, &args, (char *)0, 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: }
239: PetscCall(PetscLogStageRegister("KSPSetUp", &stages[0]));
240: PetscCall(PetscLogStageRegister("KSPSolve", &stages[1]));
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));
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: }
273: /* we reorder the indices since the element matrices are given in lexicographic order,
274: whereas the elements indices returned by DMDAGetElements follow the usual FEM ordering
275: i.e., element matrices DMDA ordering
276: 2---3 3---2
277: / / / /
278: 0---1 0---1
279: */
280: for (i = 0; i < nel; ++i) {
281: PetscInt ord[8] = {0, 1, 3, 2, 4, 5, 7, 6};
282: PetscInt j, idxs[8];
284: PetscCheck(nen <= 8, PETSC_COMM_WORLD, PETSC_ERR_SUP, "Not coded");
285: if (!e_glo) {
286: for (j = 0; j < nen; j++) idxs[j] = e_loc[i * nen + ord[j]];
287: PetscCall(MatSetValuesBlockedLocal(A, nen, idxs, nen, idxs, user.elemMat, ADD_VALUES));
288: } else {
289: for (j = 0; j < nen; j++) idxs[j] = e_glo[i * nen + ord[j]];
290: PetscCall(MatSetValuesBlocked(A, nen, idxs, nen, idxs, user.elemMat, ADD_VALUES));
291: }
292: }
293: PetscCall(DMDARestoreElements(da, &nel, &nen, &e_loc));
294: PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
295: PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
296: PetscCall(MatSetOption(A, MAT_SPD, PETSC_TRUE));
297: PetscCall(MatSetOption(A, MAT_SPD_ETERNAL, PETSC_TRUE));
299: /* Boundary conditions */
300: zero = NULL;
301: if (user.dirbc) { /* fix one side of DMDA */
302: Vec nat, glob;
303: PetscScalar *vals;
304: PetscInt n, *idx, j, st;
306: n = PetscGlobalRank ? 0 : (user.cells[1] + 1) * (user.cells[2] + 1);
307: PetscCall(ISCreateStride(PETSC_COMM_WORLD, n, 0, user.cells[0] + 1, &zero));
308: if (user.dof > 1) { /* zero all components */
309: const PetscInt *idx;
310: IS bzero;
312: PetscCall(ISGetIndices(zero, (const PetscInt **)&idx));
313: PetscCall(ISCreateBlock(PETSC_COMM_WORLD, user.dof, n, idx, PETSC_COPY_VALUES, &bzero));
314: PetscCall(ISRestoreIndices(zero, (const PetscInt **)&idx));
315: PetscCall(ISDestroy(&zero));
316: zero = bzero;
317: }
318: /* map indices from natural to global */
319: PetscCall(DMDACreateNaturalVector(da, &nat));
320: PetscCall(ISGetLocalSize(zero, &n));
321: PetscCall(PetscMalloc1(n, &vals));
322: for (i = 0; i < n; i++) vals[i] = 1.0;
323: PetscCall(ISGetIndices(zero, (const PetscInt **)&idx));
324: PetscCall(VecSetValues(nat, n, idx, vals, INSERT_VALUES));
325: PetscCall(ISRestoreIndices(zero, (const PetscInt **)&idx));
326: PetscCall(PetscFree(vals));
327: PetscCall(VecAssemblyBegin(nat));
328: PetscCall(VecAssemblyEnd(nat));
329: PetscCall(DMCreateGlobalVector(da, &glob));
330: PetscCall(DMDANaturalToGlobalBegin(da, nat, INSERT_VALUES, glob));
331: PetscCall(DMDANaturalToGlobalEnd(da, nat, INSERT_VALUES, glob));
332: PetscCall(VecDestroy(&nat));
333: PetscCall(ISDestroy(&zero));
334: PetscCall(VecGetLocalSize(glob, &n));
335: PetscCall(PetscMalloc1(n, &idx));
336: PetscCall(VecGetOwnershipRange(glob, &st, NULL));
337: PetscCall(VecGetArray(glob, &vals));
338: for (i = 0, j = 0; i < n; i++)
339: if (PetscRealPart(vals[i]) == 1.0) idx[j++] = i + st;
340: PetscCall(VecRestoreArray(glob, &vals));
341: PetscCall(VecDestroy(&glob));
342: PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, j, idx, PETSC_OWN_POINTER, &zero));
343: PetscCall(MatZeroRowsColumnsIS(A, zero, 1.0, NULL, NULL));
344: } else {
345: switch (user.pde) {
346: case PDE_POISSON:
347: PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, NULL, &nullsp));
348: break;
349: case PDE_ELASTICITY:
350: PetscCall(MatNullSpaceCreateRigidBody(xcoor, &nullsp));
351: break;
352: }
353: /* with periodic BC and Elasticity, just the displacements are in the nullspace
354: this is no harm since we eliminate all the components of the rhs */
355: PetscCall(MatSetNullSpace(A, nullsp));
356: }
358: if (user.test) {
359: Mat AA;
361: PetscCall(MatConvert(A, MATAIJ, MAT_INITIAL_MATRIX, &AA));
362: PetscCall(MatViewFromOptions(AA, NULL, "-assembled_view"));
363: PetscCall(MatDestroy(&AA));
364: }
366: /* Attach near null space for elasticity */
367: if (user.pde == PDE_ELASTICITY) {
368: MatNullSpace nearnullsp;
370: PetscCall(MatNullSpaceCreateRigidBody(xcoor, &nearnullsp));
371: PetscCall(MatSetNearNullSpace(A, nearnullsp));
372: PetscCall(MatNullSpaceDestroy(&nearnullsp));
373: }
375: /* we may want to use MG for the local solvers: attach local nearnullspace to the local matrices */
376: PetscCall(DMGetCoordinatesLocal(da, &xcoorl));
377: PetscCall(PetscObjectTypeCompare((PetscObject)A, MATIS, &ismatis));
378: if (ismatis) {
379: MatNullSpace lnullsp = NULL;
380: Mat lA;
382: PetscCall(MatISGetLocalMat(A, &lA));
383: if (user.pde == PDE_ELASTICITY) {
384: Vec lc;
385: ISLocalToGlobalMapping l2l;
386: IS is;
387: const PetscScalar *a;
388: const PetscInt *idxs;
389: PetscInt n, bs;
391: /* when using a DMDA, the local matrices have an additional local-to-local map
392: that maps from the DA local ordering to the ordering induced by the elements */
393: PetscCall(MatCreateVecs(lA, &lc, NULL));
394: PetscCall(MatGetLocalToGlobalMapping(lA, &l2l, NULL));
395: PetscCall(VecSetLocalToGlobalMapping(lc, l2l));
396: PetscCall(VecSetOption(lc, VEC_IGNORE_NEGATIVE_INDICES, PETSC_TRUE));
397: PetscCall(VecGetLocalSize(xcoorl, &n));
398: PetscCall(VecGetBlockSize(xcoorl, &bs));
399: PetscCall(ISCreateStride(PETSC_COMM_SELF, n / bs, 0, 1, &is));
400: PetscCall(ISGetIndices(is, &idxs));
401: PetscCall(VecGetArrayRead(xcoorl, &a));
402: PetscCall(VecSetValuesBlockedLocal(lc, n / bs, idxs, a, INSERT_VALUES));
403: PetscCall(VecAssemblyBegin(lc));
404: PetscCall(VecAssemblyEnd(lc));
405: PetscCall(VecRestoreArrayRead(xcoorl, &a));
406: PetscCall(ISRestoreIndices(is, &idxs));
407: PetscCall(ISDestroy(&is));
408: PetscCall(MatNullSpaceCreateRigidBody(lc, &lnullsp));
409: PetscCall(VecDestroy(&lc));
410: } else if (user.pde == PDE_POISSON) {
411: PetscCall(MatNullSpaceCreate(PETSC_COMM_SELF, PETSC_TRUE, 0, NULL, &lnullsp));
412: }
413: PetscCall(MatSetNearNullSpace(lA, lnullsp));
414: PetscCall(MatNullSpaceDestroy(&lnullsp));
415: PetscCall(MatISRestoreLocalMat(A, &lA));
416: }
418: PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
419: PetscCall(KSPSetOperators(ksp, A, A));
420: PetscCall(KSPSetType(ksp, KSPCG));
421: PetscCall(KSPGetPC(ksp, &pc));
422: if (user.use_composite_pc) {
423: PC pcksp, pcjacobi;
424: KSP ksprich;
425: PetscCall(PCSetType(pc, PCCOMPOSITE));
426: PetscCall(PCCompositeSetType(pc, PC_COMPOSITE_MULTIPLICATIVE));
427: PetscCall(PCCompositeAddPCType(pc, PCBDDC));
428: PetscCall(PCCompositeAddPCType(pc, PCKSP));
429: PetscCall(PCCompositeGetPC(pc, 1, &pcksp));
430: PetscCall(PCKSPGetKSP(pcksp, &ksprich));
431: PetscCall(KSPSetType(ksprich, KSPRICHARDSON));
432: PetscCall(KSPSetTolerances(ksprich, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 1));
433: PetscCall(KSPSetNormType(ksprich, KSP_NORM_NONE));
434: PetscCall(KSPSetConvergenceTest(ksprich, KSPConvergedSkip, NULL, NULL));
435: PetscCall(KSPGetPC(ksprich, &pcjacobi));
436: PetscCall(PCSetType(pcjacobi, PCJACOBI));
437: } else {
438: PetscCall(PCSetType(pc, PCBDDC));
439: }
440: /* PetscCall(PCBDDCSetDirichletBoundaries(pc,zero)); */
441: PetscCall(KSPSetFromOptions(ksp));
442: PetscCall(PetscLogStagePush(stages[0]));
443: PetscCall(KSPSetUp(ksp));
444: PetscCall(PetscLogStagePop());
446: PetscCall(MatCreateVecs(A, &x, &b));
447: if (user.random_initial_guess) {
448: /* Solving A x = 0 with random initial guess allows Arnoldi to run for more iterations, thereby yielding a more
449: * complete Hessenberg matrix and more accurate eigenvalues. */
450: PetscCall(VecZeroEntries(b));
451: PetscCall(VecSetRandom(x, NULL));
452: if (user.random_real) PetscCall(VecRealPart(x));
453: if (nullsp) PetscCall(MatNullSpaceRemove(nullsp, x));
454: PetscCall(KSPSetInitialGuessNonzero(ksp, PETSC_TRUE));
455: PetscCall(KSPSetComputeEigenvalues(ksp, PETSC_TRUE));
456: PetscCall(KSPGMRESSetRestart(ksp, 100));
457: } else {
458: PetscCall(VecSetRandom(b, NULL));
459: if (user.random_real) PetscCall(VecRealPart(x));
460: if (nullsp) PetscCall(MatNullSpaceRemove(nullsp, b));
461: }
462: PetscCall(PetscLogStagePush(stages[1]));
463: PetscCall(KSPSolve(ksp, b, x));
464: PetscCall(PetscLogStagePop());
466: /* cleanup */
467: PetscCall(VecDestroy(&x));
468: PetscCall(VecDestroy(&b));
469: PetscCall(ISDestroy(&zero));
470: PetscCall(PetscFree(e_glo));
471: PetscCall(MatNullSpaceDestroy(&nullsp));
472: PetscCall(KSPDestroy(&ksp));
473: PetscCall(MatDestroy(&A));
474: PetscCall(DMDestroy(&da));
475: PetscCall(PetscFinalize());
476: return 0;
477: }
479: /*TEST
481: test:
482: nsize: 8
483: filter: grep -v "variant HERMITIAN"
484: suffix: bddc_1
485: args: -pde_type Poisson -dim 3 -dirichlet 0 -ksp_view -pc_bddc_coarse_redundant_pc_type svd -ksp_error_if_not_converged
486: test:
487: nsize: 8
488: filter: grep -v "variant HERMITIAN"
489: suffix: bddc_2
490: args: -pde_type Poisson -dim 3 -dirichlet 0 -ksp_view -use_global -pc_bddc_coarse_redundant_pc_type svd -ksp_error_if_not_converged
491: test:
492: nsize: 8
493: filter: grep -v "variant HERMITIAN"
494: suffix: bddc_elast
495: 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
496: test:
497: nsize: 8
498: filter: grep -v "variant HERMITIAN"
499: suffix: bddc_elast_3lev
500: 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
501: testset:
502: nsize: 8
503: requires: hpddm slepc defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
504: # on some architectures, this test will converge in 19 or 21 iterations
505: filter: grep -v "variant HERMITIAN" | grep -v " tolerance" | sed -e "s/CONVERGED_RTOL iterations [1-2][91]\{0,1\}$/CONVERGED_RTOL iterations 20/g"
506: 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
507: test:
508: args: -pc_bddc_coarse_pc_hpddm_coarse_mat_type baij -options_left no
509: suffix: bddc_elast_3lev_hpddm_baij
510: test:
511: requires: !complex
512: suffix: bddc_elast_3lev_hpddm
513: test:
514: nsize: 8
515: requires: !single
516: filter: grep -v "variant HERMITIAN"
517: suffix: bddc_elast_4lev
518: 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
519: test:
520: nsize: 8
521: filter: grep -v "variant HERMITIAN"
522: suffix: bddc_elast_deluxe_layers
523: 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
524: test:
525: nsize: 8
526: filter: grep -v "variant HERMITIAN" | sed -e "s/iterations 1[0-9]/iterations 10/g"
527: suffix: bddc_elast_dir_approx
528: 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
529: test:
530: nsize: 8
531: filter: grep -v "variant HERMITIAN" | sed -e "s/iterations 1[0-9]/iterations 10/g"
532: suffix: bddc_elast_neu_approx
533: 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
534: test:
535: nsize: 8
536: filter: grep -v "variant HERMITIAN" | sed -e "s/iterations 1[0-9]/iterations 10/g"
537: suffix: bddc_elast_both_approx
538: 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
539: test:
540: nsize: 8
541: filter: grep -v "variant HERMITIAN"
542: suffix: fetidp_1
543: 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
544: test:
545: nsize: 8
546: filter: grep -v "variant HERMITIAN"
547: suffix: fetidp_2
548: 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
549: test:
550: nsize: 8
551: filter: grep -v "variant HERMITIAN"
552: suffix: fetidp_elast
553: 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
554: testset:
555: nsize: 8
556: requires: hpddm slepc defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
557: 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
558: test:
559: args: -matis_localmat_type {{aij baij sbaij}shared output} -pc_hpddm_coarse_mat_type {{baij sbaij}shared output}
560: suffix: hpddm
561: output_file: output/ex71_hpddm.out
562: test:
563: args: -matis_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
564: suffix: hpddm_lapack
565: output_file: output/ex71_hpddm.out
566: testset:
567: nsize: 9
568: args: -test_assembly -assembled_view -pc_bddc_coarse_redundant_pc_type svd -ksp_error_if_not_converged
569: test:
570: args: -dim 1 -cells 12 -pde_type Poisson
571: suffix: dmda_matis_poiss_1d_loc
572: output_file: output/ex71_dmda_matis_poiss_1d.out
573: test:
574: args: -dim 1 -cells 12 -pde_type Poisson -use_global
575: suffix: dmda_matis_poiss_1d_glob
576: output_file: output/ex71_dmda_matis_poiss_1d.out
577: test:
578: args: -dim 1 -cells 12 -pde_type Elasticity
579: suffix: dmda_matis_elast_1d_loc
580: output_file: output/ex71_dmda_matis_elast_1d.out
581: test:
582: args: -dim 1 -cells 12 -pde_type Elasticity -use_global
583: suffix: dmda_matis_elast_1d_glob
584: output_file: output/ex71_dmda_matis_elast_1d.out
585: test:
586: args: -dim 2 -cells 5,7 -pde_type Poisson
587: suffix: dmda_matis_poiss_2d_loc
588: output_file: output/ex71_dmda_matis_poiss_2d.out
589: test:
590: args: -dim 2 -cells 5,7 -pde_type Poisson -use_global
591: suffix: dmda_matis_poiss_2d_glob
592: output_file: output/ex71_dmda_matis_poiss_2d.out
593: test:
594: args: -dim 2 -cells 5,7 -pde_type Elasticity
595: suffix: dmda_matis_elast_2d_loc
596: output_file: output/ex71_dmda_matis_elast_2d.out
597: test:
598: args: -dim 2 -cells 5,7 -pde_type Elasticity -use_global
599: suffix: dmda_matis_elast_2d_glob
600: output_file: output/ex71_dmda_matis_elast_2d.out
601: test:
602: args: -dim 3 -cells 3,3,3 -pde_type Poisson
603: suffix: dmda_matis_poiss_3d_loc
604: output_file: output/ex71_dmda_matis_poiss_3d.out
605: test:
606: args: -dim 3 -cells 3,3,3 -pde_type Poisson -use_global
607: suffix: dmda_matis_poiss_3d_glob
608: output_file: output/ex71_dmda_matis_poiss_3d.out
609: test:
610: args: -dim 3 -cells 3,3,3 -pde_type Elasticity
611: suffix: dmda_matis_elast_3d_loc
612: output_file: output/ex71_dmda_matis_elast_3d.out
613: test:
614: args: -dim 3 -cells 3,3,3 -pde_type Elasticity -use_global
615: suffix: dmda_matis_elast_3d_glob
616: output_file: output/ex71_dmda_matis_elast_3d.out
617: test:
618: nsize: 8
619: filter: grep -v "variant HERMITIAN"
620: suffix: bddc_elast_deluxe_layers_adapt
621: requires: mumps !complex
622: 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
623: # gitlab runners have a quite old MKL (2016) which interacts badly with AMD machines (not Intel-based ones!)
624: # this is the reason behind the filtering rule
625: test:
626: nsize: 8
627: suffix: bddc_elast_deluxe_layers_adapt_mkl_pardiso
628: 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"
629: requires: mkl_pardiso !complex
630: 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
631: test:
632: nsize: 8
633: filter: grep -v "variant HERMITIAN"
634: suffix: bddc_cusparse
635: # no kokkos since it seems kokkos's resource demand is too much with 8 ranks and the test will fail on cuda related initialization.
636: requires: cuda !kokkos
637: 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
638: test:
639: nsize: 8
640: filter: grep -v "variant HERMITIAN"
641: suffix: bddc_elast_deluxe_layers_adapt_cuda
642: requires: !complex mumps cuda defined(PETSC_HAVE_CUSOLVERDNDPOTRI)
643: 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 seqaijcusparse -sub_schurs_schur_mat_type {{seqdensecuda seqdense}}
644: test:
645: nsize: 8
646: filter: grep -v "variant HERMITIAN" | grep -v "I-node routines" | sed -e "s/seqaijcusparse/seqaij/g"
647: suffix: bddc_elast_deluxe_layers_adapt_cuda_approx
648: requires: !complex mumps cuda defined(PETSC_HAVE_CUSOLVERDNDPOTRI)
649: 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 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
650: test:
651: nsize: 8
652: suffix: bddc_elast_deluxe_layers_adapt_mkl_pardiso_cuda
653: requires: !complex mkl_pardiso cuda defined(PETSC_HAVE_CUSOLVERDNDPOTRI)
654: filter: sed -e "s/CONVERGED_RTOL iterations 6/CONVERGED_RTOL iterations 5/g"
655: 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 seqaijcusparse -sub_schurs_schur_mat_type {{seqdensecuda seqdense}}
657: testset:
658: nsize: 2
659: output_file: output/ex71_aij_dmda_preall.out
660: filter: sed -e "s/CONVERGED_RTOL iterations 7/CONVERGED_RTOL iterations 6/g"
661: args: -pde_type Poisson -dim 1 -cells 6 -pc_type none -ksp_converged_reason
662: test:
663: suffix: aijviennacl_dmda_preall
664: requires: viennacl
665: args: -dm_mat_type aijviennacl -dm_preallocate_only {{0 1}} -dirichlet {{0 1}}
666: # -dm_preallocate_only 0 is broken
667: test:
668: suffix: aijcusparse_dmda_preall
669: requires: cuda
670: args: -dm_mat_type aijcusparse -dm_preallocate_only -dirichlet {{0 1}}
671: test:
672: suffix: aij_dmda_preall
673: args: -dm_mat_type aij -dm_preallocate_only {{0 1}} -dirichlet {{0 1}}
674: testset:
675: nsize: 4
676: 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
677: test:
678: args: -sub_0_pc_bddc_interface_ext_type lump
679: suffix: composite_bddc_lumped
680: test:
681: requires: !single
682: args: -sub_0_pc_bddc_interface_ext_type dirichlet
683: suffix: composite_bddc_dirichlet
685: testset:
686: nsize: 8
687: filter: grep -v "variant HERMITIAN"
688: 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
689: test:
690: suffix: gdsw_poisson
691: args: -pde_type Poisson
692: test:
693: requires: mumps !complex
694: suffix: gdsw_poisson_adaptive
695: 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
696: test:
697: suffix: gdsw_elast
698: args: -pde_type Elasticity
699: test:
700: requires: hpddm
701: suffix: gdsw_elast_hpddm
702: args: -pde_type Elasticity -mg_levels_gdsw_ksp_type hpddm -mg_levels_gdsw_ksp_hpddm_type cg
703: test:
704: requires: mumps !complex
705: suffix: gdsw_elast_adaptive
706: args: -pde_type Elasticity -mg_levels_gdsw_tolerance 0.01 -ksp_monitor_singular_value -mg_levels_gdsw_userdefined {{0 1}separate output}
708: TEST*/