xref: /petsc/src/snes/tutorials/ex56.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "3D, tri-quadratic hexahedra (Q1), displacement finite element formulation\n\
2*c4762a1bSJed Brown of linear elasticity.  E=1.0, nu=1/3.\n\
3*c4762a1bSJed Brown Unit cube domain with Dirichlet boundary\n\n";
4*c4762a1bSJed Brown 
5*c4762a1bSJed Brown #include <petscdmplex.h>
6*c4762a1bSJed Brown #include <petscsnes.h>
7*c4762a1bSJed Brown #include <petscds.h>
8*c4762a1bSJed Brown #include <petscdmforest.h>
9*c4762a1bSJed Brown 
10*c4762a1bSJed Brown static PetscReal s_soft_alpha=1.e-3;
11*c4762a1bSJed Brown static PetscReal s_mu=0.4;
12*c4762a1bSJed Brown static PetscReal s_lambda=0.4;
13*c4762a1bSJed Brown 
14*c4762a1bSJed Brown static void f0_bd_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
15*c4762a1bSJed Brown                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
16*c4762a1bSJed Brown                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
17*c4762a1bSJed Brown                        PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
18*c4762a1bSJed Brown {
19*c4762a1bSJed Brown   f0[0] = 1;     /* x direction pull */
20*c4762a1bSJed Brown   f0[1] = -x[2]; /* add a twist around x-axis */
21*c4762a1bSJed Brown   f0[2] =  x[1];
22*c4762a1bSJed Brown }
23*c4762a1bSJed Brown 
24*c4762a1bSJed Brown static void f1_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
25*c4762a1bSJed Brown                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
26*c4762a1bSJed Brown                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
27*c4762a1bSJed Brown                     PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
28*c4762a1bSJed Brown {
29*c4762a1bSJed Brown   const PetscInt Ncomp = dim;
30*c4762a1bSJed Brown   PetscInt       comp, d;
31*c4762a1bSJed Brown   for (comp = 0; comp < Ncomp; ++comp) {
32*c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
33*c4762a1bSJed Brown       f1[comp*dim+d] = 0.0;
34*c4762a1bSJed Brown     }
35*c4762a1bSJed Brown   }
36*c4762a1bSJed Brown }
37*c4762a1bSJed Brown 
38*c4762a1bSJed Brown /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */
39*c4762a1bSJed Brown static void f1_u_3d_alpha(PetscInt dim, PetscInt Nf, PetscInt NfAux,
40*c4762a1bSJed Brown                           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
41*c4762a1bSJed Brown                           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
42*c4762a1bSJed Brown                           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
43*c4762a1bSJed Brown {
44*c4762a1bSJed Brown   PetscReal trace,mu=s_mu,lambda=s_lambda,rad;
45*c4762a1bSJed Brown   PetscInt i,j;
46*c4762a1bSJed Brown   for (i=0,rad=0.;i<dim;i++) {
47*c4762a1bSJed Brown     PetscReal t=x[i];
48*c4762a1bSJed Brown     rad += t*t;
49*c4762a1bSJed Brown   }
50*c4762a1bSJed Brown   rad = PetscSqrtReal(rad);
51*c4762a1bSJed Brown   if (rad>0.25) {
52*c4762a1bSJed Brown     mu *= s_soft_alpha;
53*c4762a1bSJed Brown     lambda *= s_soft_alpha; /* we could keep the bulk the same like rubberish */
54*c4762a1bSJed Brown   }
55*c4762a1bSJed Brown   for (i=0,trace=0; i < dim; ++i) {
56*c4762a1bSJed Brown     trace += PetscRealPart(u_x[i*dim+i]);
57*c4762a1bSJed Brown   }
58*c4762a1bSJed Brown   for (i=0; i < dim; ++i) {
59*c4762a1bSJed Brown     for (j=0; j < dim; ++j) {
60*c4762a1bSJed Brown       f1[i*dim+j] = mu*(u_x[i*dim+j]+u_x[j*dim+i]);
61*c4762a1bSJed Brown     }
62*c4762a1bSJed Brown     f1[i*dim+i] += lambda*trace;
63*c4762a1bSJed Brown   }
64*c4762a1bSJed Brown }
65*c4762a1bSJed Brown 
66*c4762a1bSJed Brown /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */
67*c4762a1bSJed Brown static void f1_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
68*c4762a1bSJed Brown                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
69*c4762a1bSJed Brown                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
70*c4762a1bSJed Brown                     PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
71*c4762a1bSJed Brown {
72*c4762a1bSJed Brown   PetscReal trace,mu=s_mu,lambda=s_lambda;
73*c4762a1bSJed Brown   PetscInt i,j;
74*c4762a1bSJed Brown   for (i=0,trace=0; i < dim; ++i) {
75*c4762a1bSJed Brown     trace += PetscRealPart(u_x[i*dim+i]);
76*c4762a1bSJed Brown   }
77*c4762a1bSJed Brown   for (i=0; i < dim; ++i) {
78*c4762a1bSJed Brown     for (j=0; j < dim; ++j) {
79*c4762a1bSJed Brown       f1[i*dim+j] = mu*(u_x[i*dim+j]+u_x[j*dim+i]);
80*c4762a1bSJed Brown     }
81*c4762a1bSJed Brown     f1[i*dim+i] += lambda*trace;
82*c4762a1bSJed Brown   }
83*c4762a1bSJed Brown }
84*c4762a1bSJed Brown 
85*c4762a1bSJed Brown /* 3D elasticity */
86*c4762a1bSJed Brown #define IDX(ii,jj,kk,ll) (27*ii+9*jj+3*kk+ll)
87*c4762a1bSJed Brown 
88*c4762a1bSJed Brown void g3_uu_3d_private( PetscScalar g3[], const PetscReal mu, const PetscReal lambda)
89*c4762a1bSJed Brown {
90*c4762a1bSJed Brown   if (1) {
91*c4762a1bSJed Brown     g3[0] += lambda;
92*c4762a1bSJed Brown     g3[0] += mu;
93*c4762a1bSJed Brown     g3[0] += mu;
94*c4762a1bSJed Brown     g3[4] += lambda;
95*c4762a1bSJed Brown     g3[8] += lambda;
96*c4762a1bSJed Brown     g3[10] += mu;
97*c4762a1bSJed Brown     g3[12] += mu;
98*c4762a1bSJed Brown     g3[20] += mu;
99*c4762a1bSJed Brown     g3[24] += mu;
100*c4762a1bSJed Brown     g3[28] += mu;
101*c4762a1bSJed Brown     g3[30] += mu;
102*c4762a1bSJed Brown     g3[36] += lambda;
103*c4762a1bSJed Brown     g3[40] += lambda;
104*c4762a1bSJed Brown     g3[40] += mu;
105*c4762a1bSJed Brown     g3[40] += mu;
106*c4762a1bSJed Brown     g3[44] += lambda;
107*c4762a1bSJed Brown     g3[50] += mu;
108*c4762a1bSJed Brown     g3[52] += mu;
109*c4762a1bSJed Brown     g3[56] += mu;
110*c4762a1bSJed Brown     g3[60] += mu;
111*c4762a1bSJed Brown     g3[68] += mu;
112*c4762a1bSJed Brown     g3[70] += mu;
113*c4762a1bSJed Brown     g3[72] += lambda;
114*c4762a1bSJed Brown     g3[76] += lambda;
115*c4762a1bSJed Brown     g3[80] += lambda;
116*c4762a1bSJed Brown     g3[80] += mu;
117*c4762a1bSJed Brown     g3[80] += mu;
118*c4762a1bSJed Brown   } else {
119*c4762a1bSJed Brown     int        i,j,k,l;
120*c4762a1bSJed Brown     static int cc=-1;
121*c4762a1bSJed Brown     cc++;
122*c4762a1bSJed Brown     for (i = 0; i < 3; ++i) {
123*c4762a1bSJed Brown       for (j = 0; j < 3; ++j) {
124*c4762a1bSJed Brown         for (k = 0; k < 3; ++k) {
125*c4762a1bSJed Brown           for (l = 0; l < 3; ++l) {
126*c4762a1bSJed Brown             if (k==l && i==j) g3[IDX(i,j,k,l)] += lambda;
127*c4762a1bSJed Brown             if (i==k && j==l) g3[IDX(i,j,k,l)] += mu;
128*c4762a1bSJed Brown             if (i==l && j==k) g3[IDX(i,j,k,l)] += mu;
129*c4762a1bSJed Brown 	    if (k==l && i==j && !cc) (void) PetscPrintf(PETSC_COMM_WORLD,"g3[%d] += lambda;\n",IDX(i,j,k,l));
130*c4762a1bSJed Brown 	    if (i==k && j==l && !cc) (void) PetscPrintf(PETSC_COMM_WORLD,"g3[%d] += mu;\n",IDX(i,j,k,l));
131*c4762a1bSJed Brown 	    if (i==l && j==k && !cc) (void) PetscPrintf(PETSC_COMM_WORLD,"g3[%d] += mu;\n",IDX(i,j,k,l));
132*c4762a1bSJed Brown           }
133*c4762a1bSJed Brown         }
134*c4762a1bSJed Brown       }
135*c4762a1bSJed Brown     }
136*c4762a1bSJed Brown   }
137*c4762a1bSJed Brown }
138*c4762a1bSJed Brown 
139*c4762a1bSJed Brown static void g3_uu_3d_alpha(PetscInt dim, PetscInt Nf, PetscInt NfAux,
140*c4762a1bSJed Brown                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
141*c4762a1bSJed Brown                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
142*c4762a1bSJed Brown                            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
143*c4762a1bSJed Brown {
144*c4762a1bSJed Brown   PetscReal mu=s_mu, lambda=s_lambda,rad;
145*c4762a1bSJed Brown   PetscInt i;
146*c4762a1bSJed Brown   for (i=0,rad=0.;i<dim;i++) {
147*c4762a1bSJed Brown     PetscReal t=x[i];
148*c4762a1bSJed Brown     rad += t*t;
149*c4762a1bSJed Brown   }
150*c4762a1bSJed Brown   rad = PetscSqrtReal(rad);
151*c4762a1bSJed Brown   if (rad>0.25) {
152*c4762a1bSJed Brown     mu *= s_soft_alpha;
153*c4762a1bSJed Brown     lambda *= s_soft_alpha; /* we could keep the bulk the same like rubberish */
154*c4762a1bSJed Brown   }
155*c4762a1bSJed Brown   g3_uu_3d_private(g3,mu,lambda);
156*c4762a1bSJed Brown }
157*c4762a1bSJed Brown 
158*c4762a1bSJed Brown static void g3_uu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
159*c4762a1bSJed Brown                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
160*c4762a1bSJed Brown                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
161*c4762a1bSJed Brown                      PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
162*c4762a1bSJed Brown {
163*c4762a1bSJed Brown   g3_uu_3d_private(g3,s_mu,s_lambda);
164*c4762a1bSJed Brown }
165*c4762a1bSJed Brown 
166*c4762a1bSJed Brown static void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
167*c4762a1bSJed Brown                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
168*c4762a1bSJed Brown                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
169*c4762a1bSJed Brown                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
170*c4762a1bSJed Brown {
171*c4762a1bSJed Brown   const    PetscInt Ncomp = dim;
172*c4762a1bSJed Brown   PetscInt comp;
173*c4762a1bSJed Brown 
174*c4762a1bSJed Brown   for (comp = 0; comp < Ncomp; ++comp) f0[comp] = 0.0;
175*c4762a1bSJed Brown }
176*c4762a1bSJed Brown 
177*c4762a1bSJed Brown /* PI_i (x_i^4 - x_i^2) */
178*c4762a1bSJed Brown static void f0_u_x4(PetscInt dim, PetscInt Nf, PetscInt NfAux,
179*c4762a1bSJed Brown                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
180*c4762a1bSJed Brown                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
181*c4762a1bSJed Brown                     PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
182*c4762a1bSJed Brown {
183*c4762a1bSJed Brown   const    PetscInt Ncomp = dim;
184*c4762a1bSJed Brown   PetscInt comp,i;
185*c4762a1bSJed Brown 
186*c4762a1bSJed Brown   for (comp = 0; comp < Ncomp; ++comp) {
187*c4762a1bSJed Brown     f0[comp] = 1e5;
188*c4762a1bSJed Brown     for (i = 0; i < Ncomp; ++i) {
189*c4762a1bSJed Brown       f0[comp] *= /* (comp+1)* */(x[i]*x[i]*x[i]*x[i] - x[i]*x[i]); /* assumes (0,1]^D domain */
190*c4762a1bSJed Brown     }
191*c4762a1bSJed Brown   }
192*c4762a1bSJed Brown }
193*c4762a1bSJed Brown 
194*c4762a1bSJed Brown PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
195*c4762a1bSJed Brown {
196*c4762a1bSJed Brown   const PetscInt Ncomp = dim;
197*c4762a1bSJed Brown   PetscInt       comp;
198*c4762a1bSJed Brown 
199*c4762a1bSJed Brown   for (comp = 0; comp < Ncomp; ++comp) u[comp] = 0;
200*c4762a1bSJed Brown   return 0;
201*c4762a1bSJed Brown }
202*c4762a1bSJed Brown 
203*c4762a1bSJed Brown int main(int argc,char **args)
204*c4762a1bSJed Brown {
205*c4762a1bSJed Brown   Mat            Amat;
206*c4762a1bSJed Brown   PetscErrorCode ierr;
207*c4762a1bSJed Brown   SNES           snes;
208*c4762a1bSJed Brown   KSP            ksp;
209*c4762a1bSJed Brown   MPI_Comm       comm;
210*c4762a1bSJed Brown   PetscMPIInt    rank;
211*c4762a1bSJed Brown   PetscLogStage  stage[17];
212*c4762a1bSJed Brown   PetscBool      test_nonzero_cols=PETSC_FALSE,use_nearnullspace=PETSC_TRUE,attach_nearnullspace=PETSC_FALSE;
213*c4762a1bSJed Brown   Vec            xx,bb;
214*c4762a1bSJed Brown   PetscInt       iter,i,N,dim=3,cells[3]={1,1,1},max_conv_its,local_sizes[7],run_type=1;
215*c4762a1bSJed Brown   DM             dm,distdm,basedm;
216*c4762a1bSJed Brown   PetscBool      flg;
217*c4762a1bSJed Brown   char           convType[256];
218*c4762a1bSJed Brown   PetscReal      Lx,mdisp[10],err[10];
219*c4762a1bSJed Brown   const char * const options[10] = {"-ex56_dm_refine 0",
220*c4762a1bSJed Brown                                     "-ex56_dm_refine 1",
221*c4762a1bSJed Brown                                     "-ex56_dm_refine 2",
222*c4762a1bSJed Brown                                     "-ex56_dm_refine 3",
223*c4762a1bSJed Brown                                     "-ex56_dm_refine 4",
224*c4762a1bSJed Brown                                     "-ex56_dm_refine 5",
225*c4762a1bSJed Brown                                     "-ex56_dm_refine 6",
226*c4762a1bSJed Brown                                     "-ex56_dm_refine 7",
227*c4762a1bSJed Brown                                     "-ex56_dm_refine 8",
228*c4762a1bSJed Brown                                     "-ex56_dm_refine 9"};
229*c4762a1bSJed Brown   PetscFunctionBeginUser;
230*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
231*c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
232*c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
233*c4762a1bSJed Brown   /* options */
234*c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm,NULL,"3D bilinear Q1 elasticity options","");CHKERRQ(ierr);
235*c4762a1bSJed Brown   {
236*c4762a1bSJed Brown     i = 3;
237*c4762a1bSJed Brown     ierr = PetscOptionsIntArray("-cells", "Number of (flux tube) processor in each dimension", "ex56.c", cells, &i, NULL);CHKERRQ(ierr);
238*c4762a1bSJed Brown 
239*c4762a1bSJed Brown     Lx = 1.; /* or ne for rod */
240*c4762a1bSJed Brown     max_conv_its = 3;
241*c4762a1bSJed Brown     ierr = PetscOptionsInt("-max_conv_its","Number of iterations in convergence study","",max_conv_its,&max_conv_its,NULL);CHKERRQ(ierr);
242*c4762a1bSJed Brown     if (max_conv_its<=0 || max_conv_its>7) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_USER, "Bad number of iterations for convergence test (%D)",max_conv_its);
243*c4762a1bSJed Brown     ierr = PetscOptionsReal("-lx","Length of domain","",Lx,&Lx,NULL);CHKERRQ(ierr);
244*c4762a1bSJed Brown     ierr = PetscOptionsReal("-alpha","material coefficient inside circle","",s_soft_alpha,&s_soft_alpha,NULL);CHKERRQ(ierr);
245*c4762a1bSJed Brown     ierr = PetscOptionsBool("-test_nonzero_cols","nonzero test","",test_nonzero_cols,&test_nonzero_cols,NULL);CHKERRQ(ierr);
246*c4762a1bSJed Brown     ierr = PetscOptionsBool("-use_mat_nearnullspace","MatNearNullSpace API test","",use_nearnullspace,&use_nearnullspace,NULL);CHKERRQ(ierr);
247*c4762a1bSJed Brown     ierr = PetscOptionsBool("-attach_mat_nearnullspace","MatNearNullSpace API test (via MatSetNearNullSpace)","",attach_nearnullspace,&attach_nearnullspace,NULL);CHKERRQ(ierr);
248*c4762a1bSJed Brown     ierr = PetscOptionsInt("-run_type","0: twisting load on cantalever, 1: 3rd order accurate convergence test","",run_type,&run_type,NULL);CHKERRQ(ierr);
249*c4762a1bSJed Brown     i = 3;
250*c4762a1bSJed Brown     ierr = PetscOptionsInt("-mat_block_size","","",i,&i,&flg);CHKERRQ(ierr);
251*c4762a1bSJed Brown     if (!flg || i!=3) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_USER, "'-mat_block_size 3' must be set (%D) and = 3 (%D)",flg,flg? i : 3);
252*c4762a1bSJed Brown   }
253*c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
254*c4762a1bSJed Brown   ierr = PetscLogStageRegister("Mesh Setup", &stage[16]);CHKERRQ(ierr);
255*c4762a1bSJed Brown   for (iter=0 ; iter<max_conv_its ; iter++) {
256*c4762a1bSJed Brown     char str[] = "Solve 0";
257*c4762a1bSJed Brown     str[6] += iter;
258*c4762a1bSJed Brown     ierr = PetscLogStageRegister(str, &stage[iter]);CHKERRQ(ierr);
259*c4762a1bSJed Brown   }
260*c4762a1bSJed Brown   /* create DM, Plex calls DMSetup */
261*c4762a1bSJed Brown   ierr = PetscLogStagePush(stage[16]);CHKERRQ(ierr);
262*c4762a1bSJed Brown   ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, cells, NULL, NULL, NULL, PETSC_TRUE, &dm);CHKERRQ(ierr);
263*c4762a1bSJed Brown   {
264*c4762a1bSJed Brown     DMLabel         label;
265*c4762a1bSJed Brown     IS              is;
266*c4762a1bSJed Brown     ierr = DMCreateLabel(dm, "boundary");CHKERRQ(ierr);
267*c4762a1bSJed Brown     ierr = DMGetLabel(dm, "boundary", &label);CHKERRQ(ierr);
268*c4762a1bSJed Brown     ierr = DMPlexMarkBoundaryFaces(dm, 1, label);CHKERRQ(ierr);
269*c4762a1bSJed Brown     if (!run_type) {
270*c4762a1bSJed Brown       ierr = DMGetStratumIS(dm, "boundary", 1,  &is);CHKERRQ(ierr);
271*c4762a1bSJed Brown       ierr = DMCreateLabel(dm,"Faces");CHKERRQ(ierr);
272*c4762a1bSJed Brown       if (is) {
273*c4762a1bSJed Brown         PetscInt        d, f, Nf;
274*c4762a1bSJed Brown         const PetscInt *faces;
275*c4762a1bSJed Brown         PetscInt        csize;
276*c4762a1bSJed Brown         PetscSection    cs;
277*c4762a1bSJed Brown         Vec             coordinates ;
278*c4762a1bSJed Brown         DM              cdm;
279*c4762a1bSJed Brown         ierr = ISGetLocalSize(is, &Nf);CHKERRQ(ierr);
280*c4762a1bSJed Brown         ierr = ISGetIndices(is, &faces);CHKERRQ(ierr);
281*c4762a1bSJed Brown         ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
282*c4762a1bSJed Brown         ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
283*c4762a1bSJed Brown         ierr = DMGetLocalSection(cdm, &cs);CHKERRQ(ierr);
284*c4762a1bSJed Brown         /* Check for each boundary face if any component of its centroid is either 0.0 or 1.0 */
285*c4762a1bSJed Brown         for (f = 0; f < Nf; ++f) {
286*c4762a1bSJed Brown           PetscReal   faceCoord;
287*c4762a1bSJed Brown           PetscInt    b,v;
288*c4762a1bSJed Brown           PetscScalar *coords = NULL;
289*c4762a1bSJed Brown           PetscInt    Nv;
290*c4762a1bSJed Brown           ierr = DMPlexVecGetClosure(cdm, cs, coordinates, faces[f], &csize, &coords);CHKERRQ(ierr);
291*c4762a1bSJed Brown           Nv   = csize/dim; /* Calculate mean coordinate vector */
292*c4762a1bSJed Brown           for (d = 0; d < dim; ++d) {
293*c4762a1bSJed Brown             faceCoord = 0.0;
294*c4762a1bSJed Brown             for (v = 0; v < Nv; ++v) faceCoord += PetscRealPart(coords[v*dim+d]);
295*c4762a1bSJed Brown             faceCoord /= Nv;
296*c4762a1bSJed Brown             for (b = 0; b < 2; ++b) {
297*c4762a1bSJed Brown               if (PetscAbs(faceCoord - b) < PETSC_SMALL) { /* domain have not been set yet, still [0,1]^3 */
298*c4762a1bSJed Brown                 ierr = DMSetLabelValue(dm, "Faces", faces[f], d*2+b+1);CHKERRQ(ierr);
299*c4762a1bSJed Brown               }
300*c4762a1bSJed Brown             }
301*c4762a1bSJed Brown           }
302*c4762a1bSJed Brown           ierr = DMPlexVecRestoreClosure(cdm, cs, coordinates, faces[f], &csize, &coords);CHKERRQ(ierr);
303*c4762a1bSJed Brown         }
304*c4762a1bSJed Brown         ierr = ISRestoreIndices(is, &faces);CHKERRQ(ierr);
305*c4762a1bSJed Brown       }
306*c4762a1bSJed Brown       ierr = ISDestroy(&is);CHKERRQ(ierr);
307*c4762a1bSJed Brown       ierr = DMGetLabel(dm, "Faces", &label);CHKERRQ(ierr);
308*c4762a1bSJed Brown       ierr = DMPlexLabelComplete(dm, label);CHKERRQ(ierr);
309*c4762a1bSJed Brown     }
310*c4762a1bSJed Brown   }
311*c4762a1bSJed Brown   {
312*c4762a1bSJed Brown     PetscInt    dimEmbed, i;
313*c4762a1bSJed Brown     PetscInt    nCoords;
314*c4762a1bSJed Brown     PetscScalar *coords,bounds[] = {0,1,-.5,.5,-.5,.5,}; /* x_min,x_max,y_min,y_max */
315*c4762a1bSJed Brown     Vec         coordinates;
316*c4762a1bSJed Brown     bounds[1] = Lx;
317*c4762a1bSJed Brown     if (run_type==1) {
318*c4762a1bSJed Brown       for (i = 0; i < 2*dim; i++) bounds[i] = (i%2) ? 1 : 0;
319*c4762a1bSJed Brown     }
320*c4762a1bSJed Brown     ierr = DMGetCoordinatesLocal(dm,&coordinates);CHKERRQ(ierr);
321*c4762a1bSJed Brown     ierr = DMGetCoordinateDim(dm,&dimEmbed);CHKERRQ(ierr);
322*c4762a1bSJed Brown     if (dimEmbed != dim) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"dimEmbed != dim %D",dimEmbed);
323*c4762a1bSJed Brown     ierr = VecGetLocalSize(coordinates,&nCoords);CHKERRQ(ierr);
324*c4762a1bSJed Brown     if (nCoords % dimEmbed) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Coordinate vector the wrong size");
325*c4762a1bSJed Brown     ierr = VecGetArray(coordinates,&coords);CHKERRQ(ierr);
326*c4762a1bSJed Brown     for (i = 0; i < nCoords; i += dimEmbed) {
327*c4762a1bSJed Brown       PetscInt    j;
328*c4762a1bSJed Brown       PetscScalar *coord = &coords[i];
329*c4762a1bSJed Brown       for (j = 0; j < dimEmbed; j++) {
330*c4762a1bSJed Brown         coord[j] = bounds[2 * j] + coord[j] * (bounds[2 * j + 1] - bounds[2 * j]);
331*c4762a1bSJed Brown       }
332*c4762a1bSJed Brown     }
333*c4762a1bSJed Brown     ierr = VecRestoreArray(coordinates,&coords);CHKERRQ(ierr);
334*c4762a1bSJed Brown     ierr = DMSetCoordinatesLocal(dm,coordinates);CHKERRQ(ierr);
335*c4762a1bSJed Brown   }
336*c4762a1bSJed Brown 
337*c4762a1bSJed Brown   /* convert to p4est, and distribute */
338*c4762a1bSJed Brown 
339*c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");CHKERRQ(ierr);
340*c4762a1bSJed Brown   ierr = PetscOptionsFList("-dm_type","Convert DMPlex to another format (should not be Plex!)","ex56.c",DMList,DMPLEX,convType,256,&flg);CHKERRQ(ierr);
341*c4762a1bSJed Brown   ierr = PetscOptionsEnd();
342*c4762a1bSJed Brown   if (flg) {
343*c4762a1bSJed Brown     DM newdm;
344*c4762a1bSJed Brown     ierr = DMConvert(dm,convType,&newdm);CHKERRQ(ierr);
345*c4762a1bSJed Brown     if (newdm) {
346*c4762a1bSJed Brown       const char *prefix;
347*c4762a1bSJed Brown       PetscBool isForest;
348*c4762a1bSJed Brown       ierr = PetscObjectGetOptionsPrefix((PetscObject)dm,&prefix);CHKERRQ(ierr);
349*c4762a1bSJed Brown       ierr = PetscObjectSetOptionsPrefix((PetscObject)newdm,prefix);CHKERRQ(ierr);
350*c4762a1bSJed Brown       ierr = DMIsForest(newdm,&isForest);CHKERRQ(ierr);
351*c4762a1bSJed Brown       if (!isForest) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_USER, "Converted to non Forest?");
352*c4762a1bSJed Brown       ierr = DMDestroy(&dm);CHKERRQ(ierr);
353*c4762a1bSJed Brown       dm   = newdm;
354*c4762a1bSJed Brown     } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_USER, "Convert failed?");
355*c4762a1bSJed Brown   } else {
356*c4762a1bSJed Brown     PetscPartitioner part;
357*c4762a1bSJed Brown     /* Plex Distribute mesh over processes */
358*c4762a1bSJed Brown     ierr = DMPlexGetPartitioner(dm,&part);CHKERRQ(ierr);
359*c4762a1bSJed Brown     ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
360*c4762a1bSJed Brown     ierr = DMPlexDistribute(dm, 0, NULL, &distdm);CHKERRQ(ierr);
361*c4762a1bSJed Brown     if (distdm) {
362*c4762a1bSJed Brown       const char *prefix;
363*c4762a1bSJed Brown       ierr = PetscObjectGetOptionsPrefix((PetscObject)dm,&prefix);CHKERRQ(ierr);
364*c4762a1bSJed Brown       ierr = PetscObjectSetOptionsPrefix((PetscObject)distdm,prefix);CHKERRQ(ierr);
365*c4762a1bSJed Brown       ierr = DMDestroy(&dm);CHKERRQ(ierr);
366*c4762a1bSJed Brown       dm   = distdm;
367*c4762a1bSJed Brown     }
368*c4762a1bSJed Brown   }
369*c4762a1bSJed Brown   ierr = PetscLogStagePop();CHKERRQ(ierr);
370*c4762a1bSJed Brown   basedm = dm; dm = NULL;
371*c4762a1bSJed Brown 
372*c4762a1bSJed Brown   for (iter=0 ; iter<max_conv_its ; iter++) {
373*c4762a1bSJed Brown     ierr = PetscLogStagePush(stage[16]);CHKERRQ(ierr);
374*c4762a1bSJed Brown     /* make new DM */
375*c4762a1bSJed Brown     ierr = DMClone(basedm, &dm);CHKERRQ(ierr);
376*c4762a1bSJed Brown     ierr = PetscObjectSetOptionsPrefix((PetscObject) dm, "ex56_");CHKERRQ(ierr);
377*c4762a1bSJed Brown     ierr = PetscObjectSetName( (PetscObject)dm,"Mesh");CHKERRQ(ierr);
378*c4762a1bSJed Brown     if (max_conv_its > 1) {
379*c4762a1bSJed Brown       /* If max_conv_its == 1, then we are not doing a convergence study. In this case, do not overwrite the -ex56_dm_refine
380*c4762a1bSJed Brown        * options with the values in the options[] array; instead, use the user-specified value. */
381*c4762a1bSJed Brown       ierr = PetscOptionsClearValue(NULL,"-ex56_dm_refine");CHKERRQ(ierr);
382*c4762a1bSJed Brown       ierr = PetscOptionsInsertString(NULL,options[iter]);CHKERRQ(ierr);
383*c4762a1bSJed Brown     }
384*c4762a1bSJed Brown     ierr = DMSetFromOptions(dm);CHKERRQ(ierr); /* refinement done here in Plex, p4est */
385*c4762a1bSJed Brown     /* snes */
386*c4762a1bSJed Brown     ierr = SNESCreate(comm, &snes);CHKERRQ(ierr);
387*c4762a1bSJed Brown     ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
388*c4762a1bSJed Brown     /* fem */
389*c4762a1bSJed Brown     {
390*c4762a1bSJed Brown       const PetscInt Ncomp = dim;
391*c4762a1bSJed Brown       const PetscInt components[] = {0,1,2};
392*c4762a1bSJed Brown       const PetscInt Nfid = 1, Npid = 1;
393*c4762a1bSJed Brown       const PetscInt fid[] = {1}; /* The fixed faces (x=0) */
394*c4762a1bSJed Brown       const PetscInt pid[] = {2}; /* The faces with loading (x=L_x) */
395*c4762a1bSJed Brown       PetscFE        fe;
396*c4762a1bSJed Brown       PetscDS        prob;
397*c4762a1bSJed Brown       DM             cdm = dm;
398*c4762a1bSJed Brown 
399*c4762a1bSJed Brown       ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, dim, PETSC_FALSE, NULL, PETSC_DECIDE, &fe);CHKERRQ(ierr); /* elasticity */
400*c4762a1bSJed Brown       ierr = PetscObjectSetName((PetscObject) fe, "deformation");CHKERRQ(ierr);
401*c4762a1bSJed Brown       /* FEM prob */
402*c4762a1bSJed Brown       ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr);
403*c4762a1bSJed Brown       ierr = DMCreateDS(dm);CHKERRQ(ierr);
404*c4762a1bSJed Brown       ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
405*c4762a1bSJed Brown       /* setup problem */
406*c4762a1bSJed Brown       if (run_type==1) {
407*c4762a1bSJed Brown         ierr = PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu_3d);CHKERRQ(ierr);
408*c4762a1bSJed Brown         ierr = PetscDSSetResidual(prob, 0, f0_u_x4, f1_u_3d);CHKERRQ(ierr);
409*c4762a1bSJed Brown       } else {
410*c4762a1bSJed Brown         ierr = PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu_3d_alpha);CHKERRQ(ierr);
411*c4762a1bSJed Brown         ierr = PetscDSSetResidual(prob, 0, f0_u, f1_u_3d_alpha);CHKERRQ(ierr);
412*c4762a1bSJed Brown         ierr = PetscDSSetBdResidual(prob, 0, f0_bd_u_3d, f1_bd_u);CHKERRQ(ierr);
413*c4762a1bSJed Brown       }
414*c4762a1bSJed Brown       /* bcs */
415*c4762a1bSJed Brown       if (run_type==1) {
416*c4762a1bSJed Brown         PetscInt id = 1;
417*c4762a1bSJed Brown         ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "boundary", 0, 0, NULL, (void (*)(void)) zero, 1, &id, NULL);CHKERRQ(ierr);
418*c4762a1bSJed Brown       } else {
419*c4762a1bSJed Brown         ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "fixed", "Faces", 0, Ncomp, components, (void (*)(void)) zero, Nfid, fid, NULL);CHKERRQ(ierr);
420*c4762a1bSJed Brown         ierr = PetscDSAddBoundary(prob, DM_BC_NATURAL, "traction", "Faces", 0, Ncomp, components, NULL, Npid, pid, NULL);CHKERRQ(ierr);
421*c4762a1bSJed Brown       }
422*c4762a1bSJed Brown       while (cdm) {
423*c4762a1bSJed Brown         ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
424*c4762a1bSJed Brown         ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
425*c4762a1bSJed Brown       }
426*c4762a1bSJed Brown       ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
427*c4762a1bSJed Brown     }
428*c4762a1bSJed Brown     /* vecs & mat */
429*c4762a1bSJed Brown     ierr = DMCreateGlobalVector(dm,&xx);CHKERRQ(ierr);
430*c4762a1bSJed Brown     ierr = VecDuplicate(xx, &bb);CHKERRQ(ierr);
431*c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject) bb, "b");CHKERRQ(ierr);
432*c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject) xx, "u");CHKERRQ(ierr);
433*c4762a1bSJed Brown     ierr = DMCreateMatrix(dm, &Amat);CHKERRQ(ierr);
434*c4762a1bSJed Brown     ierr = MatSetOption(Amat,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);        /* Some matrix kernels can take advantage of symmetry if we set this. */
435*c4762a1bSJed Brown     ierr = MatSetOption(Amat,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);CHKERRQ(ierr); /* Inform PETSc that Amat is always symmetric, so info set above isn't lost. */
436*c4762a1bSJed Brown     ierr = VecGetSize(bb,&N);CHKERRQ(ierr);
437*c4762a1bSJed Brown     local_sizes[iter] = N;
438*c4762a1bSJed Brown     ierr = PetscInfo2(snes,"%D global equations, %D vertices\n",N,N/dim);CHKERRQ(ierr);
439*c4762a1bSJed Brown     if ((use_nearnullspace || attach_nearnullspace) && N/dim > 1) {
440*c4762a1bSJed Brown       /* Set up the near null space (a.k.a. rigid body modes) that will be used by the multigrid preconditioner */
441*c4762a1bSJed Brown       DM           subdm;
442*c4762a1bSJed Brown       MatNullSpace nearNullSpace;
443*c4762a1bSJed Brown       PetscInt     fields = 0;
444*c4762a1bSJed Brown       PetscObject  deformation;
445*c4762a1bSJed Brown       ierr = DMCreateSubDM(dm, 1, &fields, NULL, &subdm);CHKERRQ(ierr);
446*c4762a1bSJed Brown       ierr = DMPlexCreateRigidBody(subdm, &nearNullSpace);CHKERRQ(ierr);
447*c4762a1bSJed Brown       ierr = DMGetField(dm, 0, NULL, &deformation);CHKERRQ(ierr);
448*c4762a1bSJed Brown       ierr = PetscObjectCompose(deformation, "nearnullspace", (PetscObject) nearNullSpace);CHKERRQ(ierr);
449*c4762a1bSJed Brown       ierr = DMDestroy(&subdm);CHKERRQ(ierr);
450*c4762a1bSJed Brown       if (attach_nearnullspace) {
451*c4762a1bSJed Brown         ierr = MatSetNearNullSpace(Amat,nearNullSpace);CHKERRQ(ierr);
452*c4762a1bSJed Brown       }
453*c4762a1bSJed Brown       ierr = MatNullSpaceDestroy(&nearNullSpace);CHKERRQ(ierr); /* created by DM and destroyed by Mat */
454*c4762a1bSJed Brown     }
455*c4762a1bSJed Brown     ierr = DMPlexSetSNESLocalFEM(dm,NULL,NULL,NULL);CHKERRQ(ierr);
456*c4762a1bSJed Brown     ierr = SNESSetJacobian(snes, Amat, Amat, NULL, NULL);CHKERRQ(ierr);
457*c4762a1bSJed Brown     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
458*c4762a1bSJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr);
459*c4762a1bSJed Brown     ierr = PetscLogStagePop();CHKERRQ(ierr);
460*c4762a1bSJed Brown     ierr = PetscLogStagePush(stage[16]);CHKERRQ(ierr);
461*c4762a1bSJed Brown     /* ksp */
462*c4762a1bSJed Brown     ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
463*c4762a1bSJed Brown     ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE);CHKERRQ(ierr);
464*c4762a1bSJed Brown     /* test BCs */
465*c4762a1bSJed Brown     ierr = VecZeroEntries(xx);CHKERRQ(ierr);
466*c4762a1bSJed Brown     if (test_nonzero_cols) {
467*c4762a1bSJed Brown       if (!rank) {
468*c4762a1bSJed Brown         ierr = VecSetValue(xx,0,1.0,INSERT_VALUES);CHKERRQ(ierr);
469*c4762a1bSJed Brown       }
470*c4762a1bSJed Brown       ierr = VecAssemblyBegin(xx);CHKERRQ(ierr);
471*c4762a1bSJed Brown       ierr = VecAssemblyEnd(xx);CHKERRQ(ierr);
472*c4762a1bSJed Brown     }
473*c4762a1bSJed Brown     ierr = VecZeroEntries(bb);CHKERRQ(ierr);
474*c4762a1bSJed Brown     ierr = VecGetSize(bb,&i);CHKERRQ(ierr);
475*c4762a1bSJed Brown     local_sizes[iter] = i;
476*c4762a1bSJed Brown     ierr = PetscInfo2(snes,"%D equations in vector, %D vertices\n",i,i/dim);CHKERRQ(ierr);
477*c4762a1bSJed Brown     ierr = PetscLogStagePop();CHKERRQ(ierr);
478*c4762a1bSJed Brown     /* solve */
479*c4762a1bSJed Brown     ierr = PetscLogStagePush(stage[iter]);CHKERRQ(ierr);
480*c4762a1bSJed Brown     ierr = SNESSolve(snes, bb, xx);CHKERRQ(ierr);
481*c4762a1bSJed Brown     ierr = PetscLogStagePop();CHKERRQ(ierr);
482*c4762a1bSJed Brown     ierr = VecNorm(xx,NORM_INFINITY,&mdisp[iter]);CHKERRQ(ierr);
483*c4762a1bSJed Brown     ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
484*c4762a1bSJed Brown     {
485*c4762a1bSJed Brown       PetscViewer       viewer = NULL;
486*c4762a1bSJed Brown       PetscViewerFormat fmt;
487*c4762a1bSJed Brown       ierr = PetscOptionsGetViewer(comm,NULL,"ex56_","-vec_view",&viewer,&fmt,&flg);CHKERRQ(ierr);
488*c4762a1bSJed Brown       if (flg) {
489*c4762a1bSJed Brown         ierr = PetscViewerPushFormat(viewer,fmt);CHKERRQ(ierr);
490*c4762a1bSJed Brown         ierr = VecView(xx,viewer);CHKERRQ(ierr);
491*c4762a1bSJed Brown         ierr = VecView(bb,viewer);CHKERRQ(ierr);
492*c4762a1bSJed Brown         ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
493*c4762a1bSJed Brown       }
494*c4762a1bSJed Brown       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
495*c4762a1bSJed Brown     }
496*c4762a1bSJed Brown     /* Free work space */
497*c4762a1bSJed Brown     ierr = DMDestroy(&dm);CHKERRQ(ierr);
498*c4762a1bSJed Brown     ierr = SNESDestroy(&snes);CHKERRQ(ierr);
499*c4762a1bSJed Brown     ierr = VecDestroy(&xx);CHKERRQ(ierr);
500*c4762a1bSJed Brown     ierr = VecDestroy(&bb);CHKERRQ(ierr);
501*c4762a1bSJed Brown     ierr = MatDestroy(&Amat);CHKERRQ(ierr);
502*c4762a1bSJed Brown   }
503*c4762a1bSJed Brown   ierr = DMDestroy(&basedm);CHKERRQ(ierr);
504*c4762a1bSJed Brown   if (run_type==1) err[0] = 59.975208 - mdisp[0]; /* error with what I think is the exact solution */
505*c4762a1bSJed Brown   else             err[0] = 171.038 - mdisp[0];
506*c4762a1bSJed Brown   for (iter=1 ; iter<max_conv_its ; iter++) {
507*c4762a1bSJed Brown     if (run_type==1) err[iter] = 59.975208 - mdisp[iter];
508*c4762a1bSJed Brown     else             err[iter] = 171.038 - mdisp[iter];
509*c4762a1bSJed Brown     ierr = 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],
510*c4762a1bSJed Brown                        (double)(mdisp[iter]-mdisp[iter-1]),(double)err[iter],(double)(PetscLogReal(err[iter-1]/err[iter])/PetscLogReal(2.)));CHKERRQ(ierr);
511*c4762a1bSJed Brown   }
512*c4762a1bSJed Brown 
513*c4762a1bSJed Brown   ierr = PetscFinalize();
514*c4762a1bSJed Brown   return ierr;
515*c4762a1bSJed Brown }
516*c4762a1bSJed Brown 
517*c4762a1bSJed Brown /*TEST
518*c4762a1bSJed Brown 
519*c4762a1bSJed Brown   test:
520*c4762a1bSJed Brown     suffix: 0
521*c4762a1bSJed Brown     nsize: 4
522*c4762a1bSJed Brown     requires: !single
523*c4762a1bSJed 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-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 -matptap_via scalable -ex56_dm_view -run_type 1
524*c4762a1bSJed Brown     timeoutfactor: 2
525*c4762a1bSJed Brown 
526*c4762a1bSJed Brown   # HYPRE PtAP broken with complex numbers
527*c4762a1bSJed Brown   test:
528*c4762a1bSJed Brown     suffix: hypre
529*c4762a1bSJed Brown     requires: hypre !single !complex
530*c4762a1bSJed Brown     nsize: 4
531*c4762a1bSJed Brown     args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -ex56_dm_refine 1 -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 -mat_block_size 3 -petscpartitioner_type simple
532*c4762a1bSJed Brown 
533*c4762a1bSJed Brown   test:
534*c4762a1bSJed Brown     suffix: ml
535*c4762a1bSJed Brown     requires: ml !single
536*c4762a1bSJed Brown     nsize: 4
537*c4762a1bSJed Brown     args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -ex56_dm_refine 1 -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_esteig_ksp_max_it 10 -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 -mg_levels_pc_type sor -mg_levels_esteig_ksp_type cg -mat_block_size 3 -petscpartitioner_type simple -use_mat_nearnullspace
538*c4762a1bSJed Brown 
539*c4762a1bSJed Brown   test:
540*c4762a1bSJed Brown     suffix: hpddm
541*c4762a1bSJed Brown     requires: hpddm slepc !single
542*c4762a1bSJed Brown     nsize: 4
543*c4762a1bSJed Brown     args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -ex56_dm_refine 1 -petscspace_degree 2 -ksp_type fgmres -ksp_monitor_short -ksp_converged_reason -ksp_rtol 1.e-8 -pc_type hpddm -mat_block_size 3 -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
544*c4762a1bSJed Brown 
545*c4762a1bSJed Brown   test:
546*c4762a1bSJed Brown     nsize: 4
547*c4762a1bSJed Brown     requires: parmetis !single
548*c4762a1bSJed 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-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 10 -pc_gamg_reuse_interpolation true -pc_gamg_square_graph 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -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.05 -mg_levels_pc_type jacobi -pc_gamg_mat_partitioning_type parmetis -mat_block_size 3 -run_type 1 -pc_gamg_repartition true
549*c4762a1bSJed Brown 
550*c4762a1bSJed Brown   test:
551*c4762a1bSJed Brown     suffix: bddc
552*c4762a1bSJed Brown     nsize: 4
553*c4762a1bSJed Brown     requires: !single
554*c4762a1bSJed Brown     args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -ex56_dm_refine 1 -petscspace_degree 2 -ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-8 -ksp_converged_reason -mat_block_size 3 -petscpartitioner_type simple -ex56_dm_mat_type is -matis_localmat_type {{sbaij baij aij}} -pc_type bddc
555*c4762a1bSJed Brown 
556*c4762a1bSJed Brown   testset:
557*c4762a1bSJed Brown     nsize: 4
558*c4762a1bSJed Brown     requires: !single
559*c4762a1bSJed Brown     args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -ex56_dm_refine 1 -petscspace_degree 2 -ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-10 -ksp_converged_reason -mat_block_size 3 -petscpartitioner_type simple -ex56_dm_mat_type is -matis_localmat_type aij -pc_type bddc -attach_mat_nearnullspace {{0 1}separate output}
560*c4762a1bSJed Brown     test:
561*c4762a1bSJed Brown       suffix: bddc_approx_gamg
562*c4762a1bSJed Brown       args: -pc_bddc_switch_static -prefix_push pc_bddc_dirichlet_ -approximate -pc_type gamg -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 -mg_levels_esteig_ksp_type cg -mg_levels_esteig_ksp_max_it 10 -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 -prefix_pop -prefix_push pc_bddc_neumann_ -approximate -pc_type gamg -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 -mg_levels_esteig_ksp_type cg -mg_levels_esteig_ksp_max_it 10 -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 -prefix_pop
563*c4762a1bSJed Brown     # HYPRE PtAP broken with complex numbers
564*c4762a1bSJed Brown     test:
565*c4762a1bSJed Brown       requires: hypre !complex
566*c4762a1bSJed Brown       suffix: bddc_approx_hypre
567*c4762a1bSJed 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
568*c4762a1bSJed Brown     test:
569*c4762a1bSJed Brown       requires: ml
570*c4762a1bSJed Brown       suffix: bddc_approx_ml
571*c4762a1bSJed Brown       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 -mg_levels_esteig_ksp_type cg -mg_levels_esteig_ksp_max_it 10 -prefix_pop -prefix_push pc_bddc_neumann_ -approximate -pc_type ml -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -mg_levels_esteig_ksp_type cg -mg_levels_esteig_ksp_max_it 10 -prefix_pop
572*c4762a1bSJed Brown 
573*c4762a1bSJed Brown   test:
574*c4762a1bSJed Brown     suffix: fetidp
575*c4762a1bSJed Brown     nsize: 4
576*c4762a1bSJed Brown     requires: !single
577*c4762a1bSJed Brown     args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -ex56_dm_refine 1 -petscspace_degree 2 -ksp_type fetidp -fetidp_ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-8 -ksp_converged_reason -mat_block_size 3 -petscpartitioner_type simple -ex56_dm_mat_type is -matis_localmat_type {{sbaij baij aij}}
578*c4762a1bSJed Brown 
579*c4762a1bSJed Brown   test:
580*c4762a1bSJed Brown     suffix: bddc_elast
581*c4762a1bSJed Brown     nsize: 4
582*c4762a1bSJed Brown     requires: !single
583*c4762a1bSJed Brown     args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -ex56_dm_refine 1 -petscspace_degree 2 -ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-8 -ksp_converged_reason -mat_block_size 3 -petscpartitioner_type simple -ex56_dm_mat_type is -matis_localmat_type sbaij -pc_type bddc -pc_bddc_monolithic -attach_mat_nearnullspace
584*c4762a1bSJed Brown 
585*c4762a1bSJed Brown   test:
586*c4762a1bSJed Brown     suffix: fetidp_elast
587*c4762a1bSJed Brown     nsize: 4
588*c4762a1bSJed Brown     requires: !single
589*c4762a1bSJed Brown     args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -ex56_dm_refine 1 -petscspace_degree 2 -ksp_type fetidp -fetidp_ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-8 -ksp_converged_reason -mat_block_size 3 -petscpartitioner_type simple -ex56_dm_mat_type is -matis_localmat_type sbaij -fetidp_bddc_pc_bddc_monolithic -attach_mat_nearnullspace
590*c4762a1bSJed Brown 
591*c4762a1bSJed Brown   test:
592*c4762a1bSJed Brown     suffix: cuda
593*c4762a1bSJed Brown     nsize: 2
594*c4762a1bSJed Brown     requires: cuda !single
595*c4762a1bSJed 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-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 10 -pc_gamg_reuse_interpolation true -pc_gamg_square_graph 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -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.05 -mg_levels_pc_type jacobi -mat_block_size 3 -run_type 1 -ex56_dm_mat_type aijcusparse -ex56_dm_vec_type cuda -ksp_monitor_short
596*c4762a1bSJed Brown 
597*c4762a1bSJed Brown TEST*/
598