1*c4762a1bSJed Brown static char help[] = "Stokes Problem in 2d and 3d with simplicial finite elements.\n\ 2*c4762a1bSJed Brown We solve the Stokes problem in a rectangular\n\ 3*c4762a1bSJed Brown domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n"; 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown /* 6*c4762a1bSJed Brown The isoviscous Stokes problem, which we discretize using the finite 7*c4762a1bSJed Brown element method on an unstructured mesh. The weak form equations are 8*c4762a1bSJed Brown 9*c4762a1bSJed Brown < \nabla v, \nabla u + {\nabla u}^T > - < \nabla\cdot v, p > + < v, f > = 0 10*c4762a1bSJed Brown < q, \nabla\cdot u > = 0 11*c4762a1bSJed Brown 12*c4762a1bSJed Brown Viewing: 13*c4762a1bSJed Brown 14*c4762a1bSJed Brown To produce nice output, use 15*c4762a1bSJed Brown 16*c4762a1bSJed Brown -dm_refine 3 -show_error -dm_view hdf5:sol1.h5 -error_vec_view hdf5:sol1.h5::append -sol_vec_view hdf5:sol1.h5::append -exact_vec_view hdf5:sol1.h5::append 17*c4762a1bSJed Brown 18*c4762a1bSJed Brown You can get a LaTeX view of the mesh, with point numbering using 19*c4762a1bSJed Brown 20*c4762a1bSJed Brown -dm_view :mesh.tex:ascii_latex -dm_plex_view_scale 8.0 21*c4762a1bSJed Brown 22*c4762a1bSJed Brown The data layout can be viewed using 23*c4762a1bSJed Brown 24*c4762a1bSJed Brown -dm_petscsection_view 25*c4762a1bSJed Brown 26*c4762a1bSJed Brown Lots of information about the FEM assembly can be printed using 27*c4762a1bSJed Brown 28*c4762a1bSJed Brown -dm_plex_print_fem 2 29*c4762a1bSJed Brown 30*c4762a1bSJed Brown Field Data: 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown DMPLEX data is organized by point, and the closure operation just stacks up the 33*c4762a1bSJed Brown data from each sieve point in the closure. Thus, for a P_2-P_1 Stokes element, we 34*c4762a1bSJed Brown have 35*c4762a1bSJed Brown 36*c4762a1bSJed Brown cl{e} = {f e_0 e_1 e_2 v_0 v_1 v_2} 37*c4762a1bSJed Brown x = [u_{e_0} v_{e_0} u_{e_1} v_{e_1} u_{e_2} v_{e_2} u_{v_0} v_{v_0} p_{v_0} u_{v_1} v_{v_1} p_{v_1} u_{v_2} v_{v_2} p_{v_2}] 38*c4762a1bSJed Brown 39*c4762a1bSJed Brown The problem here is that we would like to loop over each field separately for 40*c4762a1bSJed Brown integration. Therefore, the closure visitor in DMPlexVecGetClosure() reorders 41*c4762a1bSJed Brown the data so that each field is contiguous 42*c4762a1bSJed Brown 43*c4762a1bSJed Brown x' = [u_{e_0} v_{e_0} u_{e_1} v_{e_1} u_{e_2} v_{e_2} u_{v_0} v_{v_0} u_{v_1} v_{v_1} u_{v_2} v_{v_2} p_{v_0} p_{v_1} p_{v_2}] 44*c4762a1bSJed Brown 45*c4762a1bSJed Brown Likewise, DMPlexVecSetClosure() takes data partitioned by field, and correctly 46*c4762a1bSJed Brown puts it into the Sieve ordering. 47*c4762a1bSJed Brown 48*c4762a1bSJed Brown TODO: 49*c4762a1bSJed Brown - Update FETI test output 50*c4762a1bSJed Brown - Reorder mesh 51*c4762a1bSJed Brown - Check the q1-p0 Vanka domains are correct (I think its correct) 52*c4762a1bSJed Brown - Check scaling of iterates, right now it is bad 53*c4762a1bSJed Brown - Check the q2-q1 domains since convergence is bad 54*c4762a1bSJed Brown - Ask Patrick about domains 55*c4762a1bSJed Brown - Plot residual by fields after each smoother iterate 56*c4762a1bSJed Brown - Get Diskin checks going 57*c4762a1bSJed Brown */ 58*c4762a1bSJed Brown 59*c4762a1bSJed Brown #include <petscdmplex.h> 60*c4762a1bSJed Brown #include <petscsnes.h> 61*c4762a1bSJed Brown #include <petscds.h> 62*c4762a1bSJed Brown 63*c4762a1bSJed Brown PetscInt spatialDim = 0; 64*c4762a1bSJed Brown 65*c4762a1bSJed Brown typedef enum {NEUMANN, DIRICHLET, NUM_BC_TYPES} BCType; 66*c4762a1bSJed Brown const char *bcTypes[NUM_BC_TYPES+1] = {"neumann", "dirichlet", "unknown"}; 67*c4762a1bSJed Brown typedef enum {RUN_FULL, RUN_TEST, NUM_RUN_TYPES} RunType; 68*c4762a1bSJed Brown const char *runTypes[NUM_RUN_TYPES+1] = {"full", "test", "unknown"}; 69*c4762a1bSJed Brown typedef enum {SOL_QUADRATIC, SOL_CUBIC, SOL_TRIG, NUM_SOL_TYPES} SolType; 70*c4762a1bSJed Brown const char *solTypes[NUM_SOL_TYPES+1] = {"quadratic", "cubic", "trig", "unknown"}; 71*c4762a1bSJed Brown 72*c4762a1bSJed Brown typedef struct { 73*c4762a1bSJed Brown PetscInt debug; /* The debugging level */ 74*c4762a1bSJed Brown RunType runType; /* Whether to run tests, or solve the full problem */ 75*c4762a1bSJed Brown PetscLogEvent createMeshEvent; 76*c4762a1bSJed Brown PetscBool showInitial, showError; 77*c4762a1bSJed Brown /* Domain and mesh definition */ 78*c4762a1bSJed Brown PetscInt dim; /* The topological mesh dimension */ 79*c4762a1bSJed Brown PetscBool interpolate; /* Generate intermediate mesh elements */ 80*c4762a1bSJed Brown PetscBool simplex; /* Use simplices or tensor product cells */ 81*c4762a1bSJed Brown PetscInt cells[3]; /* The initial domain division */ 82*c4762a1bSJed Brown PetscReal refinementLimit; /* The largest allowable cell volume */ 83*c4762a1bSJed Brown PetscBool testPartition; /* Use a fixed partitioning for testing */ 84*c4762a1bSJed Brown /* Problem definition */ 85*c4762a1bSJed Brown BCType bcType; 86*c4762a1bSJed Brown SolType solType; 87*c4762a1bSJed Brown PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 88*c4762a1bSJed Brown } AppCtx; 89*c4762a1bSJed Brown 90*c4762a1bSJed Brown PetscErrorCode zero_scalar(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 91*c4762a1bSJed Brown { 92*c4762a1bSJed Brown u[0] = 0.0; 93*c4762a1bSJed Brown return 0; 94*c4762a1bSJed Brown } 95*c4762a1bSJed Brown PetscErrorCode zero_vector(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 96*c4762a1bSJed Brown { 97*c4762a1bSJed Brown PetscInt d; 98*c4762a1bSJed Brown for (d = 0; d < spatialDim; ++d) u[d] = 0.0; 99*c4762a1bSJed Brown return 0; 100*c4762a1bSJed Brown } 101*c4762a1bSJed Brown 102*c4762a1bSJed Brown /* 103*c4762a1bSJed Brown In 2D we use exact solution: 104*c4762a1bSJed Brown 105*c4762a1bSJed Brown u = x^2 + y^2 106*c4762a1bSJed Brown v = 2 x^2 - 2xy 107*c4762a1bSJed Brown p = x + y - 1 108*c4762a1bSJed Brown f_x = f_y = 3 109*c4762a1bSJed Brown 110*c4762a1bSJed Brown so that 111*c4762a1bSJed Brown 112*c4762a1bSJed Brown -\Delta u + \nabla p + f = <-4, -4> + <1, 1> + <3, 3> = 0 113*c4762a1bSJed Brown \nabla \cdot u = 2x - 2x = 0 114*c4762a1bSJed Brown */ 115*c4762a1bSJed Brown PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 116*c4762a1bSJed Brown { 117*c4762a1bSJed Brown u[0] = x[0]*x[0] + x[1]*x[1]; 118*c4762a1bSJed Brown u[1] = 2.0*x[0]*x[0] - 2.0*x[0]*x[1]; 119*c4762a1bSJed Brown return 0; 120*c4762a1bSJed Brown } 121*c4762a1bSJed Brown 122*c4762a1bSJed Brown PetscErrorCode linear_p_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *p, void *ctx) 123*c4762a1bSJed Brown { 124*c4762a1bSJed Brown *p = x[0] + x[1] - 1.0; 125*c4762a1bSJed Brown return 0; 126*c4762a1bSJed Brown } 127*c4762a1bSJed Brown PetscErrorCode constant_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *p, void *ctx) 128*c4762a1bSJed Brown { 129*c4762a1bSJed Brown *p = 1.0; 130*c4762a1bSJed Brown return 0; 131*c4762a1bSJed Brown } 132*c4762a1bSJed Brown 133*c4762a1bSJed Brown /* 134*c4762a1bSJed Brown In 2D we use exact solution: 135*c4762a1bSJed Brown 136*c4762a1bSJed Brown u = x^3 + y^3 137*c4762a1bSJed Brown v = 2 x^3 - 3 x^2 y 138*c4762a1bSJed Brown p = 3/2 x^2 + 3/2 y^2 - 1 139*c4762a1bSJed Brown f_x = 6 (x + y) 140*c4762a1bSJed Brown f_y = 12 x - 3 y 141*c4762a1bSJed Brown 142*c4762a1bSJed Brown so that 143*c4762a1bSJed Brown 144*c4762a1bSJed Brown -\Delta u + \nabla p + f = <-6 x - 6 y, -12 x + 6 y> + <3 x, 3 y> + <6 (x + y), 12 x - 6 y> = 0 145*c4762a1bSJed Brown \nabla \cdot u = 3 x^2 - 3 x^2 = 0 146*c4762a1bSJed Brown */ 147*c4762a1bSJed Brown PetscErrorCode cubic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 148*c4762a1bSJed Brown { 149*c4762a1bSJed Brown u[0] = x[0]*x[0]*x[0] + x[1]*x[1]*x[1]; 150*c4762a1bSJed Brown u[1] = 2.0*x[0]*x[0]*x[0] - 3.0*x[0]*x[0]*x[1]; 151*c4762a1bSJed Brown return 0; 152*c4762a1bSJed Brown } 153*c4762a1bSJed Brown 154*c4762a1bSJed Brown PetscErrorCode quadratic_p_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *p, void *ctx) 155*c4762a1bSJed Brown { 156*c4762a1bSJed Brown *p = 3.0*x[0]*x[0]/2.0 + 3.0*x[1]*x[1]/2.0 - 1.0; 157*c4762a1bSJed Brown return 0; 158*c4762a1bSJed Brown } 159*c4762a1bSJed Brown 160*c4762a1bSJed Brown /* 161*c4762a1bSJed Brown In 2D we use exact solution: 162*c4762a1bSJed Brown 163*c4762a1bSJed Brown u = sin(n pi x) + y^2 164*c4762a1bSJed Brown v = -sin(n pi y) 165*c4762a1bSJed Brown p = 3/2 x^2 + 3/2 y^2 - 1 166*c4762a1bSJed Brown f_x = 4 - 3x - n^2 pi^2 sin (n pi x) 167*c4762a1bSJed Brown f_y = - 3y + n^2 pi^2 sin(n pi y) 168*c4762a1bSJed Brown 169*c4762a1bSJed Brown so that 170*c4762a1bSJed Brown 171*c4762a1bSJed Brown -\Delta u + \nabla p + f = <n^2 pi^2 sin (n pi x) - 4, -n^2 pi^2 sin(n pi y)> + <3 x, 3 y> + <4 - 3x - n^2 pi^2 sin (n pi x), -3y + n^2 pi^2 sin(n pi y)> = 0 172*c4762a1bSJed Brown \nabla \cdot u = n pi cos(n pi x) - n pi cos(n pi y) = 0 173*c4762a1bSJed Brown */ 174*c4762a1bSJed Brown PetscErrorCode trig_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 175*c4762a1bSJed Brown { 176*c4762a1bSJed Brown const PetscReal n = 1.0; 177*c4762a1bSJed Brown 178*c4762a1bSJed Brown u[0] = PetscSinReal(n*PETSC_PI*x[0]) + x[1]*x[1]; 179*c4762a1bSJed Brown u[1] = -PetscSinReal(n*PETSC_PI*x[1]); 180*c4762a1bSJed Brown return 0; 181*c4762a1bSJed Brown } 182*c4762a1bSJed Brown 183*c4762a1bSJed Brown void f0_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 184*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 185*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 186*c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 187*c4762a1bSJed Brown { 188*c4762a1bSJed Brown PetscInt c; 189*c4762a1bSJed Brown for (c = 0; c < dim; ++c) f0[c] = 3.0; 190*c4762a1bSJed Brown } 191*c4762a1bSJed Brown 192*c4762a1bSJed Brown void f0_cubic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 193*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 194*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 195*c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 196*c4762a1bSJed Brown { 197*c4762a1bSJed Brown f0[0] = 3.0*x[0] + 6.0*x[1]; 198*c4762a1bSJed Brown f0[1] = 12.0*x[0] - 9.0*x[1]; 199*c4762a1bSJed Brown } 200*c4762a1bSJed Brown 201*c4762a1bSJed Brown void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 202*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 203*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 204*c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 205*c4762a1bSJed Brown { 206*c4762a1bSJed Brown const PetscReal n = 1.0; 207*c4762a1bSJed Brown 208*c4762a1bSJed Brown f0[0] = 4.0 - 3.0*x[0] - PetscSqr(n*PETSC_PI)*PetscSinReal(n*PETSC_PI*x[0]); 209*c4762a1bSJed Brown f0[1] = -3.0*x[1] + PetscSqr(n*PETSC_PI)*PetscSinReal(n*PETSC_PI*x[1]); 210*c4762a1bSJed Brown } 211*c4762a1bSJed Brown 212*c4762a1bSJed Brown /* gradU[comp*dim+d] = {u_x, u_y, v_x, v_y} or {u_x, u_y, u_z, v_x, v_y, v_z, w_x, w_y, w_z} 213*c4762a1bSJed Brown u[Ncomp] = {p} */ 214*c4762a1bSJed Brown void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 215*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 216*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 217*c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 218*c4762a1bSJed Brown { 219*c4762a1bSJed Brown const PetscInt Ncomp = dim; 220*c4762a1bSJed Brown PetscInt comp, d; 221*c4762a1bSJed Brown 222*c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) { 223*c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 224*c4762a1bSJed Brown /* f1[comp*dim+d] = 0.5*(gradU[comp*dim+d] + gradU[d*dim+comp]); */ 225*c4762a1bSJed Brown f1[comp*dim+d] = u_x[comp*dim+d]; 226*c4762a1bSJed Brown } 227*c4762a1bSJed Brown f1[comp*dim+comp] -= u[Ncomp]; 228*c4762a1bSJed Brown } 229*c4762a1bSJed Brown } 230*c4762a1bSJed Brown 231*c4762a1bSJed Brown /* gradU[comp*dim+d] = {u_x, u_y, v_x, v_y} or {u_x, u_y, u_z, v_x, v_y, v_z, w_x, w_y, w_z} */ 232*c4762a1bSJed Brown void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux, 233*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 234*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 235*c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 236*c4762a1bSJed Brown { 237*c4762a1bSJed Brown PetscInt d; 238*c4762a1bSJed Brown for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d]; 239*c4762a1bSJed Brown } 240*c4762a1bSJed Brown 241*c4762a1bSJed Brown void f1_p(PetscInt dim, PetscInt Nf, PetscInt NfAux, 242*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 243*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 244*c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 245*c4762a1bSJed Brown { 246*c4762a1bSJed Brown PetscInt d; 247*c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = 0.0; 248*c4762a1bSJed Brown } 249*c4762a1bSJed Brown 250*c4762a1bSJed Brown /* < q, \nabla\cdot u > 251*c4762a1bSJed Brown NcompI = 1, NcompJ = dim */ 252*c4762a1bSJed Brown void g1_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 253*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 254*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 255*c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 256*c4762a1bSJed Brown { 257*c4762a1bSJed Brown PetscInt d; 258*c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; /* \frac{\partial\phi^{u_d}}{\partial x_d} */ 259*c4762a1bSJed Brown } 260*c4762a1bSJed Brown 261*c4762a1bSJed Brown /* -< \nabla\cdot v, p > 262*c4762a1bSJed Brown NcompI = dim, NcompJ = 1 */ 263*c4762a1bSJed Brown void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux, 264*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 265*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 266*c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 267*c4762a1bSJed Brown { 268*c4762a1bSJed Brown PetscInt d; 269*c4762a1bSJed Brown for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; /* \frac{\partial\psi^{u_d}}{\partial x_d} */ 270*c4762a1bSJed Brown } 271*c4762a1bSJed Brown 272*c4762a1bSJed Brown /* < \nabla v, \nabla u + {\nabla u}^T > 273*c4762a1bSJed Brown This just gives \nabla u, give the perdiagonal for the transpose */ 274*c4762a1bSJed Brown void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 275*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 276*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 277*c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 278*c4762a1bSJed Brown { 279*c4762a1bSJed Brown const PetscInt Ncomp = dim; 280*c4762a1bSJed Brown PetscInt compI, d; 281*c4762a1bSJed Brown 282*c4762a1bSJed Brown for (compI = 0; compI < Ncomp; ++compI) { 283*c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 284*c4762a1bSJed Brown g3[((compI*Ncomp+compI)*dim+d)*dim+d] = 1.0; 285*c4762a1bSJed Brown } 286*c4762a1bSJed Brown } 287*c4762a1bSJed Brown } 288*c4762a1bSJed Brown 289*c4762a1bSJed Brown /* 290*c4762a1bSJed Brown In 3D we use exact solution: 291*c4762a1bSJed Brown 292*c4762a1bSJed Brown u = x^2 + y^2 293*c4762a1bSJed Brown v = y^2 + z^2 294*c4762a1bSJed Brown w = x^2 + y^2 - 2(x+y)z 295*c4762a1bSJed Brown p = x + y + z - 3/2 296*c4762a1bSJed Brown f_x = f_y = f_z = 3 297*c4762a1bSJed Brown 298*c4762a1bSJed Brown so that 299*c4762a1bSJed Brown 300*c4762a1bSJed Brown -\Delta u + \nabla p + f = <-4, -4, -4> + <1, 1, 1> + <3, 3, 3> = 0 301*c4762a1bSJed Brown \nabla \cdot u = 2x + 2y - 2(x + y) = 0 302*c4762a1bSJed Brown */ 303*c4762a1bSJed Brown PetscErrorCode quadratic_u_3d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 304*c4762a1bSJed Brown { 305*c4762a1bSJed Brown u[0] = x[0]*x[0] + x[1]*x[1]; 306*c4762a1bSJed Brown u[1] = x[1]*x[1] + x[2]*x[2]; 307*c4762a1bSJed Brown u[2] = x[0]*x[0] + x[1]*x[1] - 2.0*(x[0] + x[1])*x[2]; 308*c4762a1bSJed Brown return 0; 309*c4762a1bSJed Brown } 310*c4762a1bSJed Brown 311*c4762a1bSJed Brown PetscErrorCode linear_p_3d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *p, void *ctx) 312*c4762a1bSJed Brown { 313*c4762a1bSJed Brown *p = x[0] + x[1] + x[2] - 1.5; 314*c4762a1bSJed Brown return 0; 315*c4762a1bSJed Brown } 316*c4762a1bSJed Brown 317*c4762a1bSJed Brown void pressure(PetscInt dim, PetscInt Nf, PetscInt NfAux, 318*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 319*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 320*c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar p[]) 321*c4762a1bSJed Brown { 322*c4762a1bSJed Brown p[0] = u[uOff[1]]; 323*c4762a1bSJed Brown } 324*c4762a1bSJed Brown 325*c4762a1bSJed Brown PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 326*c4762a1bSJed Brown { 327*c4762a1bSJed Brown PetscInt bc, run, sol, n; 328*c4762a1bSJed Brown PetscErrorCode ierr; 329*c4762a1bSJed Brown 330*c4762a1bSJed Brown PetscFunctionBeginUser; 331*c4762a1bSJed Brown options->debug = 0; 332*c4762a1bSJed Brown options->runType = RUN_FULL; 333*c4762a1bSJed Brown options->dim = 2; 334*c4762a1bSJed Brown options->interpolate = PETSC_FALSE; 335*c4762a1bSJed Brown options->simplex = PETSC_TRUE; 336*c4762a1bSJed Brown options->cells[0] = 3; 337*c4762a1bSJed Brown options->cells[1] = 3; 338*c4762a1bSJed Brown options->cells[2] = 3; 339*c4762a1bSJed Brown options->refinementLimit = 0.0; 340*c4762a1bSJed Brown options->testPartition = PETSC_FALSE; 341*c4762a1bSJed Brown options->bcType = DIRICHLET; 342*c4762a1bSJed Brown options->solType = SOL_QUADRATIC; 343*c4762a1bSJed Brown options->showInitial = PETSC_FALSE; 344*c4762a1bSJed Brown options->showError = PETSC_FALSE; 345*c4762a1bSJed Brown 346*c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Stokes Problem Options", "DMPLEX");CHKERRQ(ierr); 347*c4762a1bSJed Brown ierr = PetscOptionsInt("-debug", "The debugging level", "ex62.c", options->debug, &options->debug, NULL);CHKERRQ(ierr); 348*c4762a1bSJed Brown run = options->runType; 349*c4762a1bSJed Brown ierr = PetscOptionsEList("-run_type", "The run type", "ex62.c", runTypes, NUM_RUN_TYPES, runTypes[options->runType], &run, NULL);CHKERRQ(ierr); 350*c4762a1bSJed Brown options->runType = (RunType) run; 351*c4762a1bSJed Brown ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex62.c", options->dim, &options->dim, NULL);CHKERRQ(ierr); 352*c4762a1bSJed Brown spatialDim = options->dim; 353*c4762a1bSJed Brown ierr = PetscOptionsBool("-interpolate", "Generate intermediate mesh elements", "ex62.c", options->interpolate, &options->interpolate, NULL);CHKERRQ(ierr); 354*c4762a1bSJed Brown ierr = PetscOptionsBool("-simplex", "Use simplices or tensor product cells", "ex62.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr); 355*c4762a1bSJed Brown if (options->simplex) { 356*c4762a1bSJed Brown options->cells[0] = 4 - options->dim; 357*c4762a1bSJed Brown options->cells[1] = 4 - options->dim; 358*c4762a1bSJed Brown options->cells[2] = 4 - options->dim; 359*c4762a1bSJed Brown } 360*c4762a1bSJed Brown n = 3; 361*c4762a1bSJed Brown ierr = PetscOptionsIntArray("-cells", "The initial mesh division", "ex62.c", options->cells, &n, NULL);CHKERRQ(ierr); 362*c4762a1bSJed Brown ierr = PetscOptionsReal("-refinement_limit", "The largest allowable cell volume", "ex62.c", options->refinementLimit, &options->refinementLimit, NULL);CHKERRQ(ierr); 363*c4762a1bSJed Brown ierr = PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex62.c", options->testPartition, &options->testPartition, NULL);CHKERRQ(ierr); 364*c4762a1bSJed Brown bc = options->bcType; 365*c4762a1bSJed Brown ierr = PetscOptionsEList("-bc_type","Type of boundary condition","ex62.c", bcTypes, NUM_BC_TYPES, bcTypes[options->bcType], &bc, NULL);CHKERRQ(ierr); 366*c4762a1bSJed Brown options->bcType = (BCType) bc; 367*c4762a1bSJed Brown sol = options->solType; 368*c4762a1bSJed Brown ierr = PetscOptionsEList("-sol_type", "The solution type", "ex62.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL);CHKERRQ(ierr); 369*c4762a1bSJed Brown options->solType = (SolType) sol; 370*c4762a1bSJed Brown ierr = PetscOptionsBool("-show_initial", "Output the initial guess for verification", "ex62.c", options->showInitial, &options->showInitial, NULL);CHKERRQ(ierr); 371*c4762a1bSJed Brown ierr = PetscOptionsBool("-show_error", "Output the error for verification", "ex62.c", options->showError, &options->showError, NULL);CHKERRQ(ierr); 372*c4762a1bSJed Brown ierr = PetscOptionsEnd(); 373*c4762a1bSJed Brown 374*c4762a1bSJed Brown ierr = PetscLogEventRegister("CreateMesh", DM_CLASSID, &options->createMeshEvent);CHKERRQ(ierr); 375*c4762a1bSJed Brown PetscFunctionReturn(0); 376*c4762a1bSJed Brown } 377*c4762a1bSJed Brown 378*c4762a1bSJed Brown PetscErrorCode DMVecViewLocal(DM dm, Vec v) 379*c4762a1bSJed Brown { 380*c4762a1bSJed Brown Vec lv; 381*c4762a1bSJed Brown PetscErrorCode ierr; 382*c4762a1bSJed Brown 383*c4762a1bSJed Brown PetscFunctionBeginUser; 384*c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &lv);CHKERRQ(ierr); 385*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(dm, v, INSERT_VALUES, lv);CHKERRQ(ierr); 386*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(dm, v, INSERT_VALUES, lv);CHKERRQ(ierr); 387*c4762a1bSJed Brown ierr = DMPrintLocalVec(dm, "Local function", 0.0, lv);CHKERRQ(ierr); 388*c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &lv);CHKERRQ(ierr); 389*c4762a1bSJed Brown PetscFunctionReturn(0); 390*c4762a1bSJed Brown } 391*c4762a1bSJed Brown 392*c4762a1bSJed Brown PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 393*c4762a1bSJed Brown { 394*c4762a1bSJed Brown PetscInt dim = user->dim; 395*c4762a1bSJed Brown PetscBool interpolate = user->interpolate; 396*c4762a1bSJed Brown PetscReal refinementLimit = user->refinementLimit; 397*c4762a1bSJed Brown PetscErrorCode ierr; 398*c4762a1bSJed Brown 399*c4762a1bSJed Brown PetscFunctionBeginUser; 400*c4762a1bSJed Brown ierr = PetscLogEventBegin(user->createMeshEvent,0,0,0,0);CHKERRQ(ierr); 401*c4762a1bSJed Brown ierr = DMPlexCreateBoxMesh(comm, dim, user->simplex, user->cells, NULL, NULL, NULL, interpolate, dm);CHKERRQ(ierr); 402*c4762a1bSJed Brown { 403*c4762a1bSJed Brown DM refinedMesh = NULL; 404*c4762a1bSJed Brown DM distributedMesh = NULL; 405*c4762a1bSJed Brown 406*c4762a1bSJed Brown /* Refine mesh using a volume constraint */ 407*c4762a1bSJed Brown ierr = DMPlexSetRefinementLimit(*dm, refinementLimit);CHKERRQ(ierr); 408*c4762a1bSJed Brown if (user->simplex) {ierr = DMRefine(*dm, comm, &refinedMesh);CHKERRQ(ierr);} 409*c4762a1bSJed Brown if (refinedMesh) { 410*c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 411*c4762a1bSJed Brown *dm = refinedMesh; 412*c4762a1bSJed Brown } 413*c4762a1bSJed Brown /* Setup test partitioning */ 414*c4762a1bSJed Brown if (user->testPartition) { 415*c4762a1bSJed Brown PetscInt triSizes_n2[2] = {4, 4}; 416*c4762a1bSJed Brown PetscInt triPoints_n2[8] = {3, 5, 6, 7, 0, 1, 2, 4}; 417*c4762a1bSJed Brown PetscInt triSizes_n3[3] = {2, 3, 3}; 418*c4762a1bSJed Brown PetscInt triPoints_n3[8] = {3, 5, 1, 6, 7, 0, 2, 4}; 419*c4762a1bSJed Brown PetscInt triSizes_n5[5] = {1, 2, 2, 1, 2}; 420*c4762a1bSJed Brown PetscInt triPoints_n5[8] = {3, 5, 6, 4, 7, 0, 1, 2}; 421*c4762a1bSJed Brown PetscInt triSizes_ref_n2[2] = {8, 8}; 422*c4762a1bSJed Brown PetscInt triPoints_ref_n2[16] = {1, 5, 6, 7, 10, 11, 14, 15, 0, 2, 3, 4, 8, 9, 12, 13}; 423*c4762a1bSJed Brown PetscInt triSizes_ref_n3[3] = {5, 6, 5}; 424*c4762a1bSJed Brown PetscInt triPoints_ref_n3[16] = {1, 7, 10, 14, 15, 2, 6, 8, 11, 12, 13, 0, 3, 4, 5, 9}; 425*c4762a1bSJed Brown PetscInt triSizes_ref_n5[5] = {3, 4, 3, 3, 3}; 426*c4762a1bSJed Brown PetscInt triPoints_ref_n5[16] = {1, 7, 10, 2, 11, 13, 14, 5, 6, 15, 0, 8, 9, 3, 4, 12}; 427*c4762a1bSJed Brown PetscInt triSizes_ref_n5_d3[5] = {1, 1, 1, 1, 2}; 428*c4762a1bSJed Brown PetscInt triPoints_ref_n5_d3[6] = {0, 1, 2, 3, 4, 5}; 429*c4762a1bSJed Brown const PetscInt *sizes = NULL; 430*c4762a1bSJed Brown const PetscInt *points = NULL; 431*c4762a1bSJed Brown PetscPartitioner part; 432*c4762a1bSJed Brown PetscInt cEnd; 433*c4762a1bSJed Brown PetscMPIInt rank, size; 434*c4762a1bSJed Brown 435*c4762a1bSJed Brown ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 436*c4762a1bSJed Brown ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 437*c4762a1bSJed Brown ierr = DMPlexGetHeightStratum(*dm, 0, NULL, &cEnd);CHKERRQ(ierr); 438*c4762a1bSJed Brown if (!rank) { 439*c4762a1bSJed Brown if (dim == 2 && user->simplex && size == 2 && cEnd == 8) { 440*c4762a1bSJed Brown sizes = triSizes_n2; points = triPoints_n2; 441*c4762a1bSJed Brown } else if (dim == 2 && user->simplex && size == 3 && cEnd == 8) { 442*c4762a1bSJed Brown sizes = triSizes_n3; points = triPoints_n3; 443*c4762a1bSJed Brown } else if (dim == 2 && user->simplex && size == 5 && cEnd == 8) { 444*c4762a1bSJed Brown sizes = triSizes_n5; points = triPoints_n5; 445*c4762a1bSJed Brown } else if (dim == 2 && user->simplex && size == 2 && cEnd == 16) { 446*c4762a1bSJed Brown sizes = triSizes_ref_n2; points = triPoints_ref_n2; 447*c4762a1bSJed Brown } else if (dim == 2 && user->simplex && size == 3 && cEnd == 16) { 448*c4762a1bSJed Brown sizes = triSizes_ref_n3; points = triPoints_ref_n3; 449*c4762a1bSJed Brown } else if (dim == 2 && user->simplex && size == 5 && cEnd == 16) { 450*c4762a1bSJed Brown sizes = triSizes_ref_n5; points = triPoints_ref_n5; 451*c4762a1bSJed Brown } else if (dim == 3 && user->simplex && size == 5 && cEnd == 6) { 452*c4762a1bSJed Brown sizes = triSizes_ref_n5_d3; points = triPoints_ref_n5_d3; 453*c4762a1bSJed Brown } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "No stored partition matching run parameters"); 454*c4762a1bSJed Brown } 455*c4762a1bSJed Brown ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr); 456*c4762a1bSJed Brown ierr = PetscPartitionerSetType(part, PETSCPARTITIONERSHELL);CHKERRQ(ierr); 457*c4762a1bSJed Brown ierr = PetscPartitionerShellSetPartition(part, size, sizes, points);CHKERRQ(ierr); 458*c4762a1bSJed Brown } else { 459*c4762a1bSJed Brown PetscPartitioner part; 460*c4762a1bSJed Brown 461*c4762a1bSJed Brown ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr); 462*c4762a1bSJed Brown ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 463*c4762a1bSJed Brown } 464*c4762a1bSJed Brown /* Distribute mesh over processes */ 465*c4762a1bSJed Brown ierr = DMPlexDistribute(*dm, 0, NULL, &distributedMesh);CHKERRQ(ierr); 466*c4762a1bSJed Brown if (distributedMesh) { 467*c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 468*c4762a1bSJed Brown *dm = distributedMesh; 469*c4762a1bSJed Brown } 470*c4762a1bSJed Brown } 471*c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 472*c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 473*c4762a1bSJed Brown ierr = PetscLogEventEnd(user->createMeshEvent,0,0,0,0);CHKERRQ(ierr); 474*c4762a1bSJed Brown PetscFunctionReturn(0); 475*c4762a1bSJed Brown } 476*c4762a1bSJed Brown 477*c4762a1bSJed Brown PetscErrorCode SetupProblem(DM dm, AppCtx *user) 478*c4762a1bSJed Brown { 479*c4762a1bSJed Brown PetscDS prob; 480*c4762a1bSJed Brown const PetscInt id = 1; 481*c4762a1bSJed Brown PetscErrorCode ierr; 482*c4762a1bSJed Brown 483*c4762a1bSJed Brown PetscFunctionBeginUser; 484*c4762a1bSJed Brown ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 485*c4762a1bSJed Brown switch (user->solType) { 486*c4762a1bSJed Brown case SOL_QUADRATIC: 487*c4762a1bSJed Brown switch (user->dim) { 488*c4762a1bSJed Brown case 2: 489*c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 0, f0_quadratic_u, f1_u);CHKERRQ(ierr); 490*c4762a1bSJed Brown user->exactFuncs[0] = quadratic_u_2d; 491*c4762a1bSJed Brown user->exactFuncs[1] = linear_p_2d; 492*c4762a1bSJed Brown break; 493*c4762a1bSJed Brown case 3: 494*c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 0, f0_quadratic_u, f1_u);CHKERRQ(ierr); 495*c4762a1bSJed Brown user->exactFuncs[0] = quadratic_u_3d; 496*c4762a1bSJed Brown user->exactFuncs[1] = linear_p_3d; 497*c4762a1bSJed Brown break; 498*c4762a1bSJed Brown default: SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for quadratic solution", user->dim); 499*c4762a1bSJed Brown } 500*c4762a1bSJed Brown break; 501*c4762a1bSJed Brown case SOL_CUBIC: 502*c4762a1bSJed Brown switch (user->dim) { 503*c4762a1bSJed Brown case 2: 504*c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 0, f0_cubic_u, f1_u);CHKERRQ(ierr); 505*c4762a1bSJed Brown user->exactFuncs[0] = cubic_u_2d; 506*c4762a1bSJed Brown user->exactFuncs[1] = quadratic_p_2d; 507*c4762a1bSJed Brown break; 508*c4762a1bSJed Brown default: SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for quadratic solution", user->dim); 509*c4762a1bSJed Brown } 510*c4762a1bSJed Brown break; 511*c4762a1bSJed Brown case SOL_TRIG: 512*c4762a1bSJed Brown switch (user->dim) { 513*c4762a1bSJed Brown case 2: 514*c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 0, f0_trig_u, f1_u);CHKERRQ(ierr); 515*c4762a1bSJed Brown user->exactFuncs[0] = trig_u_2d; 516*c4762a1bSJed Brown user->exactFuncs[1] = quadratic_p_2d; 517*c4762a1bSJed Brown break; 518*c4762a1bSJed Brown default: SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for trigonometric solution", user->dim); 519*c4762a1bSJed Brown } 520*c4762a1bSJed Brown break; 521*c4762a1bSJed Brown default: SETERRQ2(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType); 522*c4762a1bSJed Brown } 523*c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 1, f0_p, f1_p);CHKERRQ(ierr); 524*c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr); 525*c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_up, NULL);CHKERRQ(ierr); 526*c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 1, 0, NULL, g1_pu, NULL, NULL);CHKERRQ(ierr); 527*c4762a1bSJed Brown 528*c4762a1bSJed Brown ierr = PetscDSAddBoundary(prob, user->bcType == DIRICHLET ? DM_BC_ESSENTIAL : DM_BC_NATURAL, "wall", user->bcType == NEUMANN ? "boundary" : "marker", 0, 0, NULL, (void (*)(void)) user->exactFuncs[0], 1, &id, user);CHKERRQ(ierr); 529*c4762a1bSJed Brown ierr = PetscDSSetExactSolution(prob, 0, user->exactFuncs[0], user);CHKERRQ(ierr); 530*c4762a1bSJed Brown ierr = PetscDSSetExactSolution(prob, 1, user->exactFuncs[1], user);CHKERRQ(ierr); 531*c4762a1bSJed Brown PetscFunctionReturn(0); 532*c4762a1bSJed Brown } 533*c4762a1bSJed Brown 534*c4762a1bSJed Brown PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 535*c4762a1bSJed Brown { 536*c4762a1bSJed Brown DM cdm = dm; 537*c4762a1bSJed Brown const PetscInt dim = user->dim; 538*c4762a1bSJed Brown PetscFE fe[2]; 539*c4762a1bSJed Brown MPI_Comm comm; 540*c4762a1bSJed Brown PetscErrorCode ierr; 541*c4762a1bSJed Brown 542*c4762a1bSJed Brown PetscFunctionBeginUser; 543*c4762a1bSJed Brown /* Create finite element */ 544*c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 545*c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, dim, user->simplex, "vel_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr); 546*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fe[0], "velocity");CHKERRQ(ierr); 547*c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, 1, user->simplex, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr); 548*c4762a1bSJed Brown ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr); 549*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr); 550*c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 551*c4762a1bSJed Brown ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr); 552*c4762a1bSJed Brown ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr); 553*c4762a1bSJed Brown ierr = DMCreateDS(dm);CHKERRQ(ierr); 554*c4762a1bSJed Brown ierr = SetupProblem(dm, user);CHKERRQ(ierr); 555*c4762a1bSJed Brown while (cdm) { 556*c4762a1bSJed Brown ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr); 557*c4762a1bSJed Brown ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 558*c4762a1bSJed Brown } 559*c4762a1bSJed Brown ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr); 560*c4762a1bSJed Brown ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr); 561*c4762a1bSJed Brown PetscFunctionReturn(0); 562*c4762a1bSJed Brown } 563*c4762a1bSJed Brown 564*c4762a1bSJed Brown static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt dummy, MatNullSpace *nullspace) 565*c4762a1bSJed Brown { 566*c4762a1bSJed Brown Vec vec; 567*c4762a1bSJed Brown PetscErrorCode (*funcs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void* ctx) = {zero_vector, constant_p}; 568*c4762a1bSJed Brown PetscErrorCode ierr; 569*c4762a1bSJed Brown 570*c4762a1bSJed Brown PetscFunctionBeginUser; 571*c4762a1bSJed Brown ierr = DMCreateGlobalVector(dm, &vec);CHKERRQ(ierr); 572*c4762a1bSJed Brown ierr = DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec);CHKERRQ(ierr); 573*c4762a1bSJed Brown ierr = VecNormalize(vec, NULL);CHKERRQ(ierr); 574*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) vec, "Pressure Null Space");CHKERRQ(ierr); 575*c4762a1bSJed Brown ierr = VecViewFromOptions(vec, NULL, "-pressure_nullspace_view");CHKERRQ(ierr); 576*c4762a1bSJed Brown ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_FALSE, 1, &vec, nullspace);CHKERRQ(ierr); 577*c4762a1bSJed Brown ierr = VecDestroy(&vec);CHKERRQ(ierr); 578*c4762a1bSJed Brown /* New style for field null spaces */ 579*c4762a1bSJed Brown { 580*c4762a1bSJed Brown PetscObject pressure; 581*c4762a1bSJed Brown MatNullSpace nullspacePres; 582*c4762a1bSJed Brown 583*c4762a1bSJed Brown ierr = DMGetField(dm, 1, NULL, &pressure);CHKERRQ(ierr); 584*c4762a1bSJed Brown ierr = MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres);CHKERRQ(ierr); 585*c4762a1bSJed Brown ierr = PetscObjectCompose(pressure, "nullspace", (PetscObject) nullspacePres);CHKERRQ(ierr); 586*c4762a1bSJed Brown ierr = MatNullSpaceDestroy(&nullspacePres);CHKERRQ(ierr); 587*c4762a1bSJed Brown } 588*c4762a1bSJed Brown PetscFunctionReturn(0); 589*c4762a1bSJed Brown } 590*c4762a1bSJed Brown 591*c4762a1bSJed Brown /* Add a vector in the nullspace to make the continuum integral 0. 592*c4762a1bSJed Brown 593*c4762a1bSJed Brown If int(u) = a and int(n) = b, then int(u - a/b n) = a - a/b b = 0 594*c4762a1bSJed Brown */ 595*c4762a1bSJed Brown static PetscErrorCode CorrectDiscretePressure(DM dm, MatNullSpace nullspace, Vec u, AppCtx *user) 596*c4762a1bSJed Brown { 597*c4762a1bSJed Brown PetscDS prob; 598*c4762a1bSJed Brown const Vec *nullvecs; 599*c4762a1bSJed Brown PetscScalar pintd, intc[2], intn[2]; 600*c4762a1bSJed Brown MPI_Comm comm; 601*c4762a1bSJed Brown PetscErrorCode ierr; 602*c4762a1bSJed Brown 603*c4762a1bSJed Brown PetscFunctionBeginUser; 604*c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 605*c4762a1bSJed Brown ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 606*c4762a1bSJed Brown ierr = PetscDSSetObjective(prob, 1, pressure);CHKERRQ(ierr); 607*c4762a1bSJed Brown ierr = MatNullSpaceGetVecs(nullspace, NULL, NULL, &nullvecs);CHKERRQ(ierr); 608*c4762a1bSJed Brown ierr = VecDot(nullvecs[0], u, &pintd);CHKERRQ(ierr); 609*c4762a1bSJed Brown if (PetscAbsScalar(pintd) > 1.0e-10) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Discrete integral of pressure: %g\n", (double) PetscRealPart(pintd)); 610*c4762a1bSJed Brown ierr = DMPlexComputeIntegralFEM(dm, nullvecs[0], intn, user);CHKERRQ(ierr); 611*c4762a1bSJed Brown ierr = DMPlexComputeIntegralFEM(dm, u, intc, user);CHKERRQ(ierr); 612*c4762a1bSJed Brown ierr = VecAXPY(u, -intc[1]/intn[1], nullvecs[0]);CHKERRQ(ierr); 613*c4762a1bSJed Brown ierr = DMPlexComputeIntegralFEM(dm, u, intc, user);CHKERRQ(ierr); 614*c4762a1bSJed Brown if (PetscAbsScalar(intc[1]) > 1.0e-10) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Continuum integral of pressure after correction: %g\n", (double) PetscRealPart(intc[1])); 615*c4762a1bSJed Brown PetscFunctionReturn(0); 616*c4762a1bSJed Brown } 617*c4762a1bSJed Brown 618*c4762a1bSJed Brown static PetscErrorCode SNESConvergenceCorrectPressure(SNES snes, PetscInt it, PetscReal xnorm, PetscReal gnorm, PetscReal f, SNESConvergedReason *reason, void *user) 619*c4762a1bSJed Brown { 620*c4762a1bSJed Brown PetscErrorCode ierr; 621*c4762a1bSJed Brown 622*c4762a1bSJed Brown PetscFunctionBeginUser; 623*c4762a1bSJed Brown ierr = SNESConvergedDefault(snes, it, xnorm, gnorm, f, reason, user);CHKERRQ(ierr); 624*c4762a1bSJed Brown if (*reason > 0) { 625*c4762a1bSJed Brown DM dm; 626*c4762a1bSJed Brown Mat J; 627*c4762a1bSJed Brown Vec u; 628*c4762a1bSJed Brown MatNullSpace nullspace; 629*c4762a1bSJed Brown 630*c4762a1bSJed Brown ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr); 631*c4762a1bSJed Brown ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr); 632*c4762a1bSJed Brown ierr = SNESGetJacobian(snes, &J, NULL, NULL, NULL);CHKERRQ(ierr); 633*c4762a1bSJed Brown ierr = MatGetNullSpace(J, &nullspace);CHKERRQ(ierr); 634*c4762a1bSJed Brown ierr = CorrectDiscretePressure(dm, nullspace, u, (AppCtx *) user);CHKERRQ(ierr); 635*c4762a1bSJed Brown } 636*c4762a1bSJed Brown PetscFunctionReturn(0); 637*c4762a1bSJed Brown } 638*c4762a1bSJed Brown 639*c4762a1bSJed Brown int main(int argc, char **argv) 640*c4762a1bSJed Brown { 641*c4762a1bSJed Brown SNES snes; /* nonlinear solver */ 642*c4762a1bSJed Brown DM dm; /* problem definition */ 643*c4762a1bSJed Brown Vec u, r; /* solution and residual */ 644*c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 645*c4762a1bSJed Brown PetscReal error = 0.0; /* L_2 error in the solution */ 646*c4762a1bSJed Brown PetscReal ferrors[2]; 647*c4762a1bSJed Brown PetscErrorCode ierr; 648*c4762a1bSJed Brown 649*c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 650*c4762a1bSJed Brown ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 651*c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr); 652*c4762a1bSJed Brown ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 653*c4762a1bSJed Brown ierr = SNESSetDM(snes, dm);CHKERRQ(ierr); 654*c4762a1bSJed Brown ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr); 655*c4762a1bSJed Brown 656*c4762a1bSJed Brown ierr = PetscMalloc(2 * sizeof(void (*)(const PetscReal[], PetscScalar *, void *)), &user.exactFuncs);CHKERRQ(ierr); 657*c4762a1bSJed Brown ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr); 658*c4762a1bSJed Brown ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr); 659*c4762a1bSJed Brown 660*c4762a1bSJed Brown ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 661*c4762a1bSJed Brown ierr = VecDuplicate(u, &r);CHKERRQ(ierr); 662*c4762a1bSJed Brown 663*c4762a1bSJed Brown ierr = DMSetNullSpaceConstructor(dm, 2, CreatePressureNullSpace);CHKERRQ(ierr); 664*c4762a1bSJed Brown ierr = SNESSetConvergenceTest(snes, SNESConvergenceCorrectPressure, &user, NULL);CHKERRQ(ierr); 665*c4762a1bSJed Brown 666*c4762a1bSJed Brown ierr = DMPlexSetSNESLocalFEM(dm,&user,&user,&user);CHKERRQ(ierr); 667*c4762a1bSJed Brown 668*c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 669*c4762a1bSJed Brown 670*c4762a1bSJed Brown ierr = DMProjectFunction(dm, 0.0, user.exactFuncs, NULL, INSERT_ALL_VALUES, u);CHKERRQ(ierr); 671*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) u, "Exact Solution");CHKERRQ(ierr); 672*c4762a1bSJed Brown ierr = VecViewFromOptions(u, NULL, "-exact_vec_view");CHKERRQ(ierr); 673*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) u, "Solution");CHKERRQ(ierr); 674*c4762a1bSJed Brown if (user.showInitial) {ierr = DMVecViewLocal(dm, u);CHKERRQ(ierr);} 675*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) u, "Solution");CHKERRQ(ierr); 676*c4762a1bSJed Brown if (user.runType == RUN_FULL) { 677*c4762a1bSJed Brown PetscErrorCode (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void* ctx) = {zero_vector, zero_scalar}; 678*c4762a1bSJed Brown 679*c4762a1bSJed Brown ierr = DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u);CHKERRQ(ierr); 680*c4762a1bSJed Brown if (user.debug) { 681*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n");CHKERRQ(ierr); 682*c4762a1bSJed Brown ierr = VecView(u, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 683*c4762a1bSJed Brown } 684*c4762a1bSJed Brown ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr); 685*c4762a1bSJed Brown ierr = DMComputeL2Diff(dm, 0.0, user.exactFuncs, NULL, u, &error);CHKERRQ(ierr); 686*c4762a1bSJed Brown ierr = DMComputeL2FieldDiff(dm, 0.0, user.exactFuncs, NULL, u, ferrors);CHKERRQ(ierr); 687*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g [%g, %g]\n", (double)error, (double)ferrors[0], (double)ferrors[1]);CHKERRQ(ierr); 688*c4762a1bSJed Brown if (user.showError) { 689*c4762a1bSJed Brown Vec r; 690*c4762a1bSJed Brown 691*c4762a1bSJed Brown ierr = DMGetGlobalVector(dm, &r);CHKERRQ(ierr); 692*c4762a1bSJed Brown ierr = DMProjectFunction(dm, 0.0, user.exactFuncs, NULL, INSERT_ALL_VALUES, r);CHKERRQ(ierr); 693*c4762a1bSJed Brown ierr = VecAXPY(r, -1.0, u);CHKERRQ(ierr); 694*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) r, "Solution Error");CHKERRQ(ierr); 695*c4762a1bSJed Brown ierr = VecViewFromOptions(r, NULL, "-error_vec_view");CHKERRQ(ierr); 696*c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dm, &r);CHKERRQ(ierr); 697*c4762a1bSJed Brown } 698*c4762a1bSJed Brown } else { 699*c4762a1bSJed Brown PetscReal res = 0.0; 700*c4762a1bSJed Brown 701*c4762a1bSJed Brown /* Check discretization error */ 702*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n");CHKERRQ(ierr); 703*c4762a1bSJed Brown ierr = VecView(u, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 704*c4762a1bSJed Brown ierr = DMComputeL2Diff(dm, 0.0, user.exactFuncs, NULL, u, &error);CHKERRQ(ierr); 705*c4762a1bSJed Brown if (error >= 1.0e-11) {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g\n", (double)error);CHKERRQ(ierr);} 706*c4762a1bSJed Brown else {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: < 1.0e-11\n");CHKERRQ(ierr);} 707*c4762a1bSJed Brown /* Check residual */ 708*c4762a1bSJed Brown ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr); 709*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n");CHKERRQ(ierr); 710*c4762a1bSJed Brown ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr); 711*c4762a1bSJed Brown ierr = VecView(r, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 712*c4762a1bSJed Brown ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr); 713*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", (double)res);CHKERRQ(ierr); 714*c4762a1bSJed Brown /* Check Jacobian */ 715*c4762a1bSJed Brown { 716*c4762a1bSJed Brown Mat J, M; 717*c4762a1bSJed Brown MatNullSpace nullspace; 718*c4762a1bSJed Brown Vec b; 719*c4762a1bSJed Brown PetscBool isNull; 720*c4762a1bSJed Brown 721*c4762a1bSJed Brown ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 722*c4762a1bSJed Brown ierr = SNESGetJacobian(snes, &J, &M, NULL, NULL);CHKERRQ(ierr); 723*c4762a1bSJed Brown ierr = SNESComputeJacobian(snes, u, J, M);CHKERRQ(ierr); 724*c4762a1bSJed Brown ierr = MatGetNullSpace(J, &nullspace);CHKERRQ(ierr); 725*c4762a1bSJed Brown ierr = MatNullSpaceTest(nullspace, J, &isNull);CHKERRQ(ierr); 726*c4762a1bSJed Brown ierr = VecDuplicate(u, &b);CHKERRQ(ierr); 727*c4762a1bSJed Brown ierr = VecSet(r, 0.0);CHKERRQ(ierr); 728*c4762a1bSJed Brown ierr = SNESComputeFunction(snes, r, b);CHKERRQ(ierr); 729*c4762a1bSJed Brown ierr = MatMult(J, u, r);CHKERRQ(ierr); 730*c4762a1bSJed Brown ierr = VecAXPY(r, 1.0, b);CHKERRQ(ierr); 731*c4762a1bSJed Brown ierr = VecDestroy(&b);CHKERRQ(ierr); 732*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Au - b = Au + F(0)\n");CHKERRQ(ierr); 733*c4762a1bSJed Brown ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr); 734*c4762a1bSJed Brown ierr = VecView(r, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 735*c4762a1bSJed Brown ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr); 736*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", (double)res);CHKERRQ(ierr); 737*c4762a1bSJed Brown } 738*c4762a1bSJed Brown } 739*c4762a1bSJed Brown ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr); 740*c4762a1bSJed Brown 741*c4762a1bSJed Brown ierr = VecDestroy(&u);CHKERRQ(ierr); 742*c4762a1bSJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 743*c4762a1bSJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 744*c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 745*c4762a1bSJed Brown ierr = PetscFree(user.exactFuncs);CHKERRQ(ierr); 746*c4762a1bSJed Brown ierr = PetscFinalize(); 747*c4762a1bSJed Brown return ierr; 748*c4762a1bSJed Brown } 749*c4762a1bSJed Brown 750*c4762a1bSJed Brown /*TEST 751*c4762a1bSJed Brown 752*c4762a1bSJed Brown # 2D serial P1 tests 0-3 753*c4762a1bSJed Brown test: 754*c4762a1bSJed Brown suffix: 0 755*c4762a1bSJed Brown requires: triangle 756*c4762a1bSJed Brown args: -run_type test -refinement_limit 0.0 -bc_type dirichlet -interpolate 0 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -show_initial -dm_plex_print_fem 1 757*c4762a1bSJed Brown test: 758*c4762a1bSJed Brown suffix: 1 759*c4762a1bSJed Brown requires: triangle 760*c4762a1bSJed Brown args: -run_type test -refinement_limit 0.0 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -show_initial -dm_plex_print_fem 1 761*c4762a1bSJed Brown test: 762*c4762a1bSJed Brown suffix: 2 763*c4762a1bSJed Brown requires: triangle 764*c4762a1bSJed Brown args: -run_type test -refinement_limit 0.0625 -bc_type dirichlet -interpolate 0 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -show_initial -dm_plex_print_fem 1 765*c4762a1bSJed Brown test: 766*c4762a1bSJed Brown suffix: 3 767*c4762a1bSJed Brown requires: triangle 768*c4762a1bSJed Brown args: -run_type test -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -show_initial -dm_plex_print_fem 1 769*c4762a1bSJed Brown # 2D serial P2 tests 4-5 770*c4762a1bSJed Brown test: 771*c4762a1bSJed Brown suffix: 4 772*c4762a1bSJed Brown requires: triangle 773*c4762a1bSJed Brown args: -run_type test -refinement_limit 0.0 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -show_initial -dm_plex_print_fem 1 774*c4762a1bSJed Brown test: 775*c4762a1bSJed Brown suffix: 5 776*c4762a1bSJed Brown requires: triangle 777*c4762a1bSJed Brown args: -run_type test -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -show_initial -dm_plex_print_fem 1 778*c4762a1bSJed Brown # 2D serial P3 tests 779*c4762a1bSJed Brown test: 780*c4762a1bSJed Brown suffix: 2d_p3_0 781*c4762a1bSJed Brown requires: triangle 782*c4762a1bSJed Brown args: -run_type test -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 3 -pres_petscspace_degree 2 783*c4762a1bSJed Brown test: 784*c4762a1bSJed Brown suffix: 2d_p3_1 785*c4762a1bSJed Brown requires: triangle !single 786*c4762a1bSJed Brown args: -run_type full -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 3 -pres_petscspace_degree 2 787*c4762a1bSJed Brown # Parallel tests 6-17 788*c4762a1bSJed Brown test: 789*c4762a1bSJed Brown suffix: 6 790*c4762a1bSJed Brown requires: triangle 791*c4762a1bSJed Brown nsize: 2 792*c4762a1bSJed Brown args: -run_type test -refinement_limit 0.0 -test_partition -bc_type dirichlet -interpolate 0 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -dm_plex_print_fem 1 793*c4762a1bSJed Brown test: 794*c4762a1bSJed Brown suffix: 7 795*c4762a1bSJed Brown requires: triangle 796*c4762a1bSJed Brown nsize: 3 797*c4762a1bSJed Brown args: -run_type test -refinement_limit 0.0 -test_partition -bc_type dirichlet -interpolate 0 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -dm_plex_print_fem 1 798*c4762a1bSJed Brown test: 799*c4762a1bSJed Brown suffix: 8 800*c4762a1bSJed Brown requires: triangle 801*c4762a1bSJed Brown nsize: 5 802*c4762a1bSJed Brown args: -run_type test -refinement_limit 0.0 -test_partition -bc_type dirichlet -interpolate 0 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -dm_plex_print_fem 1 803*c4762a1bSJed Brown test: 804*c4762a1bSJed Brown suffix: 9 805*c4762a1bSJed Brown requires: triangle 806*c4762a1bSJed Brown nsize: 2 807*c4762a1bSJed Brown args: -run_type test -refinement_limit 0.0 -test_partition -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -dm_plex_print_fem 1 808*c4762a1bSJed Brown test: 809*c4762a1bSJed Brown suffix: 10 810*c4762a1bSJed Brown requires: triangle 811*c4762a1bSJed Brown nsize: 3 812*c4762a1bSJed Brown args: -run_type test -refinement_limit 0.0 -test_partition -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -dm_plex_print_fem 1 813*c4762a1bSJed Brown test: 814*c4762a1bSJed Brown suffix: 11 815*c4762a1bSJed Brown requires: triangle 816*c4762a1bSJed Brown nsize: 5 817*c4762a1bSJed Brown args: -run_type test -refinement_limit 0.0 -test_partition -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -dm_plex_print_fem 1 818*c4762a1bSJed Brown test: 819*c4762a1bSJed Brown suffix: 12 820*c4762a1bSJed Brown requires: triangle 821*c4762a1bSJed Brown nsize: 2 822*c4762a1bSJed Brown args: -run_type test -refinement_limit 0.0625 -test_partition -bc_type dirichlet -interpolate 0 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -dm_plex_print_fem 1 823*c4762a1bSJed Brown test: 824*c4762a1bSJed Brown suffix: 13 825*c4762a1bSJed Brown requires: triangle 826*c4762a1bSJed Brown nsize: 3 827*c4762a1bSJed Brown args: -run_type test -refinement_limit 0.0625 -test_partition -bc_type dirichlet -interpolate 0 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -dm_plex_print_fem 1 828*c4762a1bSJed Brown test: 829*c4762a1bSJed Brown suffix: 14 830*c4762a1bSJed Brown requires: triangle 831*c4762a1bSJed Brown nsize: 5 832*c4762a1bSJed Brown args: -run_type test -refinement_limit 0.0625 -test_partition -bc_type dirichlet -interpolate 0 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -dm_plex_print_fem 1 833*c4762a1bSJed Brown test: 834*c4762a1bSJed Brown suffix: 15 835*c4762a1bSJed Brown requires: triangle 836*c4762a1bSJed Brown nsize: 2 837*c4762a1bSJed Brown args: -run_type test -refinement_limit 0.0625 -test_partition -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -dm_plex_print_fem 1 838*c4762a1bSJed Brown test: 839*c4762a1bSJed Brown suffix: 16 840*c4762a1bSJed Brown requires: triangle 841*c4762a1bSJed Brown nsize: 3 842*c4762a1bSJed Brown args: -run_type test -refinement_limit 0.0625 -test_partition -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -dm_plex_print_fem 1 843*c4762a1bSJed Brown test: 844*c4762a1bSJed Brown suffix: 17 845*c4762a1bSJed Brown requires: triangle 846*c4762a1bSJed Brown nsize: 5 847*c4762a1bSJed Brown args: -run_type test -refinement_limit 0.0625 -test_partition -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -dm_plex_print_fem 1 848*c4762a1bSJed Brown # 3D serial P1 tests 43-46 849*c4762a1bSJed Brown test: 850*c4762a1bSJed Brown suffix: 43 851*c4762a1bSJed Brown requires: ctetgen 852*c4762a1bSJed Brown args: -run_type test -dim 3 -refinement_limit 0.0 -bc_type dirichlet -interpolate 0 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -show_initial -dm_plex_print_fem 1 853*c4762a1bSJed Brown test: 854*c4762a1bSJed Brown suffix: 44 855*c4762a1bSJed Brown requires: ctetgen 856*c4762a1bSJed Brown args: -run_type test -dim 3 -refinement_limit 0.0 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -show_initial -dm_plex_print_fem 1 857*c4762a1bSJed Brown test: 858*c4762a1bSJed Brown suffix: 45 859*c4762a1bSJed Brown requires: ctetgen 860*c4762a1bSJed Brown args: -run_type test -dim 3 -refinement_limit 0.0125 -bc_type dirichlet -interpolate 0 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -show_initial -dm_plex_print_fem 1 861*c4762a1bSJed Brown test: 862*c4762a1bSJed Brown suffix: 46 863*c4762a1bSJed Brown requires: ctetgen 864*c4762a1bSJed Brown args: -run_type test -dim 3 -refinement_limit 0.0125 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -show_initial -dm_plex_print_fem 1 865*c4762a1bSJed Brown # Full solutions 18-29 866*c4762a1bSJed Brown test: 867*c4762a1bSJed Brown suffix: 18 868*c4762a1bSJed Brown requires: triangle !single 869*c4762a1bSJed Brown filter: sed -e "s/total number of linear solver iterations=11/total number of linear solver iterations=12/g" 870*c4762a1bSJed Brown args: -run_type full -refinement_limit 0.0625 -bc_type dirichlet -interpolate 0 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view 871*c4762a1bSJed Brown test: 872*c4762a1bSJed Brown suffix: 19 873*c4762a1bSJed Brown requires: triangle !single 874*c4762a1bSJed Brown nsize: 2 875*c4762a1bSJed Brown filter: sed -e "s/total number of linear solver iterations=11/total number of linear solver iterations=12/g" 876*c4762a1bSJed Brown args: -run_type full -petscpartitioner_type simple -refinement_limit 0.0625 -bc_type dirichlet -interpolate 0 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view 877*c4762a1bSJed Brown test: 878*c4762a1bSJed Brown suffix: 20 879*c4762a1bSJed Brown requires: triangle !single 880*c4762a1bSJed Brown filter: sed -e "s/total number of linear solver iterations=11/total number of linear solver iterations=12/g" 881*c4762a1bSJed Brown nsize: 3 882*c4762a1bSJed Brown args: -run_type full -petscpartitioner_type simple -refinement_limit 0.0625 -bc_type dirichlet -interpolate 0 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view 883*c4762a1bSJed Brown test: 884*c4762a1bSJed Brown suffix: 20_parmetis 885*c4762a1bSJed Brown requires: parmetis triangle !single 886*c4762a1bSJed Brown filter: sed -e "s/total number of linear solver iterations=11/total number of linear solver iterations=12/g" 887*c4762a1bSJed Brown nsize: 3 888*c4762a1bSJed Brown args: -run_type full -petscpartitioner_type parmetis -refinement_limit 0.0625 -bc_type dirichlet -interpolate 0 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view 889*c4762a1bSJed Brown test: 890*c4762a1bSJed Brown suffix: 21 891*c4762a1bSJed Brown requires: triangle !single 892*c4762a1bSJed Brown filter: sed -e "s/total number of linear solver iterations=11/total number of linear solver iterations=12/g" 893*c4762a1bSJed Brown nsize: 5 894*c4762a1bSJed Brown args: -run_type full -petscpartitioner_type simple -refinement_limit 0.0625 -bc_type dirichlet -interpolate 0 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view 895*c4762a1bSJed Brown test: 896*c4762a1bSJed Brown suffix: 22 897*c4762a1bSJed Brown requires: triangle !single 898*c4762a1bSJed Brown filter: sed -e "s/total number of linear solver iterations=11/total number of linear solver iterations=12/g" 899*c4762a1bSJed Brown args: -run_type full -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view 900*c4762a1bSJed Brown test: 901*c4762a1bSJed Brown suffix: 23 902*c4762a1bSJed Brown requires: triangle !single 903*c4762a1bSJed Brown nsize: 2 904*c4762a1bSJed Brown filter: sed -e "s/total number of linear solver iterations=11/total number of linear solver iterations=12/g" 905*c4762a1bSJed Brown args: -run_type full -petscpartitioner_type simple -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view 906*c4762a1bSJed Brown test: 907*c4762a1bSJed Brown suffix: 24 908*c4762a1bSJed Brown requires: triangle !single 909*c4762a1bSJed Brown nsize: 3 910*c4762a1bSJed Brown filter: sed -e "s/total number of linear solver iterations=11/total number of linear solver iterations=12/g" 911*c4762a1bSJed Brown args: -run_type full -petscpartitioner_type simple -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view 912*c4762a1bSJed Brown test: 913*c4762a1bSJed Brown suffix: 25 914*c4762a1bSJed Brown requires: triangle !single 915*c4762a1bSJed Brown filter: sed -e "s/total number of linear solver iterations=11/total number of linear solver iterations=12/g" 916*c4762a1bSJed Brown nsize: 5 917*c4762a1bSJed Brown args: -run_type full -petscpartitioner_type simple -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 1 -pres_petscspace_degree 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view 918*c4762a1bSJed Brown test: 919*c4762a1bSJed Brown suffix: 26 920*c4762a1bSJed Brown requires: triangle !single 921*c4762a1bSJed Brown args: -run_type full -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view 922*c4762a1bSJed Brown test: 923*c4762a1bSJed Brown suffix: 27 924*c4762a1bSJed Brown requires: triangle !single 925*c4762a1bSJed Brown nsize: 2 926*c4762a1bSJed Brown args: -run_type full -petscpartitioner_type simple -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view 927*c4762a1bSJed Brown test: 928*c4762a1bSJed Brown suffix: 28 929*c4762a1bSJed Brown requires: triangle !single 930*c4762a1bSJed Brown nsize: 3 931*c4762a1bSJed Brown args: -run_type full -petscpartitioner_type simple -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view 932*c4762a1bSJed Brown test: 933*c4762a1bSJed Brown suffix: 29 934*c4762a1bSJed Brown requires: triangle !single 935*c4762a1bSJed Brown nsize: 5 936*c4762a1bSJed Brown args: -run_type full -petscpartitioner_type simple -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view 937*c4762a1bSJed Brown # Full solutions with quads 938*c4762a1bSJed Brown # FULL Schur with LU/Jacobi 939*c4762a1bSJed Brown test: 940*c4762a1bSJed Brown suffix: quad_q2q1_full 941*c4762a1bSJed Brown requires: !single 942*c4762a1bSJed Brown args: -run_type full -simplex 0 -refinement_limit 0.00625 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_velocity_ksp_type gmres -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view 943*c4762a1bSJed Brown test: 944*c4762a1bSJed Brown suffix: quad_q2p1_full 945*c4762a1bSJed Brown requires: !single 946*c4762a1bSJed Brown args: -run_type full -simplex 0 -refinement_limit 0.00625 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -pres_petscspace_poly_tensor 0 -pres_petscdualspace_lagrange_continuity 0 -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_velocity_ksp_type gmres -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view 947*c4762a1bSJed Brown # Stokes preconditioners 30-36 948*c4762a1bSJed Brown # Jacobi 949*c4762a1bSJed Brown test: 950*c4762a1bSJed Brown suffix: 30 951*c4762a1bSJed Brown requires: triangle !single 952*c4762a1bSJed Brown filter: sed -e "s/total number of linear solver iterations=756/total number of linear solver iterations=757/g" -e "s/total number of linear solver iterations=758/total number of linear solver iterations=757/g" 953*c4762a1bSJed Brown args: -run_type full -refinement_limit 0.00625 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -ksp_gmres_restart 100 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view 954*c4762a1bSJed Brown # Block diagonal \begin{pmatrix} A & 0 \\ 0 & I \end{pmatrix} 955*c4762a1bSJed Brown test: 956*c4762a1bSJed Brown suffix: 31 957*c4762a1bSJed Brown requires: triangle !single 958*c4762a1bSJed Brown args: -run_type full -refinement_limit 0.00625 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-4 -pc_type fieldsplit -pc_fieldsplit_type additive -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view 959*c4762a1bSJed Brown # Block triangular \begin{pmatrix} A & B \\ 0 & I \end{pmatrix} 960*c4762a1bSJed Brown test: 961*c4762a1bSJed Brown suffix: 32 962*c4762a1bSJed Brown requires: triangle !single 963*c4762a1bSJed Brown args: -run_type full -refinement_limit 0.00625 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type fieldsplit -pc_fieldsplit_type multiplicative -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view 964*c4762a1bSJed Brown # Diagonal Schur complement \begin{pmatrix} A & 0 \\ 0 & S \end{pmatrix} 965*c4762a1bSJed Brown test: 966*c4762a1bSJed Brown suffix: 33 967*c4762a1bSJed Brown requires: triangle !single 968*c4762a1bSJed Brown args: -run_type full -refinement_limit 0.00625 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type diag -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_velocity_ksp_type gmres -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view 969*c4762a1bSJed Brown # Upper triangular Schur complement \begin{pmatrix} A & B \\ 0 & S \end{pmatrix} 970*c4762a1bSJed Brown test: 971*c4762a1bSJed Brown suffix: 34 972*c4762a1bSJed Brown requires: triangle !single 973*c4762a1bSJed Brown args: -run_type full -refinement_limit 0.00625 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_velocity_ksp_type gmres -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view 974*c4762a1bSJed Brown # Lower triangular Schur complement \begin{pmatrix} A & B \\ 0 & S \end{pmatrix} 975*c4762a1bSJed Brown test: 976*c4762a1bSJed Brown suffix: 35 977*c4762a1bSJed Brown requires: triangle !single 978*c4762a1bSJed Brown args: -run_type full -refinement_limit 0.00625 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type lower -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_velocity_ksp_type gmres -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view 979*c4762a1bSJed Brown # Full Schur complement \begin{pmatrix} I & 0 \\ B^T A^{-1} & I \end{pmatrix} \begin{pmatrix} A & 0 \\ 0 & S \end{pmatrix} \begin{pmatrix} I & A^{-1} B \\ 0 & I \end{pmatrix} 980*c4762a1bSJed Brown test: 981*c4762a1bSJed Brown suffix: 36 982*c4762a1bSJed Brown requires: triangle !single 983*c4762a1bSJed Brown args: -run_type full -refinement_limit 0.00625 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_velocity_ksp_type gmres -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view 984*c4762a1bSJed Brown # SIMPLE \begin{pmatrix} I & 0 \\ B^T A^{-1} & I \end{pmatrix} \begin{pmatrix} A & 0 \\ 0 & B^T diag(A)^{-1} B \end{pmatrix} \begin{pmatrix} I & diag(A)^{-1} B \\ 0 & I \end{pmatrix} 985*c4762a1bSJed Brown test: 986*c4762a1bSJed Brown suffix: pc_simple 987*c4762a1bSJed Brown requires: triangle !single 988*c4762a1bSJed Brown args: -run_type full -refinement_limit 0.00625 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_velocity_ksp_type gmres -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi -fieldsplit_pressure_inner_ksp_type preonly -fieldsplit_pressure_inner_pc_type jacobi -fieldsplit_pressure_upper_ksp_type preonly -fieldsplit_pressure_upper_pc_type jacobi -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view 989*c4762a1bSJed Brown # SIMPLEC \begin{pmatrix} I & 0 \\ B^T A^{-1} & I \end{pmatrix} \begin{pmatrix} A & 0 \\ 0 & B^T rowsum(A)^{-1} B \end{pmatrix} \begin{pmatrix} I & rowsum(A)^{-1} B \\ 0 & I \end{pmatrix} 990*c4762a1bSJed Brown test: 991*c4762a1bSJed Brown suffix: pc_simplec 992*c4762a1bSJed Brown requires: triangle 993*c4762a1bSJed Brown args: -run_type full -dm_refine 3 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -ksp_type fgmres -ksp_max_it 5 -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_ksp_max_it 10 -fieldsplit_velocity_ksp_type gmres -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi -fieldsplit_pressure_inner_ksp_type preonly -fieldsplit_pressure_inner_pc_type jacobi -fieldsplit_pressure_inner_pc_jacobi_type rowsum -fieldsplit_pressure_upper_ksp_type preonly -fieldsplit_pressure_upper_pc_type jacobi -fieldsplit_pressure_upper_pc_jacobi_type rowsum -snes_converged_reason -ksp_converged_reason -snes_view 994*c4762a1bSJed Brown # FETI-DP solvers (these solvers are quite inefficient, they are here to exercise the code) 995*c4762a1bSJed Brown test: 996*c4762a1bSJed Brown suffix: fetidp_2d_tri 997*c4762a1bSJed Brown requires: triangle mumps 998*c4762a1bSJed Brown filter: grep -v "variant HERMITIAN" | sed -e "s/linear solver iterations=41/linear solver iterations=42/g" 999*c4762a1bSJed Brown nsize: 5 1000*c4762a1bSJed Brown args: -run_type full -dm_refine 2 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_view -snes_error_if_not_converged -dm_mat_type is -ksp_type fetidp -ksp_rtol 1.0e-8 -ksp_fetidp_saddlepoint -fetidp_ksp_type cg -fetidp_fieldsplit_p_ksp_max_it 1 -fetidp_fieldsplit_p_ksp_type richardson -fetidp_fieldsplit_p_ksp_richardson_scale 200 -fetidp_fieldsplit_p_pc_type none -ksp_fetidp_saddlepoint_flip 1 -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type mumps -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type mumps -petscpartitioner_type simple -fetidp_fieldsplit_lag_ksp_type preonly 1001*c4762a1bSJed Brown test: 1002*c4762a1bSJed Brown suffix: fetidp_3d_tet_smumps 1003*c4762a1bSJed Brown output_file: output/ex62_fetidp_3d_tet.out 1004*c4762a1bSJed Brown requires: ctetgen suitesparse mumps 1005*c4762a1bSJed Brown filter: grep -v "variant HERMITIAN" | sed -e "s/linear solver iterations=10[0-9]/linear solver iterations=100/g" | sed -e "s/linear solver iterations=9[0-9]/linear solver iterations=100/g" 1006*c4762a1bSJed Brown nsize: 5 1007*c4762a1bSJed Brown args: -run_type full -dm_refine 2 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_view -snes_error_if_not_converged -dm_mat_type is -ksp_type fetidp -ksp_rtol 1.0e-8 -ksp_fetidp_saddlepoint -fetidp_ksp_type cg -fetidp_fieldsplit_p_ksp_max_it 1 -fetidp_fieldsplit_p_ksp_type richardson -fetidp_fieldsplit_p_ksp_richardson_scale 1000 -fetidp_fieldsplit_p_pc_type none -ksp_fetidp_saddlepoint_flip 1 -fetidp_bddc_pc_bddc_use_deluxe_scaling -fetidp_bddc_pc_bddc_benign_trick -fetidp_bddc_pc_bddc_deluxe_singlemat -dim 3 -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_factor_mat_solver_type petsc -fetidp_harmonic_pc_type cholesky -fetidp_bddelta_pc_factor_mat_solver_type umfpack -fetidp_fieldsplit_lag_ksp_type preonly -test_partition -fetidp_bddc_sub_schurs_mat_solver_type mumps 1008*c4762a1bSJed Brown test: 1009*c4762a1bSJed Brown suffix: fetidp_3d_tet_spetsc 1010*c4762a1bSJed Brown requires: ctetgen suitesparse 1011*c4762a1bSJed Brown filter: grep -v "variant HERMITIAN" | sed -e "s/linear solver iterations=10[0-9]/linear solver iterations=100/g" | sed -e "s/linear solver iterations=9[0-9]/linear solver iterations=100/g" 1012*c4762a1bSJed Brown nsize: 5 1013*c4762a1bSJed Brown args: -run_type full -dm_refine 2 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_view -snes_error_if_not_converged -dm_mat_type is -ksp_type fetidp -ksp_rtol 1.0e-8 -ksp_fetidp_saddlepoint -fetidp_ksp_type cg -fetidp_fieldsplit_p_ksp_max_it 1 -fetidp_fieldsplit_p_ksp_type richardson -fetidp_fieldsplit_p_ksp_richardson_scale 1000 -fetidp_fieldsplit_p_pc_type none -ksp_fetidp_saddlepoint_flip 1 -fetidp_bddc_pc_bddc_use_deluxe_scaling -fetidp_bddc_pc_bddc_benign_trick -fetidp_bddc_pc_bddc_deluxe_singlemat -dim 3 -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_factor_mat_solver_type petsc -fetidp_harmonic_pc_type cholesky -fetidp_bddelta_pc_factor_mat_solver_type umfpack -fetidp_fieldsplit_lag_ksp_type preonly -test_partition -fetidp_bddc_sub_schurs_mat_solver_type petsc -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type umfpack -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type umfpack 1014*c4762a1bSJed Brown test: 1015*c4762a1bSJed Brown suffix: fetidp_2d_quad 1016*c4762a1bSJed Brown requires: mumps double 1017*c4762a1bSJed Brown filter: grep -v "variant HERMITIAN" 1018*c4762a1bSJed Brown nsize: 5 1019*c4762a1bSJed Brown args: -run_type full -dm_refine 2 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_view -snes_error_if_not_converged -dm_mat_type is -ksp_type fetidp -ksp_rtol 1.0e-8 -ksp_fetidp_saddlepoint -fetidp_ksp_type cg -fetidp_fieldsplit_p_ksp_max_it 1 -fetidp_fieldsplit_p_ksp_type richardson -fetidp_fieldsplit_p_ksp_richardson_scale 200 -fetidp_fieldsplit_p_pc_type none -ksp_fetidp_saddlepoint_flip 1 -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type mumps -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type mumps -simplex 0 -petscpartitioner_type simple -fetidp_fieldsplit_lag_ksp_type preonly 1020*c4762a1bSJed Brown test: 1021*c4762a1bSJed Brown suffix: fetidp_3d_hex 1022*c4762a1bSJed Brown requires: suitesparse 1023*c4762a1bSJed Brown filter: grep -v "variant HERMITIAN" | sed -e "s/linear solver iterations=7[0-9]/linear solver iterations=71/g" 1024*c4762a1bSJed Brown nsize: 5 1025*c4762a1bSJed Brown args: -run_type full -dm_refine 1 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_view -snes_error_if_not_converged -dm_mat_type is -ksp_type fetidp -ksp_rtol 1.0e-8 -ksp_fetidp_saddlepoint -fetidp_ksp_type cg -fetidp_fieldsplit_p_ksp_max_it 1 -fetidp_fieldsplit_p_ksp_type richardson -fetidp_fieldsplit_p_ksp_richardson_scale 2000 -fetidp_fieldsplit_p_pc_type none -ksp_fetidp_saddlepoint_flip 1 -dim 3 -simplex 0 -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_factor_mat_solver_type petsc -fetidp_harmonic_pc_type cholesky -petscpartitioner_type simple -fetidp_fieldsplit_lag_ksp_type preonly -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type umfpack -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type umfpack 1026*c4762a1bSJed Brown # Convergence 1027*c4762a1bSJed Brown test: 1028*c4762a1bSJed Brown suffix: 2d_quad_q1_p0_conv 1029*c4762a1bSJed Brown requires: !single 1030*c4762a1bSJed Brown args: -run_type full -bc_type dirichlet -simplex 0 -interpolate 1 -dm_refine 0 -vel_petscspace_degree 1 -pres_petscspace_degree 0 \ 1031*c4762a1bSJed Brown -snes_convergence_estimate -convest_num_refine 3 -snes_error_if_not_converged \ 1032*c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 1033*c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 1034*c4762a1bSJed Brown -fieldsplit_velocity_pc_type lu \ 1035*c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 1036*c4762a1bSJed Brown test: 1037*c4762a1bSJed Brown suffix: 2d_tri_p2_p1_conv 1038*c4762a1bSJed Brown requires: triangle !single 1039*c4762a1bSJed Brown args: -run_type full -sol_type cubic -bc_type dirichlet -interpolate 1 -dm_refine 0 \ 1040*c4762a1bSJed Brown -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 1041*c4762a1bSJed Brown -snes_convergence_estimate -convest_num_refine 3 -snes_error_if_not_converged \ 1042*c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 1043*c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 1044*c4762a1bSJed Brown -fieldsplit_velocity_pc_type lu \ 1045*c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 1046*c4762a1bSJed Brown test: 1047*c4762a1bSJed Brown suffix: 2d_quad_q2_q1_conv 1048*c4762a1bSJed Brown requires: !single 1049*c4762a1bSJed Brown args: -run_type full -sol_type cubic -bc_type dirichlet -simplex 0 -interpolate 1 -dm_refine 0 \ 1050*c4762a1bSJed Brown -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 1051*c4762a1bSJed Brown -snes_convergence_estimate -convest_num_refine 3 -snes_error_if_not_converged \ 1052*c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 1053*c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 1054*c4762a1bSJed Brown -fieldsplit_velocity_pc_type lu \ 1055*c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 1056*c4762a1bSJed Brown test: 1057*c4762a1bSJed Brown suffix: 2d_quad_q2_p1_conv 1058*c4762a1bSJed Brown requires: !single 1059*c4762a1bSJed Brown args: -run_type full -sol_type cubic -bc_type dirichlet -simplex 0 -interpolate 1 -dm_refine 0 \ 1060*c4762a1bSJed Brown -vel_petscspace_degree 2 -pres_petscspace_degree 1 -pres_petscspace_poly_tensor 0 -pres_petscdualspace_lagrange_continuity 0 \ 1061*c4762a1bSJed Brown -snes_convergence_estimate -convest_num_refine 3 -snes_error_if_not_converged \ 1062*c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 1063*c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 1064*c4762a1bSJed Brown -fieldsplit_velocity_pc_type lu \ 1065*c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 1066*c4762a1bSJed Brown # GMG solver 1067*c4762a1bSJed Brown test: 1068*c4762a1bSJed Brown suffix: 2d_tri_p2_p1_gmg_vcycle 1069*c4762a1bSJed Brown requires: triangle 1070*c4762a1bSJed Brown args: -run_type full -sol_type cubic -bc_type dirichlet -interpolate 1 -cells 2,2 -dm_refine_hierarchy 1 \ 1071*c4762a1bSJed Brown -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 1072*c4762a1bSJed Brown -snes_convergence_estimate -convest_num_refine 1 -snes_error_if_not_converged \ 1073*c4762a1bSJed Brown -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 1074*c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 1075*c4762a1bSJed Brown -fieldsplit_velocity_pc_type mg \ 1076*c4762a1bSJed Brown -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 1077*c4762a1bSJed Brown # Vanka solver 1078*c4762a1bSJed Brown test: 1079*c4762a1bSJed Brown suffix: 2d_quad_q1_p0_vanka_add 1080*c4762a1bSJed Brown requires: double !complex 1081*c4762a1bSJed Brown filter: sed -e "s/linear solver iterations=[0-9][0-9]*""/linear solver iterations=49/g" -e "s/Linear solve converged due to CONVERGED_RTOL iterations [0-9][0-9]*""/Linear solve converged due to CONVERGED_RTOL iterations 49/g" 1082*c4762a1bSJed Brown args: -run_type full -bc_type dirichlet -simplex 0 -dm_refine 1 -interpolate 1 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \ 1083*c4762a1bSJed Brown -snes_rtol 1.0e-4 -snes_error_if_not_converged -snes_view -snes_monitor -snes_converged_reason \ 1084*c4762a1bSJed Brown -ksp_type gmres -ksp_rtol 1.0e-5 -ksp_error_if_not_converged -ksp_converged_reason \ 1085*c4762a1bSJed Brown -pc_type patch -pc_patch_partition_of_unity 0 -pc_patch_construct_codim 0 -pc_patch_construct_type vanka \ 1086*c4762a1bSJed Brown -sub_ksp_type preonly -sub_pc_type lu 1087*c4762a1bSJed Brown # Vanka solver, forming dense inverses on patches 1088*c4762a1bSJed Brown test: 1089*c4762a1bSJed Brown suffix: 2d_quad_q1_p0_vanka_add_dense_inverse 1090*c4762a1bSJed Brown requires: double !complex 1091*c4762a1bSJed Brown filter: sed -e "s/linear solver iterations=[0-9][0-9]*""/linear solver iterations=49/g" -e "s/Linear solve converged due to CONVERGED_RTOL iterations [0-9][0-9]*""/Linear solve converged due to CONVERGED_RTOL iterations 49/g" 1092*c4762a1bSJed Brown args: -run_type full -bc_type dirichlet -simplex 0 -dm_refine 1 -interpolate 1 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \ 1093*c4762a1bSJed Brown -snes_rtol 1.0e-4 -snes_error_if_not_converged -snes_view -snes_monitor -snes_converged_reason \ 1094*c4762a1bSJed Brown -ksp_type gmres -ksp_rtol 1.0e-5 -ksp_error_if_not_converged -ksp_converged_reason \ 1095*c4762a1bSJed Brown -pc_type patch -pc_patch_partition_of_unity 0 -pc_patch_construct_codim 0 -pc_patch_construct_type vanka \ 1096*c4762a1bSJed Brown -pc_patch_dense_inverse -pc_patch_sub_mat_type seqdense 1097*c4762a1bSJed Brown test: 1098*c4762a1bSJed Brown suffix: 2d_quad_q1_p0_vanka_add_unity 1099*c4762a1bSJed Brown requires: double !complex 1100*c4762a1bSJed Brown filter: sed -e "s/linear solver iterations=[0-9][0-9]*""/linear solver iterations=45/g" -e "s/Linear solve converged due to CONVERGED_RTOL iterations [0-9][0-9]*""/Linear solve converged due to CONVERGED_RTOL iterations 45/g" 1101*c4762a1bSJed Brown args: -run_type full -bc_type dirichlet -simplex 0 -dm_refine 1 -interpolate 1 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \ 1102*c4762a1bSJed Brown -snes_rtol 1.0e-4 -snes_error_if_not_converged -snes_view -snes_monitor -snes_converged_reason \ 1103*c4762a1bSJed Brown -ksp_type gmres -ksp_rtol 1.0e-5 -ksp_error_if_not_converged -ksp_converged_reason \ 1104*c4762a1bSJed Brown -pc_type patch -pc_patch_partition_of_unity 1 -pc_patch_construct_codim 0 -pc_patch_construct_type vanka \ 1105*c4762a1bSJed Brown -sub_ksp_type preonly -sub_pc_type lu 1106*c4762a1bSJed Brown test: 1107*c4762a1bSJed Brown suffix: 2d_quad_q2_q1_vanka_add 1108*c4762a1bSJed Brown requires: double !complex 1109*c4762a1bSJed Brown filter: sed -e "s/linear solver iterations=[0-9][0-9][0-9]*""/linear solver iterations=489/g" 1110*c4762a1bSJed Brown args: -run_type full -bc_type dirichlet -simplex 0 -dm_refine 0 -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \ 1111*c4762a1bSJed Brown -snes_rtol 1.0e-4 -snes_error_if_not_converged -snes_view -snes_monitor -snes_converged_reason \ 1112*c4762a1bSJed Brown -ksp_type gmres -ksp_rtol 1.0e-5 -ksp_error_if_not_converged \ 1113*c4762a1bSJed Brown -pc_type patch -pc_patch_partition_of_unity 0 -pc_patch_construct_dim 0 -pc_patch_construct_type vanka \ 1114*c4762a1bSJed Brown -sub_ksp_type preonly -sub_pc_type lu 1115*c4762a1bSJed Brown test: 1116*c4762a1bSJed Brown suffix: 2d_quad_q2_q1_vanka_add_unity 1117*c4762a1bSJed Brown requires: double !complex 1118*c4762a1bSJed Brown filter: sed -e "s/linear solver iterations=[0-9][0-9][0-9]*""/linear solver iterations=795/g" 1119*c4762a1bSJed Brown args: -run_type full -bc_type dirichlet -simplex 0 -dm_refine 0 -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \ 1120*c4762a1bSJed Brown -snes_rtol 1.0e-4 -snes_error_if_not_converged -snes_view -snes_monitor -snes_converged_reason \ 1121*c4762a1bSJed Brown -ksp_type gmres -ksp_rtol 1.0e-5 -ksp_error_if_not_converged \ 1122*c4762a1bSJed Brown -pc_type patch -pc_patch_partition_of_unity 1 -pc_patch_construct_dim 0 -pc_patch_construct_type vanka \ 1123*c4762a1bSJed Brown -sub_ksp_type preonly -sub_pc_type lu 1124*c4762a1bSJed Brown # Vanka smoother 1125*c4762a1bSJed Brown test: 1126*c4762a1bSJed Brown suffix: 2d_quad_q1_p0_gmg_vanka_add 1127*c4762a1bSJed Brown requires: double !complex long_runtime 1128*c4762a1bSJed Brown args: -run_type full -bc_type dirichlet -simplex 0 -dm_refine_hierarchy 3 -interpolate 1 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \ 1129*c4762a1bSJed Brown -snes_rtol 1.0e-4 -snes_error_if_not_converged -snes_view -snes_monitor -snes_converged_reason \ 1130*c4762a1bSJed Brown -ksp_type gmres -ksp_rtol 1.0e-5 -ksp_error_if_not_converged -ksp_monitor_true_residual \ 1131*c4762a1bSJed Brown -pc_type mg -pc_mg_levels 3 \ 1132*c4762a1bSJed Brown -mg_levels_ksp_type gmres -mg_levels_ksp_max_it 30 -mg_levels_ksp_monitor_true_residual_no \ 1133*c4762a1bSJed Brown -mg_levels_pc_type patch -mg_levels_pc_patch_partition_of_unity 0 -mg_levels_pc_patch_construct_codim 0 -mg_levels_pc_patch_construct_type vanka \ 1134*c4762a1bSJed Brown -mg_levels_sub_ksp_type preonly -mg_levels_sub_pc_type lu \ 1135*c4762a1bSJed Brown -mg_coarse_pc_type svd 1136*c4762a1bSJed Brown 1137*c4762a1bSJed Brown test: 1138*c4762a1bSJed Brown requires: !single 1139*c4762a1bSJed Brown suffix: bddc_quad 1140*c4762a1bSJed Brown nsize: 2 1141*c4762a1bSJed Brown args: -run_type full -dm_refine 1 -bc_type dirichlet -interpolate 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_view -snes_error_if_not_converged -dm_mat_type is -ksp_type gmres -ksp_rtol 1.e-8 -pc_type bddc -pc_bddc_corner_selection -pc_bddc_dirichlet_pc_type svd -pc_bddc_neumann_pc_type svd -pc_bddc_coarse_redundant_pc_type svd -simplex 0 -petscpartitioner_type simple -ksp_monitor_short -pc_bddc_symmetric 0 1142*c4762a1bSJed Brown 1143*c4762a1bSJed Brown TEST*/ 1144