xref: /petsc/src/snes/tutorials/ex56.c (revision e0b7e82fd3cf27fce84cc3e37e8d70a5c36a2d4e)
147d993e7Ssuyashtn /* Portions of this code are under:
247d993e7Ssuyashtn    Copyright (C) 2022 Advanced Micro Devices, Inc. All rights reserved.
347d993e7Ssuyashtn */
43667b0a5SMark Adams static char help[] = "3D tensor hexahedra & 3D Laplacian displacement finite element formulation\n\
5c4762a1bSJed Brown of linear elasticity.  E=1.0, nu=1/3.\n\
6c4762a1bSJed Brown Unit cube domain with Dirichlet boundary\n\n";
7c4762a1bSJed Brown 
8c4762a1bSJed Brown #include <petscdmplex.h>
9c4762a1bSJed Brown #include <petscsnes.h>
10c4762a1bSJed Brown #include <petscds.h>
11c4762a1bSJed Brown #include <petscdmforest.h>
12c4762a1bSJed Brown 
133667b0a5SMark Adams static PetscReal s_soft_alpha = 0.01;
14c4762a1bSJed Brown static PetscReal s_mu         = 0.4;
15c4762a1bSJed Brown static PetscReal s_lambda     = 0.4;
16c4762a1bSJed Brown 
17d71ae5a4SJacob Faibussowitsch static void f0_bd_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
18d71ae5a4SJacob Faibussowitsch {
19c4762a1bSJed Brown   f0[0] = 1;     /* x direction pull */
20c4762a1bSJed Brown   f0[1] = -x[2]; /* add a twist around x-axis */
21c4762a1bSJed Brown   f0[2] = x[1];
22c4762a1bSJed Brown }
23c4762a1bSJed Brown 
24d71ae5a4SJacob Faibussowitsch static void f1_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
25d71ae5a4SJacob Faibussowitsch {
26c4762a1bSJed Brown   const PetscInt Ncomp = dim;
27c4762a1bSJed Brown   PetscInt       comp, d;
28c4762a1bSJed Brown   for (comp = 0; comp < Ncomp; ++comp) {
29ad540459SPierre Jolivet     for (d = 0; d < dim; ++d) f1[comp * dim + d] = 0.0;
30c4762a1bSJed Brown   }
31c4762a1bSJed Brown }
32c4762a1bSJed Brown 
33c4762a1bSJed Brown /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */
34d71ae5a4SJacob Faibussowitsch static void f1_u_3d_alpha(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
35d71ae5a4SJacob Faibussowitsch {
36c4762a1bSJed Brown   PetscReal trace, mu = s_mu, lambda = s_lambda, rad;
37c4762a1bSJed Brown   PetscInt  i, j;
38c4762a1bSJed Brown   for (i = 0, rad = 0.; i < dim; i++) {
39c4762a1bSJed Brown     PetscReal t = x[i];
40c4762a1bSJed Brown     rad += t * t;
41c4762a1bSJed Brown   }
42c4762a1bSJed Brown   rad = PetscSqrtReal(rad);
43c4762a1bSJed Brown   if (rad > 0.25) {
44c4762a1bSJed Brown     mu *= s_soft_alpha;
45c4762a1bSJed Brown     lambda *= s_soft_alpha; /* we could keep the bulk the same like rubberish */
46c4762a1bSJed Brown   }
47ad540459SPierre Jolivet   for (i = 0, trace = 0; i < dim; ++i) trace += PetscRealPart(u_x[i * dim + i]);
48c4762a1bSJed Brown   for (i = 0; i < dim; ++i) {
49ad540459SPierre Jolivet     for (j = 0; j < dim; ++j) f1[i * dim + j] = mu * (u_x[i * dim + j] + u_x[j * dim + i]);
50c4762a1bSJed Brown     f1[i * dim + i] += lambda * trace;
51c4762a1bSJed Brown   }
52c4762a1bSJed Brown }
53c4762a1bSJed Brown 
54c4762a1bSJed Brown /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */
55d71ae5a4SJacob Faibussowitsch static void f1_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
56d71ae5a4SJacob Faibussowitsch {
57c4762a1bSJed Brown   PetscReal trace, mu = s_mu, lambda = s_lambda;
58c4762a1bSJed Brown   PetscInt  i, j;
59ad540459SPierre Jolivet   for (i = 0, trace = 0; i < dim; ++i) trace += PetscRealPart(u_x[i * dim + i]);
60c4762a1bSJed Brown   for (i = 0; i < dim; ++i) {
61ad540459SPierre Jolivet     for (j = 0; j < dim; ++j) f1[i * dim + j] = mu * (u_x[i * dim + j] + u_x[j * dim + i]);
62c4762a1bSJed Brown     f1[i * dim + i] += lambda * trace;
63c4762a1bSJed Brown   }
64c4762a1bSJed Brown }
65c4762a1bSJed Brown 
663667b0a5SMark Adams static void f1_u_lap(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
673667b0a5SMark Adams {
683667b0a5SMark Adams   PetscInt d;
693667b0a5SMark Adams   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
703667b0a5SMark Adams }
713667b0a5SMark Adams 
72c4762a1bSJed Brown /* 3D elasticity */
73c4762a1bSJed Brown #define IDX(ii, jj, kk, ll) (27 * ii + 9 * jj + 3 * kk + ll)
74c4762a1bSJed Brown 
75d71ae5a4SJacob Faibussowitsch void g3_uu_3d_private(PetscScalar g3[], const PetscReal mu, const PetscReal lambda)
76d71ae5a4SJacob Faibussowitsch {
77c4762a1bSJed Brown   if (1) {
78c4762a1bSJed Brown     g3[0] += lambda;
79c4762a1bSJed Brown     g3[0] += mu;
80c4762a1bSJed Brown     g3[0] += mu;
81c4762a1bSJed Brown     g3[4] += lambda;
82c4762a1bSJed Brown     g3[8] += lambda;
83c4762a1bSJed Brown     g3[10] += mu;
84c4762a1bSJed Brown     g3[12] += mu;
85c4762a1bSJed Brown     g3[20] += mu;
86c4762a1bSJed Brown     g3[24] += mu;
87c4762a1bSJed Brown     g3[28] += mu;
88c4762a1bSJed Brown     g3[30] += mu;
89c4762a1bSJed Brown     g3[36] += lambda;
90c4762a1bSJed Brown     g3[40] += lambda;
91c4762a1bSJed Brown     g3[40] += mu;
92c4762a1bSJed Brown     g3[40] += mu;
93c4762a1bSJed Brown     g3[44] += lambda;
94c4762a1bSJed Brown     g3[50] += mu;
95c4762a1bSJed Brown     g3[52] += mu;
96c4762a1bSJed Brown     g3[56] += mu;
97c4762a1bSJed Brown     g3[60] += mu;
98c4762a1bSJed Brown     g3[68] += mu;
99c4762a1bSJed Brown     g3[70] += mu;
100c4762a1bSJed Brown     g3[72] += lambda;
101c4762a1bSJed Brown     g3[76] += lambda;
102c4762a1bSJed Brown     g3[80] += lambda;
103c4762a1bSJed Brown     g3[80] += mu;
104c4762a1bSJed Brown     g3[80] += mu;
105c4762a1bSJed Brown   } else {
106c4762a1bSJed Brown     int        i, j, k, l;
107c4762a1bSJed Brown     static int cc = -1;
108c4762a1bSJed Brown     cc++;
109c4762a1bSJed Brown     for (i = 0; i < 3; ++i) {
110c4762a1bSJed Brown       for (j = 0; j < 3; ++j) {
111c4762a1bSJed Brown         for (k = 0; k < 3; ++k) {
112c4762a1bSJed Brown           for (l = 0; l < 3; ++l) {
113c4762a1bSJed Brown             if (k == l && i == j) g3[IDX(i, j, k, l)] += lambda;
114c4762a1bSJed Brown             if (i == k && j == l) g3[IDX(i, j, k, l)] += mu;
115c4762a1bSJed Brown             if (i == l && j == k) g3[IDX(i, j, k, l)] += mu;
116c4762a1bSJed Brown             if (k == l && i == j && !cc) (void)PetscPrintf(PETSC_COMM_WORLD, "g3[%d] += lambda;\n", IDX(i, j, k, l));
117c4762a1bSJed Brown             if (i == k && j == l && !cc) (void)PetscPrintf(PETSC_COMM_WORLD, "g3[%d] += mu;\n", IDX(i, j, k, l));
118c4762a1bSJed Brown             if (i == l && j == k && !cc) (void)PetscPrintf(PETSC_COMM_WORLD, "g3[%d] += mu;\n", IDX(i, j, k, l));
119c4762a1bSJed Brown           }
120c4762a1bSJed Brown         }
121c4762a1bSJed Brown       }
122c4762a1bSJed Brown     }
123c4762a1bSJed Brown   }
124c4762a1bSJed Brown }
125c4762a1bSJed Brown 
126d71ae5a4SJacob Faibussowitsch static void g3_uu_3d_alpha(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
127d71ae5a4SJacob Faibussowitsch {
128c4762a1bSJed Brown   PetscReal mu = s_mu, lambda = s_lambda, rad;
129c4762a1bSJed Brown   PetscInt  i;
130c4762a1bSJed Brown   for (i = 0, rad = 0.; i < dim; i++) {
131c4762a1bSJed Brown     PetscReal t = x[i];
132c4762a1bSJed Brown     rad += t * t;
133c4762a1bSJed Brown   }
134c4762a1bSJed Brown   rad = PetscSqrtReal(rad);
135c4762a1bSJed Brown   if (rad > 0.25) {
136c4762a1bSJed Brown     mu *= s_soft_alpha;
137c4762a1bSJed Brown     lambda *= s_soft_alpha; /* we could keep the bulk the same like rubberish */
138c4762a1bSJed Brown   }
139c4762a1bSJed Brown   g3_uu_3d_private(g3, mu, lambda);
140c4762a1bSJed Brown }
141c4762a1bSJed Brown 
142d71ae5a4SJacob Faibussowitsch static void g3_uu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
143d71ae5a4SJacob Faibussowitsch {
144c4762a1bSJed Brown   g3_uu_3d_private(g3, s_mu, s_lambda);
145c4762a1bSJed Brown }
146c4762a1bSJed Brown 
1473667b0a5SMark Adams static void g3_lap(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
1483667b0a5SMark Adams {
1493667b0a5SMark Adams   PetscInt d;
1503667b0a5SMark Adams   for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
1513667b0a5SMark Adams }
1523667b0a5SMark Adams 
1533667b0a5SMark Adams static void g3_lap_alpha(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
1543667b0a5SMark Adams {
1553667b0a5SMark Adams   PetscReal lambda = 1, rad;
1563667b0a5SMark Adams   PetscInt  i;
1573667b0a5SMark Adams   for (i = 0, rad = 0.; i < dim; i++) {
1583667b0a5SMark Adams     PetscReal t = x[i];
1593667b0a5SMark Adams     rad += t * t;
1603667b0a5SMark Adams   }
1613667b0a5SMark Adams   rad = PetscSqrtReal(rad);
1623667b0a5SMark Adams   if (rad > 0.25) { lambda *= s_soft_alpha; /* we could keep the bulk the same like rubberish */ }
1633667b0a5SMark Adams   for (int d = 0; d < dim; ++d) g3[d * dim + d] = lambda;
1643667b0a5SMark Adams }
1653667b0a5SMark Adams 
166d71ae5a4SJacob Faibussowitsch static void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
167d71ae5a4SJacob Faibussowitsch {
168c4762a1bSJed Brown   const PetscInt Ncomp = dim;
169c4762a1bSJed Brown   PetscInt       comp;
170c4762a1bSJed Brown 
171c4762a1bSJed Brown   for (comp = 0; comp < Ncomp; ++comp) f0[comp] = 0.0;
172c4762a1bSJed Brown }
173c4762a1bSJed Brown 
174c4762a1bSJed Brown /* PI_i (x_i^4 - x_i^2) */
175d71ae5a4SJacob Faibussowitsch static void f0_u_x4(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
176d71ae5a4SJacob Faibussowitsch {
1773667b0a5SMark Adams   for (int comp = 0; comp < Nf; ++comp) {
178c4762a1bSJed Brown     f0[comp] = 1e5;
1793667b0a5SMark Adams     for (int i = 0; i < dim; ++i) { f0[comp] *= /* (comp+1)* */ (x[i] * x[i] * x[i] * x[i] - x[i] * x[i]); /* assumes (0,1]^D domain */ }
180c4762a1bSJed Brown   }
181c4762a1bSJed Brown }
182c4762a1bSJed Brown 
183d71ae5a4SJacob Faibussowitsch PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
184d71ae5a4SJacob Faibussowitsch {
185c4762a1bSJed Brown   const PetscInt Ncomp = dim;
186c4762a1bSJed Brown   PetscInt       comp;
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   for (comp = 0; comp < Ncomp; ++comp) u[comp] = 0;
1893ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
190c4762a1bSJed Brown }
191c4762a1bSJed Brown 
192d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
193d71ae5a4SJacob Faibussowitsch {
194c4762a1bSJed Brown   Mat           Amat;
195c4762a1bSJed Brown   SNES          snes;
196c4762a1bSJed Brown   KSP           ksp;
197c4762a1bSJed Brown   MPI_Comm      comm;
198c4762a1bSJed Brown   PetscMPIInt   rank;
199c4762a1bSJed Brown   PetscLogStage stage[17];
200c4762a1bSJed Brown   PetscBool     test_nonzero_cols = PETSC_FALSE, use_nearnullspace = PETSC_TRUE, attach_nearnullspace = PETSC_FALSE;
201c4762a1bSJed Brown   Vec           xx, bb;
2023667b0a5SMark Adams   PetscInt      iter, i, N, dim = 3, max_conv_its, sizes[7], run_type = 1, Ncomp = dim;
2033667b0a5SMark Adams   DM            dm;
204c4762a1bSJed Brown   PetscBool     flg;
205c4762a1bSJed Brown   PetscReal     Lx, mdisp[10], err[10];
2064d86920dSPierre Jolivet 
207327415f7SBarry Smith   PetscFunctionBeginUser;
2089566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
209c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
2109566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
211c4762a1bSJed Brown   /* options */
212d0609cedSBarry Smith   PetscOptionsBegin(comm, NULL, "3D bilinear Q1 elasticity options", "");
213c4762a1bSJed Brown   {
214c4762a1bSJed Brown     Lx           = 1.; /* or ne for rod */
215c4762a1bSJed Brown     max_conv_its = 3;
2169566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-max_conv_its", "Number of iterations in convergence study", "", max_conv_its, &max_conv_its, NULL));
2173667b0a5SMark Adams     PetscCheck(max_conv_its > 0 && max_conv_its < 8, PETSC_COMM_WORLD, PETSC_ERR_USER, "Bad number of iterations for convergence test (%" PetscInt_FMT ")", max_conv_its);
2189566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-lx", "Length of domain", "", Lx, &Lx, NULL));
2199566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-alpha", "material coefficient inside circle", "", s_soft_alpha, &s_soft_alpha, NULL));
2209566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-test_nonzero_cols", "nonzero test", "", test_nonzero_cols, &test_nonzero_cols, NULL));
2219566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-use_mat_nearnullspace", "MatNearNullSpace API test", "", use_nearnullspace, &use_nearnullspace, NULL));
2229566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-attach_mat_nearnullspace", "MatNearNullSpace API test (via MatSetNearNullSpace)", "", attach_nearnullspace, &attach_nearnullspace, NULL));
2233667b0a5SMark Adams     PetscCall(PetscOptionsInt("-run_type", "0: twisting load on cantalever, 1: Elasticty convergence test on cube, 2: Laplacian, 3: soft core Laplacian", "", run_type, &run_type, NULL));
224c4762a1bSJed Brown   }
225d0609cedSBarry Smith   PetscOptionsEnd();
2269566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("Mesh Setup", &stage[16]));
227c4762a1bSJed Brown   for (iter = 0; iter < max_conv_its; iter++) {
228c4762a1bSJed Brown     char str[] = "Solve 0";
229c4762a1bSJed Brown     str[6] += iter;
2309566063dSJacob Faibussowitsch     PetscCall(PetscLogStageRegister(str, &stage[iter]));
231c4762a1bSJed Brown   }
232c4762a1bSJed Brown   /* create DM, Plex calls DMSetup */
2339566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(stage[16]));
2343667b0a5SMark Adams   PetscCall(DMCreate(comm, &dm));
2353667b0a5SMark Adams   PetscCall(DMSetType(dm, DMPLEX));
2363667b0a5SMark Adams   PetscCall(PetscObjectSetName((PetscObject)dm, "Mesh"));
2373667b0a5SMark Adams   PetscCall(DMSetFromOptions(dm));
2383667b0a5SMark Adams   PetscCall(DMPlexDistributeSetDefault(dm, PETSC_FALSE));
2393667b0a5SMark Adams   PetscCall(DMGetDimension(dm, &dim));
240c4762a1bSJed Brown   {
241c4762a1bSJed Brown     DMLabel label;
242c4762a1bSJed Brown     IS      is;
2439566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(dm, "boundary"));
2449566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dm, "boundary", &label));
2459566063dSJacob Faibussowitsch     PetscCall(DMPlexMarkBoundaryFaces(dm, 1, label));
2462a5e0493SMark Adams     if (run_type == 0) {
2479566063dSJacob Faibussowitsch       PetscCall(DMGetStratumIS(dm, "boundary", 1, &is));
2489566063dSJacob Faibussowitsch       PetscCall(DMCreateLabel(dm, "Faces"));
249c4762a1bSJed Brown       if (is) {
250c4762a1bSJed Brown         PetscInt        d, f, Nf;
251c4762a1bSJed Brown         const PetscInt *faces;
252c4762a1bSJed Brown         PetscInt        csize;
253c4762a1bSJed Brown         PetscSection    cs;
254c4762a1bSJed Brown         Vec             coordinates;
255c4762a1bSJed Brown         DM              cdm;
2569566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(is, &Nf));
2579566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(is, &faces));
2589566063dSJacob Faibussowitsch         PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
2599566063dSJacob Faibussowitsch         PetscCall(DMGetCoordinateDM(dm, &cdm));
2609566063dSJacob Faibussowitsch         PetscCall(DMGetLocalSection(cdm, &cs));
261c4762a1bSJed Brown         /* Check for each boundary face if any component of its centroid is either 0.0 or 1.0 */
262c4762a1bSJed Brown         for (f = 0; f < Nf; ++f) {
263c4762a1bSJed Brown           PetscReal    faceCoord;
264c4762a1bSJed Brown           PetscInt     b, v;
265c4762a1bSJed Brown           PetscScalar *coords = NULL;
266c4762a1bSJed Brown           PetscInt     Nv;
2679566063dSJacob Faibussowitsch           PetscCall(DMPlexVecGetClosure(cdm, cs, coordinates, faces[f], &csize, &coords));
268c4762a1bSJed Brown           Nv = csize / dim; /* Calculate mean coordinate vector */
269c4762a1bSJed Brown           for (d = 0; d < dim; ++d) {
270c4762a1bSJed Brown             faceCoord = 0.0;
271c4762a1bSJed Brown             for (v = 0; v < Nv; ++v) faceCoord += PetscRealPart(coords[v * dim + d]);
272c4762a1bSJed Brown             faceCoord /= Nv;
273c4762a1bSJed Brown             for (b = 0; b < 2; ++b) {
274c4762a1bSJed Brown               if (PetscAbs(faceCoord - b) < PETSC_SMALL) { /* domain have not been set yet, still [0,1]^3 */
2759566063dSJacob Faibussowitsch                 PetscCall(DMSetLabelValue(dm, "Faces", faces[f], d * 2 + b + 1));
276c4762a1bSJed Brown               }
277c4762a1bSJed Brown             }
278c4762a1bSJed Brown           }
2799566063dSJacob Faibussowitsch           PetscCall(DMPlexVecRestoreClosure(cdm, cs, coordinates, faces[f], &csize, &coords));
280c4762a1bSJed Brown         }
2819566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(is, &faces));
282c4762a1bSJed Brown       }
2839566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&is));
2849566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(dm, "Faces", &label));
2859566063dSJacob Faibussowitsch       PetscCall(DMPlexLabelComplete(dm, label));
286c4762a1bSJed Brown     }
287c4762a1bSJed Brown   }
2889566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
289c4762a1bSJed Brown   for (iter = 0; iter < max_conv_its; iter++) {
2909566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(stage[16]));
291c4762a1bSJed Brown     /* snes */
2929566063dSJacob Faibussowitsch     PetscCall(SNESCreate(comm, &snes));
2939566063dSJacob Faibussowitsch     PetscCall(SNESSetDM(snes, dm));
2943667b0a5SMark Adams     PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
295c4762a1bSJed Brown     /* fem */
296c4762a1bSJed Brown     {
297c4762a1bSJed Brown       const PetscInt components[] = {0, 1, 2};
298c4762a1bSJed Brown       const PetscInt Nfid = 1, Npid = 1;
2993667b0a5SMark Adams       PetscInt       fid[] = {1}; /* The fixed faces (x=0) */
300c4762a1bSJed Brown       const PetscInt pid[] = {2}; /* The faces with loading (x=L_x) */
301c4762a1bSJed Brown       PetscFE        fe;
302c4762a1bSJed Brown       PetscDS        prob;
30345480ffeSMatthew G. Knepley       DMLabel        label;
304c4762a1bSJed Brown 
3053667b0a5SMark Adams       if (run_type == 2 || run_type == 3) Ncomp = 1;
3063667b0a5SMark Adams       else Ncomp = dim;
3073667b0a5SMark Adams       PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, Ncomp, PETSC_FALSE, NULL, PETSC_DECIDE, &fe));
3089566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)fe, "deformation"));
309c4762a1bSJed Brown       /* FEM prob */
3109566063dSJacob Faibussowitsch       PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
3119566063dSJacob Faibussowitsch       PetscCall(DMCreateDS(dm));
3129566063dSJacob Faibussowitsch       PetscCall(DMGetDS(dm, &prob));
313c4762a1bSJed Brown       /* setup problem */
3143667b0a5SMark Adams       if (run_type == 1) { // elast
3159566063dSJacob Faibussowitsch         PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu_3d));
3169566063dSJacob Faibussowitsch         PetscCall(PetscDSSetResidual(prob, 0, f0_u_x4, f1_u_3d));
3173667b0a5SMark Adams       } else if (run_type == 0) { //twisted not maintained
31845480ffeSMatthew G. Knepley         PetscWeakForm wf;
31945480ffeSMatthew G. Knepley         PetscInt      bd, i;
3209566063dSJacob Faibussowitsch         PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu_3d_alpha));
3219566063dSJacob Faibussowitsch         PetscCall(PetscDSSetResidual(prob, 0, f0_u, f1_u_3d_alpha));
3229566063dSJacob Faibussowitsch         PetscCall(DMGetLabel(dm, "Faces", &label));
3239566063dSJacob Faibussowitsch         PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "traction", label, Npid, pid, 0, Ncomp, components, NULL, NULL, NULL, &bd));
3249566063dSJacob Faibussowitsch         PetscCall(PetscDSGetBoundary(prob, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
3259566063dSJacob Faibussowitsch         for (i = 0; i < Npid; ++i) PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, pid[i], 0, 0, 0, f0_bd_u_3d, 0, f1_bd_u));
3263667b0a5SMark Adams       } else if (run_type == 2) { // Laplacian
3273667b0a5SMark Adams         PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_lap));
3283667b0a5SMark Adams         PetscCall(PetscDSSetResidual(prob, 0, f0_u_x4, f1_u_lap));
3293667b0a5SMark Adams       } else if (run_type == 3) { // soft core Laplacian
3303667b0a5SMark Adams         PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_lap_alpha));
3313667b0a5SMark Adams         PetscCall(PetscDSSetResidual(prob, 0, f0_u_x4, f1_u_lap));
332c4762a1bSJed Brown       }
333c4762a1bSJed Brown       /* bcs */
3343667b0a5SMark Adams       if (run_type != 0) {
335c4762a1bSJed Brown         PetscInt id = 1;
3369566063dSJacob Faibussowitsch         PetscCall(DMGetLabel(dm, "boundary", &label));
3379566063dSJacob Faibussowitsch         PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))zero, NULL, NULL, NULL));
338c4762a1bSJed Brown       } else {
3399566063dSJacob Faibussowitsch         PetscCall(DMGetLabel(dm, "Faces", &label));
3409566063dSJacob Faibussowitsch         PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "fixed", label, Nfid, fid, 0, Ncomp, components, (void (*)(void))zero, NULL, NULL, NULL));
341c4762a1bSJed Brown       }
3429566063dSJacob Faibussowitsch       PetscCall(PetscFEDestroy(&fe));
343c4762a1bSJed Brown     }
344c4762a1bSJed Brown     /* vecs & mat */
3459566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(dm, &xx));
3469566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx, &bb));
3479566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)bb, "b"));
3489566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)xx, "u"));
3499566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix(dm, &Amat));
3509566063dSJacob Faibussowitsch     PetscCall(MatSetOption(Amat, MAT_SYMMETRIC, PETSC_TRUE));        /* Some matrix kernels can take advantage of symmetry if we set this. */
3519566063dSJacob Faibussowitsch     PetscCall(MatSetOption(Amat, MAT_SYMMETRY_ETERNAL, PETSC_TRUE)); /* Inform PETSc that Amat is always symmetric, so info set above isn't lost. */
3523667b0a5SMark Adams     PetscCall(MatSetBlockSize(Amat, Ncomp));
3539566063dSJacob Faibussowitsch     PetscCall(MatSetOption(Amat, MAT_SPD, PETSC_TRUE));
354b94d7dedSBarry Smith     PetscCall(MatSetOption(Amat, MAT_SPD_ETERNAL, PETSC_TRUE));
3559566063dSJacob Faibussowitsch     PetscCall(VecGetSize(bb, &N));
3563667b0a5SMark Adams     sizes[iter] = N;
35763a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "%" PetscInt_FMT " global equations, %" PetscInt_FMT " vertices\n", N, N / dim));
3583667b0a5SMark Adams     if ((use_nearnullspace || attach_nearnullspace) && N / dim > 1 && Ncomp > 1) {
359c4762a1bSJed Brown       /* Set up the near null space (a.k.a. rigid body modes) that will be used by the multigrid preconditioner */
360c4762a1bSJed Brown       DM           subdm;
361c4762a1bSJed Brown       MatNullSpace nearNullSpace;
362c4762a1bSJed Brown       PetscInt     fields = 0;
363c4762a1bSJed Brown       PetscObject  deformation;
3649566063dSJacob Faibussowitsch       PetscCall(DMCreateSubDM(dm, 1, &fields, NULL, &subdm));
3659566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateRigidBody(subdm, 0, &nearNullSpace));
3669566063dSJacob Faibussowitsch       PetscCall(DMGetField(dm, 0, NULL, &deformation));
3679566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose(deformation, "nearnullspace", (PetscObject)nearNullSpace));
3689566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&subdm));
3691baa6e33SBarry Smith       if (attach_nearnullspace) PetscCall(MatSetNearNullSpace(Amat, nearNullSpace));
3709566063dSJacob Faibussowitsch       PetscCall(MatNullSpaceDestroy(&nearNullSpace)); /* created by DM and destroyed by Mat */
371c4762a1bSJed Brown     }
3726493148fSStefano Zampini     PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, NULL));
3739566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, Amat, Amat, NULL, NULL));
3749566063dSJacob Faibussowitsch     PetscCall(SNESSetFromOptions(snes));
3759566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dm));
3769566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
3779566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(stage[16]));
378c4762a1bSJed Brown     /* ksp */
3799566063dSJacob Faibussowitsch     PetscCall(SNESGetKSP(snes, &ksp));
3809566063dSJacob Faibussowitsch     PetscCall(KSPSetComputeSingularValues(ksp, PETSC_TRUE));
3818926f930SMark Adams     if (!use_nearnullspace) {
3828926f930SMark Adams       PC pc;
3838926f930SMark Adams       PetscCall(KSPGetPC(ksp, &pc));
3848926f930SMark Adams       PetscCall(PCGAMGASMSetHEM(pc, 3)); // code coverage
3858926f930SMark Adams     }
386c4762a1bSJed Brown     /* test BCs */
3879566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(xx));
388c4762a1bSJed Brown     if (test_nonzero_cols) {
38948a46eb9SPierre Jolivet       if (rank == 0) PetscCall(VecSetValue(xx, 0, 1.0, INSERT_VALUES));
3909566063dSJacob Faibussowitsch       PetscCall(VecAssemblyBegin(xx));
3919566063dSJacob Faibussowitsch       PetscCall(VecAssemblyEnd(xx));
392c4762a1bSJed Brown     }
3939566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(bb));
3949566063dSJacob Faibussowitsch     PetscCall(VecGetSize(bb, &i));
3953667b0a5SMark Adams     sizes[iter] = i;
39663a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "%" PetscInt_FMT " equations in vector, %" PetscInt_FMT " vertices\n", i, i / dim));
3979566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
398c4762a1bSJed Brown     /* solve */
3993667b0a5SMark Adams     PetscCall(SNESComputeJacobian(snes, xx, Amat, Amat));
4003667b0a5SMark Adams     PetscCall(MatViewFromOptions(Amat, NULL, "-my_mat_view"));
4019566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(stage[iter]));
4029566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes, bb, xx));
4039566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
4049566063dSJacob Faibussowitsch     PetscCall(VecNorm(xx, NORM_INFINITY, &mdisp[iter]));
405c4762a1bSJed Brown     {
406c4762a1bSJed Brown       PetscViewer       viewer = NULL;
407c4762a1bSJed Brown       PetscViewerFormat fmt;
4083667b0a5SMark Adams       PetscCall(PetscOptionsGetViewer(comm, NULL, "", "-vec_view", &viewer, &fmt, &flg));
409c4762a1bSJed Brown       if (flg) {
4109566063dSJacob Faibussowitsch         PetscCall(PetscViewerPushFormat(viewer, fmt));
4119566063dSJacob Faibussowitsch         PetscCall(VecView(xx, viewer));
4129566063dSJacob Faibussowitsch         PetscCall(VecView(bb, viewer));
4139566063dSJacob Faibussowitsch         PetscCall(PetscViewerPopFormat(viewer));
414c4762a1bSJed Brown       }
415cd791dc2SBarry Smith       PetscCall(PetscOptionsRestoreViewer(&viewer));
416c4762a1bSJed Brown     }
417c4762a1bSJed Brown     /* Free work space */
4189566063dSJacob Faibussowitsch     PetscCall(SNESDestroy(&snes));
4199566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xx));
4209566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bb));
4219566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Amat));
4223667b0a5SMark Adams     if (iter + 1 < max_conv_its) {
4233667b0a5SMark Adams       DM newdm;
4243667b0a5SMark Adams       PetscCall(DMViewFromOptions(dm, NULL, "-my_dm_view"));
4253667b0a5SMark Adams       PetscCall(DMRefine(dm, comm, &newdm));
4263667b0a5SMark Adams       if (rank == -1) {
4273667b0a5SMark Adams         PetscDS prob;
4283667b0a5SMark Adams         PetscCall(DMGetDS(dm, &prob));
4293667b0a5SMark Adams         PetscCall(PetscDSViewFromOptions(prob, NULL, "-ds_view"));
4303667b0a5SMark Adams         PetscCall(DMGetDS(newdm, &prob));
4313667b0a5SMark Adams         PetscCall(PetscDSViewFromOptions(prob, NULL, "-ds_view"));
432c4762a1bSJed Brown       }
4333667b0a5SMark Adams       PetscCall(DMDestroy(&dm));
4343667b0a5SMark Adams       dm = newdm;
4353667b0a5SMark Adams       PetscCall(PetscObjectSetName((PetscObject)dm, "Mesh"));
4363667b0a5SMark Adams       PetscCall(DMViewFromOptions(dm, NULL, "-my_dm_view"));
4373667b0a5SMark Adams       PetscCall(DMSetFromOptions(dm));
4383667b0a5SMark Adams     }
4393667b0a5SMark Adams   }
4403667b0a5SMark Adams   PetscCall(DMDestroy(&dm));
441cdb485ddSMark Adams   if (run_type == 1) err[0] = 5.97537599375e+01 - mdisp[0]; /* error with what I think is the exact solution */
4423667b0a5SMark Adams   else if (run_type == 0) err[0] = 0;
4433667b0a5SMark Adams   else if (run_type == 2) err[0] = 3.527795e+01 - mdisp[0];
4443667b0a5SMark Adams   else err[0] = 0;
4453667b0a5SMark Adams   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[%d] %d) N=%12" PetscInt_FMT ", max displ=%9.7e, error=%4.3e\n", rank, 0, sizes[0], (double)mdisp[0], (double)err[0]));
446c4762a1bSJed Brown   for (iter = 1; iter < max_conv_its; iter++) {
447cdb485ddSMark Adams     if (run_type == 1) err[iter] = 5.97537599375e+01 - mdisp[iter];
4483667b0a5SMark Adams     else if (run_type == 0) err[iter] = 0;
4493667b0a5SMark Adams     else if (run_type == 2) err[iter] = 3.527795e+01 - mdisp[iter];
4503667b0a5SMark Adams     else err[iter] = 0;
4513667b0a5SMark Adams     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[%d] %" PetscInt_FMT ") N=%12" PetscInt_FMT ", max displ=%9.7e, disp diff=%9.2e, error=%4.3e, rate=%3.2g\n", rank, iter, sizes[iter], (double)mdisp[iter], (double)(mdisp[iter] - mdisp[iter - 1]), (double)err[iter], (double)(PetscLogReal(PetscAbs(err[iter - 1] / err[iter])) / PetscLogReal(2.))));
452c4762a1bSJed Brown   }
453c4762a1bSJed Brown 
4549566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
455b122ec5aSJacob Faibussowitsch   return 0;
456c4762a1bSJed Brown }
457c4762a1bSJed Brown 
458c4762a1bSJed Brown /*TEST
459c4762a1bSJed Brown 
4603667b0a5SMark Adams   testset:
461c4762a1bSJed Brown     nsize: 4
462c4762a1bSJed Brown     requires: !single
463d4adc10fSMark Adams     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 2,2,1 -petscspace_degree 2 -snes_max_it 1 -ksp_max_it 100 -ksp_type cg -ksp_rtol 1.e-10 -ksp_norm_type unpreconditioned -pc_type gamg -pc_gamg_coarse_eq_limit 10 -pc_gamg_reuse_interpolation true -pc_gamg_aggressive_coarsening 1 -pc_gamg_threshold 0.001 -ksp_converged_reason -use_mat_nearnullspace true -mg_levels_ksp_max_it 2 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.2,0,1.1 -mg_levels_pc_type jacobi -petscpartitioner_type simple -my_dm_view -snes_lag_jacobian -2 -snes_type ksponly -pc_gamg_mis_k_minimum_degree_ordering true -pc_gamg_low_memory_threshold_filter
464c4762a1bSJed Brown     timeoutfactor: 2
4653667b0a5SMark Adams     test:
4663667b0a5SMark Adams       suffix: 0
467*e0b7e82fSBarry Smith       args: -run_type 1 -max_conv_its 3 -pc_gamg_mat_coarsen_type hem -pc_gamg_mat_coarsen_max_it 5 -pc_gamg_asm_hem_aggs 4 -ksp_rtol 1.e-6
4688926f930SMark Adams       filter: sed -e "s/Linear solve converged due to CONVERGED_RTOL iterations 7/Linear solve converged due to CONVERGED_RTOL iterations 8/g"
4693667b0a5SMark Adams     test:
4703667b0a5SMark Adams       suffix: 1
4718926f930SMark Adams       filter: grep -v HERMITIAN
4728926f930SMark Adams       args: -run_type 2 -max_conv_its 2 -use_mat_nearnullspace false -snes_view
473d529f056Smarkadams4 
474d529f056Smarkadams4   test:
475d529f056Smarkadams4     nsize: 1
476d529f056Smarkadams4     requires: !single
477d529f056Smarkadams4     suffix: 2
478d4adc10fSMark Adams     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 2,2,1 -max_conv_its 2 -petscspace_degree 2 -snes_max_it 1 -ksp_type cg -ksp_norm_type unpreconditioned -pc_type gamg -pc_gamg_coarse_eq_limit 10 -pc_gamg_aggressive_coarsening 1 -ksp_converged_reason -use_mat_nearnullspace true -my_dm_view -snes_type ksponly
479d529f056Smarkadams4     timeoutfactor: 2
480c4762a1bSJed Brown 
481c4762a1bSJed Brown   # HYPRE PtAP broken with complex numbers
482c4762a1bSJed Brown   test:
483c4762a1bSJed Brown     suffix: hypre
484263f2b91SStefano Zampini     requires: hypre !single !complex !defined(PETSC_HAVE_HYPRE_DEVICE)
485c4762a1bSJed Brown     nsize: 4
4863667b0a5SMark Adams     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-8 -pc_type hypre -pc_hypre_type boomeramg -pc_hypre_boomeramg_no_CF true -pc_hypre_boomeramg_agg_nl 1 -pc_hypre_boomeramg_coarsen_type HMIS -pc_hypre_boomeramg_interp_type ext+i -ksp_converged_reason -use_mat_nearnullspace true -petscpartitioner_type simple
487c4762a1bSJed Brown 
488c4762a1bSJed Brown   test:
489c4762a1bSJed Brown     suffix: ml
490c4762a1bSJed Brown     requires: ml !single
491c4762a1bSJed Brown     nsize: 4
4923667b0a5SMark Adams     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type cg -ksp_monitor_short -ksp_converged_reason -ksp_rtol 1.e-8 -pc_type ml -mg_levels_ksp_type chebyshev -mg_levels_ksp_max_it 3 -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 -mg_levels_pc_type sor -petscpartitioner_type simple -use_mat_nearnullspace
493c4762a1bSJed Brown 
494c4762a1bSJed Brown   test:
495c4762a1bSJed Brown     suffix: hpddm
496dfd57a17SPierre Jolivet     requires: hpddm slepc !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
497c4762a1bSJed Brown     nsize: 4
4983667b0a5SMark Adams     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type fgmres -ksp_monitor_short -ksp_converged_reason -ksp_rtol 1.e-8 -pc_type hpddm -petscpartitioner_type simple -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_nev 6 -pc_hpddm_coarse_p 1 -pc_hpddm_coarse_pc_type svd
499c4762a1bSJed Brown 
500c4762a1bSJed Brown   test:
50163b77682SMark Adams     suffix: repart
502c4762a1bSJed Brown     nsize: 4
503c4762a1bSJed Brown     requires: parmetis !single
5043667b0a5SMark Adams     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 2,2,1 -max_conv_its 2 -petscspace_degree 2 -snes_max_it 4 -ksp_max_it 100 -ksp_type cg -ksp_rtol 1.e-2 -ksp_norm_type unpreconditioned -snes_rtol 1.e-3 -pc_type gamg -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -pc_gamg_aggressive_coarsening 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -use_mat_nearnullspace true -mg_levels_ksp_max_it 2 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 -mg_levels_pc_type jacobi -pc_gamg_mat_partitioning_type parmetis -pc_gamg_repartition true -pc_gamg_process_eq_limit 20 -pc_gamg_coarse_eq_limit 10 -ksp_converged_reason -pc_gamg_reuse_interpolation true -petscpartitioner_type simple
505c4762a1bSJed Brown 
506c4762a1bSJed Brown   test:
507c4762a1bSJed Brown     suffix: bddc
508c4762a1bSJed Brown     nsize: 4
509c4762a1bSJed Brown     requires: !single
5104f58015eSStefano Zampini     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-8 -ksp_converged_reason -petscpartitioner_type simple -dm_mat_type is -mat_is_localmat_type {{sbaij baij aij}} -pc_type bddc
511c4762a1bSJed Brown 
512c4762a1bSJed Brown   testset:
513c4762a1bSJed Brown     nsize: 4
514c4762a1bSJed Brown     requires: !single
5154f58015eSStefano Zampini     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-10 -ksp_converged_reason -petscpartitioner_type simple -dm_mat_type is -mat_is_localmat_type aij -pc_type bddc -attach_mat_nearnullspace {{0 1}separate output}
516c4762a1bSJed Brown     test:
517c4762a1bSJed Brown       suffix: bddc_approx_gamg
518bae903cbSmarkadams4       args: -pc_bddc_switch_static -prefix_push pc_bddc_dirichlet_ -approximate -pc_type gamg -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -pc_gamg_reuse_interpolation true -pc_gamg_aggressive_coarsening 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -prefix_pop -prefix_push pc_bddc_neumann_ -approximate -pc_type gamg -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -pc_gamg_coarse_eq_limit 10 -pc_gamg_reuse_interpolation true -pc_gamg_aggressive_coarsening 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -prefix_pop
519c4762a1bSJed Brown     # HYPRE PtAP broken with complex numbers
520c4762a1bSJed Brown     test:
521263f2b91SStefano Zampini       requires: hypre !complex !defined(PETSC_HAVE_HYPRE_DEVICE)
522c4762a1bSJed Brown       suffix: bddc_approx_hypre
523c4762a1bSJed Brown       args: -pc_bddc_switch_static -prefix_push pc_bddc_dirichlet_ -pc_type hypre -pc_hypre_boomeramg_no_CF true -pc_hypre_boomeramg_strong_threshold 0.75 -pc_hypre_boomeramg_agg_nl 1 -pc_hypre_boomeramg_coarsen_type HMIS -pc_hypre_boomeramg_interp_type ext+i -prefix_pop -prefix_push pc_bddc_neumann_ -pc_type hypre -pc_hypre_boomeramg_no_CF true -pc_hypre_boomeramg_strong_threshold 0.75 -pc_hypre_boomeramg_agg_nl 1 -pc_hypre_boomeramg_coarsen_type HMIS -pc_hypre_boomeramg_interp_type ext+i -prefix_pop
524c4762a1bSJed Brown     test:
525c4762a1bSJed Brown       requires: ml
526c4762a1bSJed Brown       suffix: bddc_approx_ml
5270e75e42fSMark       args: -pc_bddc_switch_static -prefix_push pc_bddc_dirichlet_ -approximate -pc_type ml -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -prefix_pop -prefix_push pc_bddc_neumann_ -approximate -pc_type ml -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -prefix_pop
528c4762a1bSJed Brown 
529c4762a1bSJed Brown   test:
530c4762a1bSJed Brown     suffix: fetidp
531c4762a1bSJed Brown     nsize: 4
532c4762a1bSJed Brown     requires: !single
5334f58015eSStefano Zampini     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type fetidp -fetidp_ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-8 -ksp_converged_reason -petscpartitioner_type simple -dm_mat_type is -mat_is_localmat_type {{sbaij baij aij}}
534c4762a1bSJed Brown 
535c4762a1bSJed Brown   test:
536c4762a1bSJed Brown     suffix: bddc_elast
537c4762a1bSJed Brown     nsize: 4
538c4762a1bSJed Brown     requires: !single
5394f58015eSStefano Zampini     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-8 -ksp_converged_reason -petscpartitioner_type simple -dm_mat_type is -mat_is_localmat_type sbaij -pc_type bddc -pc_bddc_monolithic -attach_mat_nearnullspace
540c4762a1bSJed Brown 
541c4762a1bSJed Brown   test:
542c4762a1bSJed Brown     suffix: fetidp_elast
543c4762a1bSJed Brown     nsize: 4
544c4762a1bSJed Brown     requires: !single
5454f58015eSStefano Zampini     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type fetidp -fetidp_ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-8 -ksp_converged_reason -petscpartitioner_type simple -dm_mat_type is -mat_is_localmat_type sbaij -fetidp_bddc_pc_bddc_monolithic -attach_mat_nearnullspace
546c4762a1bSJed Brown 
547908b9b43SStefano Zampini   test:
548908b9b43SStefano Zampini     suffix: gdsw
549908b9b43SStefano Zampini     nsize: 4
550908b9b43SStefano Zampini     requires: !single
5513667b0a5SMark Adams     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-8 -ksp_converged_reason -petscpartitioner_type simple -dm_mat_type is -attach_mat_nearnullspace \
552908b9b43SStefano Zampini           -pc_type mg -pc_mg_galerkin -pc_mg_adapt_interp_coarse_space gdsw -pc_mg_levels 2 -mg_levels_pc_type bjacobi -mg_levels_sub_pc_type icc
553908b9b43SStefano Zampini 
55435990778SJunchao Zhang   testset:
55535990778SJunchao Zhang     nsize: 4
55635990778SJunchao Zhang     requires: !single
557d4adc10fSMark Adams     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 2,2,1 -max_conv_its 2 -petscspace_degree 2 -snes_max_it 2 -ksp_max_it 100 -ksp_type cg -ksp_rtol 1.e-10 -ksp_norm_type unpreconditioned -snes_rtol 1.e-10 -pc_type gamg -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -pc_gamg_coarse_eq_limit 10 -pc_gamg_reuse_interpolation true -pc_gamg_aggressive_coarsening 0 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -use_mat_nearnullspace true -mg_levels_ksp_max_it 2 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 -mg_levels_pc_type jacobi -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -dm_view -petscpartitioner_type simple -pc_gamg_process_eq_limit 20 -pc_gamg_coarse_eq_limit 40
55835990778SJunchao Zhang     output_file: output/ex56_cuda.out
55935990778SJunchao Zhang 
560c4762a1bSJed Brown     test:
561c4762a1bSJed Brown       suffix: cuda
56235990778SJunchao Zhang       requires: cuda
5633667b0a5SMark Adams       args: -dm_mat_type aijcusparse -dm_vec_type cuda
5646cb74228SMark Adams 
5656cb74228SMark Adams     test:
56647d993e7Ssuyashtn       suffix: hip
56747d993e7Ssuyashtn       requires: hip
5683667b0a5SMark Adams       args: -dm_mat_type aijhipsparse -dm_vec_type hip
56947d993e7Ssuyashtn 
57047d993e7Ssuyashtn     test:
5716cb74228SMark Adams       suffix: viennacl
57235990778SJunchao Zhang       requires: viennacl
5733667b0a5SMark Adams       args: -dm_mat_type aijviennacl -dm_vec_type viennacl
5746cb74228SMark Adams 
57535990778SJunchao Zhang     test:
57635990778SJunchao Zhang       suffix: kokkos
577dcfd994dSJunchao Zhang       requires: kokkos_kernels
5783667b0a5SMark Adams       args: -dm_mat_type aijkokkos -dm_vec_type kokkos
579dea3b165SRichard Tran Mills   # Don't run AIJMKL caes with complex scalars because of convergence issues.
580dea3b165SRichard Tran Mills   # Note that we need to test both single and multiple MPI rank cases, because these use different sparse MKL routines to implement the PtAP operation.
581a8736bf8SRichard Tran Mills   test:
582a8736bf8SRichard Tran Mills     suffix: seqaijmkl
583a8736bf8SRichard Tran Mills     nsize: 1
584dfd57a17SPierre Jolivet     requires: defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE) !single !complex
585a25662f8SPablo Brubeck     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 2,2,1 -max_conv_its 2 -petscspace_degree 2 -snes_max_it 2 -ksp_max_it 100 -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned -snes_rtol 1.e-10 -pc_type gamg -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -pc_gamg_coarse_eq_limit 1000 -pc_gamg_reuse_interpolation true -pc_gamg_aggressive_coarsening 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -ksp_converged_reason -use_mat_nearnullspace true -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.1 -mg_levels_pc_type jacobi -petscpartitioner_type simple -mat_block_size 3 -dm_view -mat_seqaij_type seqaijmkl
586a8736bf8SRichard Tran Mills     timeoutfactor: 2
587a8736bf8SRichard Tran Mills 
588dea3b165SRichard Tran Mills   test:
589dea3b165SRichard Tran Mills     suffix: mpiaijmkl
5903667b0a5SMark Adams     nsize: 4
591dfd57a17SPierre Jolivet     requires: defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE) !single !complex
592a25662f8SPablo Brubeck     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 2,2,1 -max_conv_its 2 -petscspace_degree 2 -snes_max_it 2 -ksp_max_it 100 -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned -snes_rtol 1.e-10 -pc_type gamg -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -pc_gamg_coarse_eq_limit 1000 -pc_gamg_reuse_interpolation true -pc_gamg_aggressive_coarsening 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -ksp_converged_reason -use_mat_nearnullspace true -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.1 -mg_levels_pc_type jacobi -petscpartitioner_type simple -mat_block_size 3 -dm_view -mat_seqaij_type seqaijmkl
593dea3b165SRichard Tran Mills     timeoutfactor: 2
594dea3b165SRichard Tran Mills 
595c4762a1bSJed Brown TEST*/
596