xref: /petsc/src/snes/tutorials/ex56.c (revision 2a5e0493c6ebd032c99d5ce63f71fde6987c27ef)
1*2a5e0493SMark Adams static char help[] = "3D, tensor hexahedra (Q1-K), displacement finite element formulation\n\
2c4762a1bSJed Brown of linear elasticity.  E=1.0, nu=1/3.\n\
3c4762a1bSJed Brown Unit cube domain with Dirichlet boundary\n\n";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown #include <petscdmplex.h>
6c4762a1bSJed Brown #include <petscsnes.h>
7c4762a1bSJed Brown #include <petscds.h>
8c4762a1bSJed Brown #include <petscdmforest.h>
9c4762a1bSJed Brown 
10c4762a1bSJed Brown static PetscReal s_soft_alpha=1.e-3;
11c4762a1bSJed Brown static PetscReal s_mu=0.4;
12c4762a1bSJed Brown static PetscReal s_lambda=0.4;
13c4762a1bSJed Brown 
14c4762a1bSJed Brown static void f0_bd_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
15c4762a1bSJed Brown                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
16c4762a1bSJed Brown                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
17c4762a1bSJed Brown                        PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
18c4762a1bSJed Brown {
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 
24c4762a1bSJed Brown static void f1_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
25c4762a1bSJed Brown                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
26c4762a1bSJed Brown                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
27c4762a1bSJed Brown                     PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
28c4762a1bSJed Brown {
29c4762a1bSJed Brown   const PetscInt Ncomp = dim;
30c4762a1bSJed Brown   PetscInt       comp, d;
31c4762a1bSJed Brown   for (comp = 0; comp < Ncomp; ++comp) {
32c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
33c4762a1bSJed Brown       f1[comp*dim+d] = 0.0;
34c4762a1bSJed Brown     }
35c4762a1bSJed Brown   }
36c4762a1bSJed Brown }
37c4762a1bSJed Brown 
38c4762a1bSJed Brown /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */
39c4762a1bSJed Brown static void f1_u_3d_alpha(PetscInt dim, PetscInt Nf, PetscInt NfAux,
40c4762a1bSJed Brown                           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
41c4762a1bSJed Brown                           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
42c4762a1bSJed Brown                           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
43c4762a1bSJed Brown {
44c4762a1bSJed Brown   PetscReal trace,mu=s_mu,lambda=s_lambda,rad;
45c4762a1bSJed Brown   PetscInt i,j;
46c4762a1bSJed Brown   for (i=0,rad=0.;i<dim;i++) {
47c4762a1bSJed Brown     PetscReal t=x[i];
48c4762a1bSJed Brown     rad += t*t;
49c4762a1bSJed Brown   }
50c4762a1bSJed Brown   rad = PetscSqrtReal(rad);
51c4762a1bSJed Brown   if (rad>0.25) {
52c4762a1bSJed Brown     mu *= s_soft_alpha;
53c4762a1bSJed Brown     lambda *= s_soft_alpha; /* we could keep the bulk the same like rubberish */
54c4762a1bSJed Brown   }
55c4762a1bSJed Brown   for (i=0,trace=0; i < dim; ++i) {
56c4762a1bSJed Brown     trace += PetscRealPart(u_x[i*dim+i]);
57c4762a1bSJed Brown   }
58c4762a1bSJed Brown   for (i=0; i < dim; ++i) {
59c4762a1bSJed Brown     for (j=0; j < dim; ++j) {
60c4762a1bSJed Brown       f1[i*dim+j] = mu*(u_x[i*dim+j]+u_x[j*dim+i]);
61c4762a1bSJed Brown     }
62c4762a1bSJed Brown     f1[i*dim+i] += lambda*trace;
63c4762a1bSJed Brown   }
64c4762a1bSJed Brown }
65c4762a1bSJed Brown 
66c4762a1bSJed Brown /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */
67c4762a1bSJed Brown static void f1_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
68c4762a1bSJed Brown                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
69c4762a1bSJed Brown                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
70c4762a1bSJed Brown                     PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
71c4762a1bSJed Brown {
72c4762a1bSJed Brown   PetscReal trace,mu=s_mu,lambda=s_lambda;
73c4762a1bSJed Brown   PetscInt i,j;
74c4762a1bSJed Brown   for (i=0,trace=0; i < dim; ++i) {
75c4762a1bSJed Brown     trace += PetscRealPart(u_x[i*dim+i]);
76c4762a1bSJed Brown   }
77c4762a1bSJed Brown   for (i=0; i < dim; ++i) {
78c4762a1bSJed Brown     for (j=0; j < dim; ++j) {
79c4762a1bSJed Brown       f1[i*dim+j] = mu*(u_x[i*dim+j]+u_x[j*dim+i]);
80c4762a1bSJed Brown     }
81c4762a1bSJed Brown     f1[i*dim+i] += lambda*trace;
82c4762a1bSJed Brown   }
83c4762a1bSJed Brown }
84c4762a1bSJed Brown 
85c4762a1bSJed Brown /* 3D elasticity */
86c4762a1bSJed Brown #define IDX(ii,jj,kk,ll) (27*ii+9*jj+3*kk+ll)
87c4762a1bSJed Brown 
88c4762a1bSJed Brown void g3_uu_3d_private( PetscScalar g3[], const PetscReal mu, const PetscReal lambda)
89c4762a1bSJed Brown {
90c4762a1bSJed Brown   if (1) {
91c4762a1bSJed Brown     g3[0] += lambda;
92c4762a1bSJed Brown     g3[0] += mu;
93c4762a1bSJed Brown     g3[0] += mu;
94c4762a1bSJed Brown     g3[4] += lambda;
95c4762a1bSJed Brown     g3[8] += lambda;
96c4762a1bSJed Brown     g3[10] += mu;
97c4762a1bSJed Brown     g3[12] += mu;
98c4762a1bSJed Brown     g3[20] += mu;
99c4762a1bSJed Brown     g3[24] += mu;
100c4762a1bSJed Brown     g3[28] += mu;
101c4762a1bSJed Brown     g3[30] += mu;
102c4762a1bSJed Brown     g3[36] += lambda;
103c4762a1bSJed Brown     g3[40] += lambda;
104c4762a1bSJed Brown     g3[40] += mu;
105c4762a1bSJed Brown     g3[40] += mu;
106c4762a1bSJed Brown     g3[44] += lambda;
107c4762a1bSJed Brown     g3[50] += mu;
108c4762a1bSJed Brown     g3[52] += mu;
109c4762a1bSJed Brown     g3[56] += mu;
110c4762a1bSJed Brown     g3[60] += mu;
111c4762a1bSJed Brown     g3[68] += mu;
112c4762a1bSJed Brown     g3[70] += mu;
113c4762a1bSJed Brown     g3[72] += lambda;
114c4762a1bSJed Brown     g3[76] += lambda;
115c4762a1bSJed Brown     g3[80] += lambda;
116c4762a1bSJed Brown     g3[80] += mu;
117c4762a1bSJed Brown     g3[80] += mu;
118c4762a1bSJed Brown   } else {
119c4762a1bSJed Brown     int        i,j,k,l;
120c4762a1bSJed Brown     static int cc=-1;
121c4762a1bSJed Brown     cc++;
122c4762a1bSJed Brown     for (i = 0; i < 3; ++i) {
123c4762a1bSJed Brown       for (j = 0; j < 3; ++j) {
124c4762a1bSJed Brown         for (k = 0; k < 3; ++k) {
125c4762a1bSJed Brown           for (l = 0; l < 3; ++l) {
126c4762a1bSJed Brown             if (k==l && i==j) g3[IDX(i,j,k,l)] += lambda;
127c4762a1bSJed Brown             if (i==k && j==l) g3[IDX(i,j,k,l)] += mu;
128c4762a1bSJed Brown             if (i==l && j==k) g3[IDX(i,j,k,l)] += mu;
129c4762a1bSJed Brown             if (k==l && i==j && !cc) (void) PetscPrintf(PETSC_COMM_WORLD,"g3[%d] += lambda;\n",IDX(i,j,k,l));
130c4762a1bSJed Brown             if (i==k && j==l && !cc) (void) PetscPrintf(PETSC_COMM_WORLD,"g3[%d] += mu;\n",IDX(i,j,k,l));
131c4762a1bSJed Brown             if (i==l && j==k && !cc) (void) PetscPrintf(PETSC_COMM_WORLD,"g3[%d] += mu;\n",IDX(i,j,k,l));
132c4762a1bSJed Brown           }
133c4762a1bSJed Brown         }
134c4762a1bSJed Brown       }
135c4762a1bSJed Brown     }
136c4762a1bSJed Brown   }
137c4762a1bSJed Brown }
138c4762a1bSJed Brown 
139c4762a1bSJed Brown static void g3_uu_3d_alpha(PetscInt dim, PetscInt Nf, PetscInt NfAux,
140c4762a1bSJed Brown                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
141c4762a1bSJed Brown                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
142c4762a1bSJed Brown                            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
143c4762a1bSJed Brown {
144c4762a1bSJed Brown   PetscReal mu=s_mu, lambda=s_lambda,rad;
145c4762a1bSJed Brown   PetscInt i;
146c4762a1bSJed Brown   for (i=0,rad=0.;i<dim;i++) {
147c4762a1bSJed Brown     PetscReal t=x[i];
148c4762a1bSJed Brown     rad += t*t;
149c4762a1bSJed Brown   }
150c4762a1bSJed Brown   rad = PetscSqrtReal(rad);
151c4762a1bSJed Brown   if (rad>0.25) {
152c4762a1bSJed Brown     mu *= s_soft_alpha;
153c4762a1bSJed Brown     lambda *= s_soft_alpha; /* we could keep the bulk the same like rubberish */
154c4762a1bSJed Brown   }
155c4762a1bSJed Brown   g3_uu_3d_private(g3,mu,lambda);
156c4762a1bSJed Brown }
157c4762a1bSJed Brown 
158c4762a1bSJed Brown static void g3_uu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
159c4762a1bSJed Brown                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
160c4762a1bSJed Brown                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
161c4762a1bSJed Brown                      PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
162c4762a1bSJed Brown {
163c4762a1bSJed Brown   g3_uu_3d_private(g3,s_mu,s_lambda);
164c4762a1bSJed Brown }
165c4762a1bSJed Brown 
166c4762a1bSJed Brown static void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
167c4762a1bSJed Brown                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
168c4762a1bSJed Brown                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
169c4762a1bSJed Brown                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
170c4762a1bSJed Brown {
171c4762a1bSJed Brown   const    PetscInt Ncomp = dim;
172c4762a1bSJed Brown   PetscInt comp;
173c4762a1bSJed Brown 
174c4762a1bSJed Brown   for (comp = 0; comp < Ncomp; ++comp) f0[comp] = 0.0;
175c4762a1bSJed Brown }
176c4762a1bSJed Brown 
177c4762a1bSJed Brown /* PI_i (x_i^4 - x_i^2) */
178c4762a1bSJed Brown static void f0_u_x4(PetscInt dim, PetscInt Nf, PetscInt NfAux,
179c4762a1bSJed Brown                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
180c4762a1bSJed Brown                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
181c4762a1bSJed Brown                     PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
182c4762a1bSJed Brown {
183c4762a1bSJed Brown   const    PetscInt Ncomp = dim;
184c4762a1bSJed Brown   PetscInt comp,i;
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   for (comp = 0; comp < Ncomp; ++comp) {
187c4762a1bSJed Brown     f0[comp] = 1e5;
188c4762a1bSJed Brown     for (i = 0; i < Ncomp; ++i) {
189c4762a1bSJed Brown       f0[comp] *= /* (comp+1)* */(x[i]*x[i]*x[i]*x[i] - x[i]*x[i]); /* assumes (0,1]^D domain */
190c4762a1bSJed Brown     }
191c4762a1bSJed Brown   }
192c4762a1bSJed Brown }
193c4762a1bSJed Brown 
194c4762a1bSJed Brown PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
195c4762a1bSJed Brown {
196c4762a1bSJed Brown   const PetscInt Ncomp = dim;
197c4762a1bSJed Brown   PetscInt       comp;
198c4762a1bSJed Brown 
199c4762a1bSJed Brown   for (comp = 0; comp < Ncomp; ++comp) u[comp] = 0;
200c4762a1bSJed Brown   return 0;
201c4762a1bSJed Brown }
202c4762a1bSJed Brown 
203c4762a1bSJed Brown int main(int argc,char **args)
204c4762a1bSJed Brown {
205c4762a1bSJed Brown   Mat                Amat;
206c4762a1bSJed Brown   SNES               snes;
207c4762a1bSJed Brown   KSP                ksp;
208c4762a1bSJed Brown   MPI_Comm           comm;
209c4762a1bSJed Brown   PetscMPIInt        rank;
210956f8c0dSBarry Smith #if defined(PETSC_USE_LOG)
211c4762a1bSJed Brown   PetscLogStage      stage[17];
212956f8c0dSBarry Smith #endif
213*2a5e0493SMark Adams   PetscErrorCode     ierr;
214c4762a1bSJed Brown   PetscBool          test_nonzero_cols = PETSC_FALSE,use_nearnullspace = PETSC_TRUE,attach_nearnullspace = PETSC_FALSE;
215c4762a1bSJed Brown   Vec                xx,bb;
216c4762a1bSJed Brown   PetscInt           iter,i,N,dim = 3,cells[3] = {1,1,1},max_conv_its,local_sizes[7],run_type = 1;
217c4762a1bSJed Brown   DM                 dm,distdm,basedm;
218c4762a1bSJed Brown   PetscBool          flg;
219c4762a1bSJed Brown   char               convType[256];
220c4762a1bSJed Brown   PetscReal          Lx,mdisp[10],err[10];
221c4762a1bSJed Brown   const char * const options[10] = {"-ex56_dm_refine 0",
222c4762a1bSJed Brown                                     "-ex56_dm_refine 1",
223c4762a1bSJed Brown                                     "-ex56_dm_refine 2",
224c4762a1bSJed Brown                                     "-ex56_dm_refine 3",
225c4762a1bSJed Brown                                     "-ex56_dm_refine 4",
226c4762a1bSJed Brown                                     "-ex56_dm_refine 5",
227c4762a1bSJed Brown                                     "-ex56_dm_refine 6",
228c4762a1bSJed Brown                                     "-ex56_dm_refine 7",
229c4762a1bSJed Brown                                     "-ex56_dm_refine 8",
230c4762a1bSJed Brown                                     "-ex56_dm_refine 9"};
231c4762a1bSJed Brown   PetscFunctionBeginUser;
2329566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
233c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
2349566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
235c4762a1bSJed Brown   /* options */
2369566063dSJacob Faibussowitsch   ierr = PetscOptionsBegin(comm,NULL,"3D bilinear Q1 elasticity options","");PetscCall(ierr);
237c4762a1bSJed Brown   {
238c4762a1bSJed Brown     i = 3;
2399566063dSJacob Faibussowitsch     PetscCall(PetscOptionsIntArray("-cells", "Number of (flux tube) processor in each dimension", "ex56.c", cells, &i, NULL));
240c4762a1bSJed Brown 
241c4762a1bSJed Brown     Lx = 1.; /* or ne for rod */
242c4762a1bSJed Brown     max_conv_its = 3;
2439566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-max_conv_its","Number of iterations in convergence study","",max_conv_its,&max_conv_its,NULL));
244e00437b9SBarry Smith     PetscCheck(max_conv_its > 0 && max_conv_its < 7,PETSC_COMM_WORLD, PETSC_ERR_USER, "Bad number of iterations for convergence test (%D)",max_conv_its);
2459566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-lx","Length of domain","",Lx,&Lx,NULL));
2469566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-alpha","material coefficient inside circle","",s_soft_alpha,&s_soft_alpha,NULL));
2479566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-test_nonzero_cols","nonzero test","",test_nonzero_cols,&test_nonzero_cols,NULL));
2489566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-use_mat_nearnullspace","MatNearNullSpace API test","",use_nearnullspace,&use_nearnullspace,NULL));
2499566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-attach_mat_nearnullspace","MatNearNullSpace API test (via MatSetNearNullSpace)","",attach_nearnullspace,&attach_nearnullspace,NULL));
2509566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-run_type","0: twisting load on cantalever, 1: 3rd order accurate convergence test","",run_type,&run_type,NULL));
251c4762a1bSJed Brown   }
2529566063dSJacob Faibussowitsch   ierr = PetscOptionsEnd();PetscCall(ierr);
2539566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("Mesh Setup", &stage[16]));
254c4762a1bSJed Brown   for (iter=0 ; iter<max_conv_its ; iter++) {
255c4762a1bSJed Brown     char str[] = "Solve 0";
256c4762a1bSJed Brown     str[6] += iter;
2579566063dSJacob Faibussowitsch     PetscCall(PetscLogStageRegister(str, &stage[iter]));
258c4762a1bSJed Brown   }
259c4762a1bSJed Brown   /* create DM, Plex calls DMSetup */
2609566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(stage[16]));
2619566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, cells, NULL, NULL, NULL, PETSC_TRUE, &dm));
262c4762a1bSJed Brown   {
263c4762a1bSJed Brown     DMLabel         label;
264c4762a1bSJed Brown     IS              is;
2659566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(dm, "boundary"));
2669566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dm, "boundary", &label));
2679566063dSJacob Faibussowitsch     PetscCall(DMPlexMarkBoundaryFaces(dm, 1, label));
268*2a5e0493SMark Adams     if (run_type == 0) {
2699566063dSJacob Faibussowitsch       PetscCall(DMGetStratumIS(dm, "boundary", 1,  &is));
2709566063dSJacob Faibussowitsch       PetscCall(DMCreateLabel(dm,"Faces"));
271c4762a1bSJed Brown       if (is) {
272c4762a1bSJed Brown         PetscInt        d, f, Nf;
273c4762a1bSJed Brown         const PetscInt *faces;
274c4762a1bSJed Brown         PetscInt        csize;
275c4762a1bSJed Brown         PetscSection    cs;
276c4762a1bSJed Brown         Vec             coordinates ;
277c4762a1bSJed Brown         DM              cdm;
2789566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(is, &Nf));
2799566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(is, &faces));
2809566063dSJacob Faibussowitsch         PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
2819566063dSJacob Faibussowitsch         PetscCall(DMGetCoordinateDM(dm, &cdm));
2829566063dSJacob Faibussowitsch         PetscCall(DMGetLocalSection(cdm, &cs));
283c4762a1bSJed Brown         /* Check for each boundary face if any component of its centroid is either 0.0 or 1.0 */
284c4762a1bSJed Brown         for (f = 0; f < Nf; ++f) {
285c4762a1bSJed Brown           PetscReal   faceCoord;
286c4762a1bSJed Brown           PetscInt    b,v;
287c4762a1bSJed Brown           PetscScalar *coords = NULL;
288c4762a1bSJed Brown           PetscInt    Nv;
2899566063dSJacob Faibussowitsch           PetscCall(DMPlexVecGetClosure(cdm, cs, coordinates, faces[f], &csize, &coords));
290c4762a1bSJed Brown           Nv   = csize/dim; /* Calculate mean coordinate vector */
291c4762a1bSJed Brown           for (d = 0; d < dim; ++d) {
292c4762a1bSJed Brown             faceCoord = 0.0;
293c4762a1bSJed Brown             for (v = 0; v < Nv; ++v) faceCoord += PetscRealPart(coords[v*dim+d]);
294c4762a1bSJed Brown             faceCoord /= Nv;
295c4762a1bSJed Brown             for (b = 0; b < 2; ++b) {
296c4762a1bSJed Brown               if (PetscAbs(faceCoord - b) < PETSC_SMALL) { /* domain have not been set yet, still [0,1]^3 */
2979566063dSJacob Faibussowitsch                 PetscCall(DMSetLabelValue(dm, "Faces", faces[f], d*2+b+1));
298c4762a1bSJed Brown               }
299c4762a1bSJed Brown             }
300c4762a1bSJed Brown           }
3019566063dSJacob Faibussowitsch           PetscCall(DMPlexVecRestoreClosure(cdm, cs, coordinates, faces[f], &csize, &coords));
302c4762a1bSJed Brown         }
3039566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(is, &faces));
304c4762a1bSJed Brown       }
3059566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&is));
3069566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(dm, "Faces", &label));
3079566063dSJacob Faibussowitsch       PetscCall(DMPlexLabelComplete(dm, label));
308c4762a1bSJed Brown     }
309c4762a1bSJed Brown   }
310c4762a1bSJed Brown   {
311c4762a1bSJed Brown     PetscInt    dimEmbed, i;
312c4762a1bSJed Brown     PetscInt    nCoords;
313c4762a1bSJed Brown     PetscScalar *coords,bounds[] = {0,1,-.5,.5,-.5,.5,}; /* x_min,x_max,y_min,y_max */
314c4762a1bSJed Brown     Vec         coordinates;
315c4762a1bSJed Brown     bounds[1] = Lx;
316c4762a1bSJed Brown     if (run_type==1) {
317c4762a1bSJed Brown       for (i = 0; i < 2*dim; i++) bounds[i] = (i%2) ? 1 : 0;
318c4762a1bSJed Brown     }
3199566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dm,&coordinates));
3209566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDim(dm,&dimEmbed));
321e00437b9SBarry Smith     PetscCheck(dimEmbed == dim,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"dimEmbed != dim %D",dimEmbed);
3229566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(coordinates,&nCoords));
323e00437b9SBarry Smith     PetscCheck((nCoords % dimEmbed) == 0,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Coordinate vector the wrong size");
3249566063dSJacob Faibussowitsch     PetscCall(VecGetArray(coordinates,&coords));
325c4762a1bSJed Brown     for (i = 0; i < nCoords; i += dimEmbed) {
326c4762a1bSJed Brown       PetscInt    j;
327c4762a1bSJed Brown       PetscScalar *coord = &coords[i];
328c4762a1bSJed Brown       for (j = 0; j < dimEmbed; j++) {
329c4762a1bSJed Brown         coord[j] = bounds[2 * j] + coord[j] * (bounds[2 * j + 1] - bounds[2 * j]);
330c4762a1bSJed Brown       }
331c4762a1bSJed Brown     }
3329566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(coordinates,&coords));
3339566063dSJacob Faibussowitsch     PetscCall(DMSetCoordinatesLocal(dm,coordinates));
334c4762a1bSJed Brown   }
335c4762a1bSJed Brown 
336c4762a1bSJed Brown   /* convert to p4est, and distribute */
3379566063dSJacob Faibussowitsch   ierr = PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");PetscCall(ierr);
3389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-dm_type","Convert DMPlex to another format (should not be Plex!)","ex56.c",DMList,DMPLEX,convType,256,&flg));
3399566063dSJacob Faibussowitsch   ierr = PetscOptionsEnd();PetscCall(ierr);
340c4762a1bSJed Brown   if (flg) {
341c4762a1bSJed Brown     DM newdm;
3429566063dSJacob Faibussowitsch     PetscCall(DMConvert(dm,convType,&newdm));
343c4762a1bSJed Brown     if (newdm) {
344c4762a1bSJed Brown       const char *prefix;
345c4762a1bSJed Brown       PetscBool isForest;
3469566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm,&prefix));
3479566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)newdm,prefix));
3489566063dSJacob Faibussowitsch       PetscCall(DMIsForest(newdm,&isForest));
34928b400f6SJacob Faibussowitsch       PetscCheck(isForest,PETSC_COMM_WORLD, PETSC_ERR_USER, "Converted to non Forest?");
3509566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&dm));
351c4762a1bSJed Brown       dm   = newdm;
352c4762a1bSJed Brown     } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_USER, "Convert failed?");
353c4762a1bSJed Brown   } else {
354c4762a1bSJed Brown     PetscPartitioner part;
355c4762a1bSJed Brown     /* Plex Distribute mesh over processes */
3569566063dSJacob Faibussowitsch     PetscCall(DMPlexGetPartitioner(dm,&part));
3579566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerSetFromOptions(part));
3589566063dSJacob Faibussowitsch     PetscCall(DMPlexDistribute(dm, 0, NULL, &distdm));
359c4762a1bSJed Brown     if (distdm) {
360c4762a1bSJed Brown       const char *prefix;
3619566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm,&prefix));
3629566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)distdm,prefix));
3639566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&dm));
364c4762a1bSJed Brown       dm   = distdm;
365c4762a1bSJed Brown     }
366c4762a1bSJed Brown   }
3679566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
368c4762a1bSJed Brown   basedm = dm; dm = NULL;
369c4762a1bSJed Brown 
370c4762a1bSJed Brown   for (iter=0 ; iter<max_conv_its ; iter++) {
3719566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(stage[16]));
372c4762a1bSJed Brown     /* make new DM */
3739566063dSJacob Faibussowitsch     PetscCall(DMClone(basedm, &dm));
3749566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject) dm, "ex56_"));
3759566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName( (PetscObject)dm,"Mesh"));
376c4762a1bSJed Brown     if (max_conv_its > 1) {
3770e75e42fSMark       /* If max_conv_its == 1, then we are not doing a convergence study. */
3789566063dSJacob Faibussowitsch       PetscCall(PetscOptionsInsertString(NULL,options[iter]));
379c4762a1bSJed Brown     }
3809566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(dm)); /* refinement done here in Plex, p4est */
381c4762a1bSJed Brown     /* snes */
3829566063dSJacob Faibussowitsch     PetscCall(SNESCreate(comm, &snes));
3839566063dSJacob Faibussowitsch     PetscCall(SNESSetDM(snes, dm));
384c4762a1bSJed Brown     /* fem */
385c4762a1bSJed Brown     {
386c4762a1bSJed Brown       const PetscInt Ncomp = dim;
387c4762a1bSJed Brown       const PetscInt components[] = {0,1,2};
388c4762a1bSJed Brown       const PetscInt Nfid = 1, Npid = 1;
389c4762a1bSJed Brown       const PetscInt fid[] = {1}; /* The fixed faces (x=0) */
390c4762a1bSJed Brown       const PetscInt pid[] = {2}; /* The faces with loading (x=L_x) */
391c4762a1bSJed Brown       PetscFE        fe;
392c4762a1bSJed Brown       PetscDS        prob;
39345480ffeSMatthew G. Knepley       DMLabel        label;
394c4762a1bSJed Brown       DM             cdm = dm;
395c4762a1bSJed Brown 
3969566063dSJacob Faibussowitsch       PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, dim, PETSC_FALSE, NULL, PETSC_DECIDE, &fe)); /* elasticity */
3979566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject) fe, "deformation"));
398c4762a1bSJed Brown       /* FEM prob */
3999566063dSJacob Faibussowitsch       PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe));
4009566063dSJacob Faibussowitsch       PetscCall(DMCreateDS(dm));
4019566063dSJacob Faibussowitsch       PetscCall(DMGetDS(dm, &prob));
402c4762a1bSJed Brown       /* setup problem */
403c4762a1bSJed Brown       if (run_type==1) {
4049566063dSJacob Faibussowitsch         PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu_3d));
4059566063dSJacob Faibussowitsch         PetscCall(PetscDSSetResidual(prob, 0, f0_u_x4, f1_u_3d));
406c4762a1bSJed Brown       } else {
40745480ffeSMatthew G. Knepley         PetscWeakForm wf;
40845480ffeSMatthew G. Knepley         PetscInt      bd, i;
40945480ffeSMatthew G. Knepley 
4109566063dSJacob Faibussowitsch         PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu_3d_alpha));
4119566063dSJacob Faibussowitsch         PetscCall(PetscDSSetResidual(prob, 0, f0_u, f1_u_3d_alpha));
41245480ffeSMatthew G. Knepley 
4139566063dSJacob Faibussowitsch         PetscCall(DMGetLabel(dm, "Faces", &label));
4149566063dSJacob Faibussowitsch         PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "traction", label, Npid, pid, 0, Ncomp, components, NULL, NULL, NULL, &bd));
4159566063dSJacob Faibussowitsch         PetscCall(PetscDSGetBoundary(prob, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
4169566063dSJacob Faibussowitsch         for (i = 0; i < Npid; ++i) PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, pid[i], 0, 0, 0, f0_bd_u_3d, 0, f1_bd_u));
417c4762a1bSJed Brown       }
418c4762a1bSJed Brown       /* bcs */
419c4762a1bSJed Brown       if (run_type==1) {
420c4762a1bSJed Brown         PetscInt id = 1;
4219566063dSJacob Faibussowitsch         PetscCall(DMGetLabel(dm, "boundary", &label));
4229566063dSJacob Faibussowitsch         PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) zero, NULL, NULL, NULL));
423c4762a1bSJed Brown       } else {
4249566063dSJacob Faibussowitsch         PetscCall(DMGetLabel(dm, "Faces", &label));
4259566063dSJacob Faibussowitsch         PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "fixed", label, Nfid, fid, 0, Ncomp, components, (void (*)(void)) zero, NULL, NULL, NULL));
426c4762a1bSJed Brown       }
427c4762a1bSJed Brown       while (cdm) {
4289566063dSJacob Faibussowitsch         PetscCall(DMCopyDisc(dm, cdm));
4299566063dSJacob Faibussowitsch         PetscCall(DMGetCoarseDM(cdm, &cdm));
430c4762a1bSJed Brown       }
4319566063dSJacob Faibussowitsch       PetscCall(PetscFEDestroy(&fe));
432c4762a1bSJed Brown     }
433c4762a1bSJed Brown     /* vecs & mat */
4349566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(dm,&xx));
4359566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx, &bb));
4369566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) bb, "b"));
4379566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) xx, "u"));
4389566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix(dm, &Amat));
4399566063dSJacob Faibussowitsch     PetscCall(MatSetOption(Amat,MAT_SYMMETRIC,PETSC_TRUE));        /* Some matrix kernels can take advantage of symmetry if we set this. */
4409566063dSJacob Faibussowitsch     PetscCall(MatSetOption(Amat,MAT_SYMMETRY_ETERNAL,PETSC_TRUE)); /* Inform PETSc that Amat is always symmetric, so info set above isn't lost. */
4419566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(Amat,3));
4429566063dSJacob Faibussowitsch     PetscCall(MatSetOption(Amat,MAT_SPD,PETSC_TRUE));
4439566063dSJacob Faibussowitsch     PetscCall(VecGetSize(bb,&N));
444c4762a1bSJed Brown     local_sizes[iter] = N;
4459566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes,"%D global equations, %D vertices\n",N,N/dim));
446c4762a1bSJed Brown     if ((use_nearnullspace || attach_nearnullspace) && N/dim > 1) {
447c4762a1bSJed Brown       /* Set up the near null space (a.k.a. rigid body modes) that will be used by the multigrid preconditioner */
448c4762a1bSJed Brown       DM           subdm;
449c4762a1bSJed Brown       MatNullSpace nearNullSpace;
450c4762a1bSJed Brown       PetscInt     fields = 0;
451c4762a1bSJed Brown       PetscObject  deformation;
4529566063dSJacob Faibussowitsch       PetscCall(DMCreateSubDM(dm, 1, &fields, NULL, &subdm));
4539566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateRigidBody(subdm, 0, &nearNullSpace));
4549566063dSJacob Faibussowitsch       PetscCall(DMGetField(dm, 0, NULL, &deformation));
4559566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose(deformation, "nearnullspace", (PetscObject) nearNullSpace));
4569566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&subdm));
457c4762a1bSJed Brown       if (attach_nearnullspace) {
4589566063dSJacob Faibussowitsch         PetscCall(MatSetNearNullSpace(Amat,nearNullSpace));
459c4762a1bSJed Brown       }
4609566063dSJacob Faibussowitsch       PetscCall(MatNullSpaceDestroy(&nearNullSpace)); /* created by DM and destroyed by Mat */
461c4762a1bSJed Brown     }
4629566063dSJacob Faibussowitsch     PetscCall(DMPlexSetSNESLocalFEM(dm,NULL,NULL,NULL));
4639566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, Amat, Amat, NULL, NULL));
4649566063dSJacob Faibussowitsch     PetscCall(SNESSetFromOptions(snes));
4659566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dm));
4669566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
4679566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(stage[16]));
468c4762a1bSJed Brown     /* ksp */
4699566063dSJacob Faibussowitsch     PetscCall(SNESGetKSP(snes, &ksp));
4709566063dSJacob Faibussowitsch     PetscCall(KSPSetComputeSingularValues(ksp,PETSC_TRUE));
471c4762a1bSJed Brown     /* test BCs */
4729566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(xx));
473c4762a1bSJed Brown     if (test_nonzero_cols) {
474dd400576SPatrick Sanan       if (rank == 0) {
4759566063dSJacob Faibussowitsch         PetscCall(VecSetValue(xx,0,1.0,INSERT_VALUES));
476c4762a1bSJed Brown       }
4779566063dSJacob Faibussowitsch       PetscCall(VecAssemblyBegin(xx));
4789566063dSJacob Faibussowitsch       PetscCall(VecAssemblyEnd(xx));
479c4762a1bSJed Brown     }
4809566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(bb));
4819566063dSJacob Faibussowitsch     PetscCall(VecGetSize(bb,&i));
482c4762a1bSJed Brown     local_sizes[iter] = i;
4839566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes,"%D equations in vector, %D vertices\n",i,i/dim));
4849566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
485c4762a1bSJed Brown     /* solve */
4869566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(stage[iter]));
4879566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes, bb, xx));
4889566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
4899566063dSJacob Faibussowitsch     PetscCall(VecNorm(xx,NORM_INFINITY,&mdisp[iter]));
4909566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
491c4762a1bSJed Brown     {
492c4762a1bSJed Brown       PetscViewer       viewer = NULL;
493c4762a1bSJed Brown       PetscViewerFormat fmt;
4949566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetViewer(comm,NULL,"ex56_","-vec_view",&viewer,&fmt,&flg));
495c4762a1bSJed Brown       if (flg) {
4969566063dSJacob Faibussowitsch         PetscCall(PetscViewerPushFormat(viewer,fmt));
4979566063dSJacob Faibussowitsch         PetscCall(VecView(xx,viewer));
4989566063dSJacob Faibussowitsch         PetscCall(VecView(bb,viewer));
4999566063dSJacob Faibussowitsch         PetscCall(PetscViewerPopFormat(viewer));
500c4762a1bSJed Brown       }
5019566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
502c4762a1bSJed Brown     }
503c4762a1bSJed Brown     /* Free work space */
5049566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dm));
5059566063dSJacob Faibussowitsch     PetscCall(SNESDestroy(&snes));
5069566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xx));
5079566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bb));
5089566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Amat));
509c4762a1bSJed Brown   }
5109566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&basedm));
511c4762a1bSJed Brown   if (run_type==1) err[0] = 59.975208 - mdisp[0]; /* error with what I think is the exact solution */
512c4762a1bSJed Brown   else             err[0] = 171.038 - mdisp[0];
513c4762a1bSJed Brown   for (iter=1 ; iter<max_conv_its ; iter++) {
514c4762a1bSJed Brown     if (run_type==1) err[iter] = 59.975208 - mdisp[iter];
515c4762a1bSJed Brown     else             err[iter] = 171.038 - mdisp[iter];
516*2a5e0493SMark Adams     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"[%d] %D) N=%12D, max displ=%9.7e, disp diff=%9.2e, error=%4.3e, rate=%3.2g\n",rank,iter,local_sizes[iter],(double)mdisp[iter],
517*2a5e0493SMark Adams                           (double)(mdisp[iter]-mdisp[iter-1]),(double)err[iter],(double)(PetscLogReal(err[iter-1]/err[iter])/PetscLogReal(2.))));
518c4762a1bSJed Brown   }
519c4762a1bSJed Brown 
5209566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
521b122ec5aSJacob Faibussowitsch   return 0;
522c4762a1bSJed Brown }
523c4762a1bSJed Brown 
524c4762a1bSJed Brown /*TEST
525c4762a1bSJed Brown 
526c4762a1bSJed Brown   test:
527c4762a1bSJed Brown     suffix: 0
528c4762a1bSJed Brown     nsize: 4
529c4762a1bSJed Brown     requires: !single
530*2a5e0493SMark Adams     args: -cells 2,2,1 -max_conv_its 2 -petscspace_degree 3 -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_square_graph 0 -pc_gamg_threshold 0.001 -ksp_converged_reason -snes_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 -ex56_dm_view -snes_lag_jacobian -2 -snes_type ksponly -use_gpu_aware_mpi true
531c4762a1bSJed Brown     timeoutfactor: 2
532c4762a1bSJed Brown 
533c4762a1bSJed Brown   # HYPRE PtAP broken with complex numbers
534c4762a1bSJed Brown   test:
535c4762a1bSJed Brown     suffix: hypre
536263f2b91SStefano Zampini     requires: hypre !single !complex !defined(PETSC_HAVE_HYPRE_DEVICE)
537c4762a1bSJed Brown     nsize: 4
5380e75e42fSMark     args: -cells 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
539c4762a1bSJed Brown 
540c4762a1bSJed Brown   test:
541c4762a1bSJed Brown     suffix: ml
542c4762a1bSJed Brown     requires: ml !single
543c4762a1bSJed Brown     nsize: 4
5440e75e42fSMark     args: -cells 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
545c4762a1bSJed Brown 
546c4762a1bSJed Brown   test:
547c4762a1bSJed Brown     suffix: hpddm
548dfd57a17SPierre Jolivet     requires: hpddm slepc !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
549c4762a1bSJed Brown     nsize: 4
5500e75e42fSMark     args: -cells 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
551c4762a1bSJed Brown 
552c4762a1bSJed Brown   test:
55363b77682SMark Adams     suffix: repart
554c4762a1bSJed Brown     nsize: 4
555c4762a1bSJed Brown     requires: parmetis !single
55673f7197eSJed Brown     args: -cells 8,2,2 -max_conv_its 1 -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_square_graph 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 -snes_converged_reason -pc_gamg_process_eq_limit 20 -pc_gamg_coarse_eq_limit 10 -ksp_converged_reason -snes_converged_reason -pc_gamg_reuse_interpolation true
557c4762a1bSJed Brown 
558c4762a1bSJed Brown   test:
559c4762a1bSJed Brown     suffix: bddc
560c4762a1bSJed Brown     nsize: 4
561c4762a1bSJed Brown     requires: !single
5620e75e42fSMark     args: -cells 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 -ex56_dm_mat_type is -matis_localmat_type {{sbaij baij aij}} -pc_type bddc
563c4762a1bSJed Brown 
564c4762a1bSJed Brown   testset:
565c4762a1bSJed Brown     nsize: 4
566c4762a1bSJed Brown     requires: !single
5670e75e42fSMark     args: -cells 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 -ex56_dm_mat_type is -matis_localmat_type aij -pc_type bddc -attach_mat_nearnullspace {{0 1}separate output}
568c4762a1bSJed Brown     test:
569c4762a1bSJed Brown       suffix: bddc_approx_gamg
57073f7197eSJed Brown       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_square_graph 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_square_graph 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -prefix_pop
571c4762a1bSJed Brown     # HYPRE PtAP broken with complex numbers
572c4762a1bSJed Brown     test:
573263f2b91SStefano Zampini       requires: hypre !complex !defined(PETSC_HAVE_HYPRE_DEVICE)
574c4762a1bSJed Brown       suffix: bddc_approx_hypre
575c4762a1bSJed 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
576c4762a1bSJed Brown     test:
577c4762a1bSJed Brown       requires: ml
578c4762a1bSJed Brown       suffix: bddc_approx_ml
5790e75e42fSMark       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
580c4762a1bSJed Brown 
581c4762a1bSJed Brown   test:
582c4762a1bSJed Brown     suffix: fetidp
583c4762a1bSJed Brown     nsize: 4
584c4762a1bSJed Brown     requires: !single
5850e75e42fSMark     args: -cells 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 -ex56_dm_mat_type is -matis_localmat_type {{sbaij baij aij}}
586c4762a1bSJed Brown 
587c4762a1bSJed Brown   test:
588c4762a1bSJed Brown     suffix: bddc_elast
589c4762a1bSJed Brown     nsize: 4
590c4762a1bSJed Brown     requires: !single
5910e75e42fSMark     args: -cells 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 -ex56_dm_mat_type is -matis_localmat_type sbaij -pc_type bddc -pc_bddc_monolithic -attach_mat_nearnullspace
592c4762a1bSJed Brown 
593c4762a1bSJed Brown   test:
594c4762a1bSJed Brown     suffix: fetidp_elast
595c4762a1bSJed Brown     nsize: 4
596c4762a1bSJed Brown     requires: !single
5970e75e42fSMark     args: -cells 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 -ex56_dm_mat_type is -matis_localmat_type sbaij -fetidp_bddc_pc_bddc_monolithic -attach_mat_nearnullspace
598c4762a1bSJed Brown 
59935990778SJunchao Zhang   testset:
60035990778SJunchao Zhang     nsize: 4
60135990778SJunchao Zhang     requires: !single
60273f7197eSJed Brown     args: -cells 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_square_graph 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 -ksp_monitor_short -ksp_converged_reason -snes_converged_reason -snes_monitor_short -ex56_dm_view -petscpartitioner_type simple -pc_gamg_process_eq_limit 20
60335990778SJunchao Zhang     output_file: output/ex56_cuda.out
60435990778SJunchao Zhang 
605c4762a1bSJed Brown     test:
606c4762a1bSJed Brown       suffix: cuda
60735990778SJunchao Zhang       requires: cuda
60835990778SJunchao Zhang       args: -ex56_dm_mat_type aijcusparse -ex56_dm_vec_type cuda
6096cb74228SMark Adams 
6106cb74228SMark Adams     test:
6116cb74228SMark Adams       suffix: viennacl
61235990778SJunchao Zhang       requires: viennacl
61335990778SJunchao Zhang       args: -ex56_dm_mat_type aijviennacl -ex56_dm_vec_type viennacl
6146cb74228SMark Adams 
61535990778SJunchao Zhang     test:
61635990778SJunchao Zhang       suffix: kokkos
6173078479eSJunchao Zhang       requires: !sycl kokkos_kernels
61835990778SJunchao Zhang       args: -ex56_dm_mat_type aijkokkos -ex56_dm_vec_type kokkos
619dea3b165SRichard Tran Mills   # Don't run AIJMKL caes with complex scalars because of convergence issues.
620dea3b165SRichard 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.
621a8736bf8SRichard Tran Mills   test:
622a8736bf8SRichard Tran Mills     suffix: seqaijmkl
623a8736bf8SRichard Tran Mills     nsize: 1
624dfd57a17SPierre Jolivet     requires: defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE) !single !complex
625a8736bf8SRichard Tran Mills     args: -cells 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_square_graph 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -ksp_converged_reason -snes_monitor_short -ksp_monitor_short -snes_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 -ex56_dm_view -run_type 1 -mat_seqaij_type seqaijmkl
626a8736bf8SRichard Tran Mills     timeoutfactor: 2
627a8736bf8SRichard Tran Mills 
628dea3b165SRichard Tran Mills   test:
629dea3b165SRichard Tran Mills     suffix: mpiaijmkl
630dea3b165SRichard Tran Mills     nsize: 2
631dfd57a17SPierre Jolivet     requires: defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE) !single !complex
632dea3b165SRichard Tran Mills     args: -cells 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_square_graph 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -ksp_converged_reason -snes_monitor_short -ksp_monitor_short -snes_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 -ex56_dm_view -run_type 1 -mat_seqaij_type seqaijmkl
633dea3b165SRichard Tran Mills     timeoutfactor: 2
634dea3b165SRichard Tran Mills 
635c4762a1bSJed Brown TEST*/
636