1c4762a1bSJed Brown static char help[] = "Poisson Problem in 2d and 3d with simplicial finite elements.\n\ 2c4762a1bSJed Brown We solve the Poisson problem in a rectangular\n\ 3c4762a1bSJed Brown domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\ 4c4762a1bSJed Brown This example supports discretized auxiliary fields (conductivity) as well as\n\ 5c4762a1bSJed Brown multilevel nonlinear solvers.\n\n\n"; 6c4762a1bSJed Brown 7c4762a1bSJed Brown /* 8c4762a1bSJed Brown A visualization of the adaptation can be accomplished using: 9c4762a1bSJed Brown 10c4762a1bSJed Brown -dm_adapt_view hdf5:$PWD/adapt.h5 -sol_adapt_view hdf5:$PWD/adapt.h5::append -dm_adapt_pre_view hdf5:$PWD/orig.h5 -sol_adapt_pre_view hdf5:$PWD/orig.h5::append 11c4762a1bSJed Brown 12c4762a1bSJed Brown Information on refinement: 13c4762a1bSJed Brown 14c20d7725SJed Brown -info :~sys,vec,is,mat,ksp,snes,ts 15c4762a1bSJed Brown */ 16c4762a1bSJed Brown 17c4762a1bSJed Brown #include <petscdmplex.h> 18c4762a1bSJed Brown #include <petscdmadaptor.h> 19c4762a1bSJed Brown #include <petscsnes.h> 20c4762a1bSJed Brown #include <petscds.h> 21c4762a1bSJed Brown #include <petscviewerhdf5.h> 22c4762a1bSJed Brown 23c4762a1bSJed Brown typedef enum {NEUMANN, DIRICHLET, NONE} BCType; 24c4762a1bSJed Brown typedef enum {RUN_FULL, RUN_EXACT, RUN_TEST, RUN_PERF} RunType; 258d1b37daSJoe Wallwork typedef enum {COEFF_NONE, COEFF_ANALYTIC, COEFF_FIELD, COEFF_NONLINEAR, COEFF_BALL, COEFF_CROSS, COEFF_CHECKERBOARD_0, COEFF_CHECKERBOARD_1} CoeffType; 26c4762a1bSJed Brown 27c4762a1bSJed Brown typedef struct { 28c4762a1bSJed Brown RunType runType; /* Whether to run tests, or solve the full problem */ 29c4762a1bSJed Brown PetscBool jacobianMF; /* Whether to calculate the Jacobian action on the fly */ 30c4762a1bSJed Brown PetscBool showInitial, showSolution, restart, quiet, nonzInit; 31c4762a1bSJed Brown /* Problem definition */ 32c4762a1bSJed Brown BCType bcType; 33c4762a1bSJed Brown CoeffType variableCoefficient; 34c4762a1bSJed Brown PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx); 35c4762a1bSJed Brown PetscBool fieldBC; 36c4762a1bSJed Brown void (**exactFields)(PetscInt, PetscInt, PetscInt, 37c4762a1bSJed Brown const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 38c4762a1bSJed Brown const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 39c4762a1bSJed Brown PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]); 40c4762a1bSJed Brown PetscBool bdIntegral; /* Compute the integral of the solution on the boundary */ 41d6837840SMatthew G. Knepley /* Reproducing tests from SISC 40(3), pp. A1473-A1493, 2018 */ 42d6837840SMatthew G. Knepley PetscInt div; /* Number of divisions */ 43d6837840SMatthew G. Knepley PetscInt k; /* Parameter for checkerboard coefficient */ 44d6837840SMatthew G. Knepley PetscInt *kgrid; /* Random parameter grid */ 4530602db0SMatthew G. Knepley PetscBool rand; /* Make random assignments */ 46c4762a1bSJed Brown /* Solver */ 47c4762a1bSJed Brown PC pcmg; /* This is needed for error monitoring */ 48c4762a1bSJed Brown PetscBool checkksp; /* Whether to check the KSPSolve for runType == RUN_TEST */ 49c4762a1bSJed Brown } AppCtx; 50c4762a1bSJed Brown 51c4762a1bSJed Brown static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 52c4762a1bSJed Brown { 53c4762a1bSJed Brown u[0] = 0.0; 54c4762a1bSJed Brown return 0; 55c4762a1bSJed Brown } 56c4762a1bSJed Brown 57c4762a1bSJed Brown static PetscErrorCode ecks(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 58c4762a1bSJed Brown { 59c4762a1bSJed Brown u[0] = x[0]; 60c4762a1bSJed Brown return 0; 61c4762a1bSJed Brown } 62c4762a1bSJed Brown 63c4762a1bSJed Brown /* 64c4762a1bSJed Brown In 2D for Dirichlet conditions, we use exact solution: 65c4762a1bSJed Brown 66c4762a1bSJed Brown u = x^2 + y^2 67c4762a1bSJed Brown f = 4 68c4762a1bSJed Brown 69c4762a1bSJed Brown so that 70c4762a1bSJed Brown 71c4762a1bSJed Brown -\Delta u + f = -4 + 4 = 0 72c4762a1bSJed Brown 73c4762a1bSJed Brown For Neumann conditions, we have 74c4762a1bSJed Brown 75c4762a1bSJed Brown -\nabla u \cdot -\hat y |_{y=0} = (2y)|_{y=0} = 0 (bottom) 76c4762a1bSJed Brown -\nabla u \cdot \hat y |_{y=1} = -(2y)|_{y=1} = -2 (top) 77c4762a1bSJed Brown -\nabla u \cdot -\hat x |_{x=0} = (2x)|_{x=0} = 0 (left) 78c4762a1bSJed Brown -\nabla u \cdot \hat x |_{x=1} = -(2x)|_{x=1} = -2 (right) 79c4762a1bSJed Brown 80c4762a1bSJed Brown Which we can express as 81c4762a1bSJed Brown 82c4762a1bSJed Brown \nabla u \cdot \hat n|_\Gamma = {2 x, 2 y} \cdot \hat n = 2 (x + y) 83c4762a1bSJed Brown 84c4762a1bSJed Brown The boundary integral of this solution is (assuming we are not orienting the edges) 85c4762a1bSJed Brown 86c4762a1bSJed Brown \int^1_0 x^2 dx + \int^1_0 (1 + y^2) dy + \int^1_0 (x^2 + 1) dx + \int^1_0 y^2 dy = 1/3 + 4/3 + 4/3 + 1/3 = 3 1/3 87c4762a1bSJed Brown */ 88c4762a1bSJed Brown static PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 89c4762a1bSJed Brown { 90c4762a1bSJed Brown *u = x[0]*x[0] + x[1]*x[1]; 91c4762a1bSJed Brown return 0; 92c4762a1bSJed Brown } 93c4762a1bSJed Brown 94c4762a1bSJed Brown static void quadratic_u_field_2d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 95c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 96c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 97c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar uexact[]) 98c4762a1bSJed Brown { 99c4762a1bSJed Brown uexact[0] = a[0]; 100c4762a1bSJed Brown } 101c4762a1bSJed Brown 1028d1b37daSJoe Wallwork static PetscErrorCode ball_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 103c4762a1bSJed Brown { 104c4762a1bSJed Brown const PetscReal alpha = 500.; 105c4762a1bSJed Brown const PetscReal radius2 = PetscSqr(0.15); 106c4762a1bSJed Brown const PetscReal r2 = PetscSqr(x[0] - 0.5) + PetscSqr(x[1] - 0.5); 107c4762a1bSJed Brown const PetscReal xi = alpha*(radius2 - r2); 108c4762a1bSJed Brown 109c4762a1bSJed Brown *u = PetscTanhScalar(xi) + 1.0; 110c4762a1bSJed Brown return 0; 111c4762a1bSJed Brown } 112c4762a1bSJed Brown 113c4762a1bSJed Brown static PetscErrorCode cross_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 114c4762a1bSJed Brown { 115c4762a1bSJed Brown const PetscReal alpha = 50*4; 116c4762a1bSJed Brown const PetscReal xy = (x[0]-0.5)*(x[1]-0.5); 117c4762a1bSJed Brown 118c4762a1bSJed Brown *u = PetscSinReal(alpha*xy) * (alpha*PetscAbsReal(xy) < 2*PETSC_PI ? 1 : 0.01); 119c4762a1bSJed Brown return 0; 120c4762a1bSJed Brown } 121c4762a1bSJed Brown 122c4762a1bSJed Brown static void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 123c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 124c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 125c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 126c4762a1bSJed Brown { 127c4762a1bSJed Brown f0[0] = 4.0; 128c4762a1bSJed Brown } 129c4762a1bSJed Brown 1308d1b37daSJoe Wallwork static void f0_ball_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 131c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 132c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 133c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 134c4762a1bSJed Brown { 1358d1b37daSJoe Wallwork PetscInt d; 1368d1b37daSJoe Wallwork const PetscReal alpha = 500., radius2 = PetscSqr(0.15); 1378d1b37daSJoe Wallwork PetscReal r2, xi; 138c4762a1bSJed Brown 1398d1b37daSJoe Wallwork for (d = 0, r2 = 0.0; d < dim; ++d) r2 += PetscSqr(x[d] - 0.5); 1408d1b37daSJoe Wallwork xi = alpha*(radius2 - r2); 1418d1b37daSJoe Wallwork f0[0] = (-2.0*dim*alpha - 8.0*PetscSqr(alpha)*r2*PetscTanhReal(xi)) * PetscSqr(1.0/PetscCoshReal(xi)); 142c4762a1bSJed Brown } 143c4762a1bSJed Brown 1448d1b37daSJoe Wallwork static void f0_cross_u_2d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 145c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 146c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 147c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 148c4762a1bSJed Brown { 149c4762a1bSJed Brown const PetscReal alpha = 50*4; 150c4762a1bSJed Brown const PetscReal xy = (x[0]-0.5)*(x[1]-0.5); 151c4762a1bSJed Brown 152c4762a1bSJed Brown f0[0] = PetscSinReal(alpha*xy) * (alpha*PetscAbsReal(xy) < 2*PETSC_PI ? 1 : 0.01); 153c4762a1bSJed Brown } 154c4762a1bSJed Brown 155d6837840SMatthew G. Knepley static void f0_checkerboard_0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 156d6837840SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 157d6837840SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 158d6837840SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 159d6837840SMatthew G. Knepley { 160d6837840SMatthew G. Knepley f0[0] = -20.0*PetscExpReal(-(PetscSqr(x[0] - 0.5) + PetscSqr(x[1] - 0.5))); 161d6837840SMatthew G. Knepley } 162d6837840SMatthew G. Knepley 163c4762a1bSJed Brown static void f0_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 164c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 165c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 166c4762a1bSJed Brown PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 167c4762a1bSJed Brown { 168c4762a1bSJed Brown PetscInt d; 169c4762a1bSJed Brown for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += -n[d]*2.0*x[d]; 170c4762a1bSJed Brown } 171c4762a1bSJed Brown 172c4762a1bSJed Brown /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */ 173c4762a1bSJed Brown static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 174c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 175c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 176c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 177c4762a1bSJed Brown { 178c4762a1bSJed Brown PetscInt d; 179c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 180c4762a1bSJed Brown } 181c4762a1bSJed Brown 182c4762a1bSJed Brown /* < \nabla v, \nabla u + {\nabla u}^T > 183c4762a1bSJed Brown This just gives \nabla u, give the perdiagonal for the transpose */ 184c4762a1bSJed Brown static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 185c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 186c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 187c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 188c4762a1bSJed Brown { 189c4762a1bSJed Brown PetscInt d; 190c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 191c4762a1bSJed Brown } 192c4762a1bSJed Brown 193c4762a1bSJed Brown /* 194c4762a1bSJed Brown In 2D for x periodicity and y Dirichlet conditions, we use exact solution: 195c4762a1bSJed Brown 196c4762a1bSJed Brown u = sin(2 pi x) 197c4762a1bSJed Brown f = -4 pi^2 sin(2 pi x) 198c4762a1bSJed Brown 199c4762a1bSJed Brown so that 200c4762a1bSJed Brown 201c4762a1bSJed Brown -\Delta u + f = 4 pi^2 sin(2 pi x) - 4 pi^2 sin(2 pi x) = 0 202c4762a1bSJed Brown */ 203c4762a1bSJed Brown static PetscErrorCode xtrig_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 204c4762a1bSJed Brown { 205c4762a1bSJed Brown *u = PetscSinReal(2.0*PETSC_PI*x[0]); 206c4762a1bSJed Brown return 0; 207c4762a1bSJed Brown } 208c4762a1bSJed Brown 209c4762a1bSJed Brown static void f0_xtrig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 210c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 211c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 212c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 213c4762a1bSJed Brown { 214c4762a1bSJed Brown f0[0] = -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[0]); 215c4762a1bSJed Brown } 216c4762a1bSJed Brown 217c4762a1bSJed Brown /* 218c4762a1bSJed Brown In 2D for x-y periodicity, we use exact solution: 219c4762a1bSJed Brown 220c4762a1bSJed Brown u = sin(2 pi x) sin(2 pi y) 221c4762a1bSJed Brown f = -8 pi^2 sin(2 pi x) 222c4762a1bSJed Brown 223c4762a1bSJed Brown so that 224c4762a1bSJed Brown 225c4762a1bSJed Brown -\Delta u + f = 4 pi^2 sin(2 pi x) sin(2 pi y) + 4 pi^2 sin(2 pi x) sin(2 pi y) - 8 pi^2 sin(2 pi x) = 0 226c4762a1bSJed Brown */ 227c4762a1bSJed Brown static PetscErrorCode xytrig_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 228c4762a1bSJed Brown { 229c4762a1bSJed Brown *u = PetscSinReal(2.0*PETSC_PI*x[0])*PetscSinReal(2.0*PETSC_PI*x[1]); 230c4762a1bSJed Brown return 0; 231c4762a1bSJed Brown } 232c4762a1bSJed Brown 233c4762a1bSJed Brown static void f0_xytrig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 234c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 235c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 236c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 237c4762a1bSJed Brown { 238c4762a1bSJed Brown f0[0] = -8.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[0]); 239c4762a1bSJed Brown } 240c4762a1bSJed Brown 241c4762a1bSJed Brown /* 242c4762a1bSJed Brown In 2D for Dirichlet conditions with a variable coefficient, we use exact solution: 243c4762a1bSJed Brown 244c4762a1bSJed Brown u = x^2 + y^2 245c4762a1bSJed Brown f = 6 (x + y) 246c4762a1bSJed Brown nu = (x + y) 247c4762a1bSJed Brown 248c4762a1bSJed Brown so that 249c4762a1bSJed Brown 250c4762a1bSJed Brown -\div \nu \grad u + f = -6 (x + y) + 6 (x + y) = 0 251c4762a1bSJed Brown */ 252c4762a1bSJed Brown static PetscErrorCode nu_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 253c4762a1bSJed Brown { 254c4762a1bSJed Brown *u = x[0] + x[1]; 255c4762a1bSJed Brown return 0; 256c4762a1bSJed Brown } 257c4762a1bSJed Brown 258d6837840SMatthew G. Knepley static PetscErrorCode checkerboardCoeff(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 259d6837840SMatthew G. Knepley { 260d6837840SMatthew G. Knepley AppCtx *user = (AppCtx *) ctx; 261d6837840SMatthew G. Knepley PetscInt div = user->div; 262d6837840SMatthew G. Knepley PetscInt k = user->k; 263d6837840SMatthew G. Knepley PetscInt mask = 0, ind = 0, d; 264d6837840SMatthew G. Knepley 265d6837840SMatthew G. Knepley PetscFunctionBeginUser; 266d6837840SMatthew G. Knepley for (d = 0; d < dim; ++d) mask = (mask + (PetscInt) (x[d]*div)) % 2; 267d6837840SMatthew G. Knepley if (user->kgrid) { 268d6837840SMatthew G. Knepley for (d = 0; d < dim; ++d) { 269d6837840SMatthew G. Knepley if (d > 0) ind *= dim; 270d6837840SMatthew G. Knepley ind += (PetscInt) (x[d]*div); 271d6837840SMatthew G. Knepley } 272d6837840SMatthew G. Knepley k = user->kgrid[ind]; 273d6837840SMatthew G. Knepley } 274d6837840SMatthew G. Knepley u[0] = mask ? 1.0 : PetscPowRealInt(10.0, -k); 275d6837840SMatthew G. Knepley PetscFunctionReturn(0); 276d6837840SMatthew G. Knepley } 277d6837840SMatthew G. Knepley 278c4762a1bSJed Brown void f0_analytic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 279c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 280c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 281c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 282c4762a1bSJed Brown { 283c4762a1bSJed Brown f0[0] = 6.0*(x[0] + x[1]); 284c4762a1bSJed Brown } 285c4762a1bSJed Brown 286c4762a1bSJed Brown /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */ 287c4762a1bSJed Brown void f1_analytic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 288c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 289c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 290c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 291c4762a1bSJed Brown { 292c4762a1bSJed Brown PetscInt d; 293c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = (x[0] + x[1])*u_x[d]; 294c4762a1bSJed Brown } 295c4762a1bSJed Brown 296c4762a1bSJed Brown void f1_field_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 297c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 298c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 299c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 300c4762a1bSJed Brown { 301c4762a1bSJed Brown PetscInt d; 302c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = a[0]*u_x[d]; 303c4762a1bSJed Brown } 304c4762a1bSJed Brown 305c4762a1bSJed Brown /* < \nabla v, \nabla u + {\nabla u}^T > 306c4762a1bSJed Brown This just gives \nabla u, give the perdiagonal for the transpose */ 307c4762a1bSJed Brown void g3_analytic_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 308c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 309c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 310c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 311c4762a1bSJed Brown { 312c4762a1bSJed Brown PetscInt d; 313c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d*dim+d] = x[0] + x[1]; 314c4762a1bSJed Brown } 315c4762a1bSJed Brown 316c4762a1bSJed Brown void g3_field_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 317c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 318c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 319c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 320c4762a1bSJed Brown { 321c4762a1bSJed Brown PetscInt d; 322c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d*dim+d] = a[0]; 323c4762a1bSJed Brown } 324c4762a1bSJed Brown 325c4762a1bSJed Brown /* 326c4762a1bSJed Brown In 2D for Dirichlet conditions with a nonlinear coefficient (p-Laplacian with p = 4), we use exact solution: 327c4762a1bSJed Brown 328c4762a1bSJed Brown u = x^2 + y^2 329c4762a1bSJed Brown f = 16 (x^2 + y^2) 330c4762a1bSJed Brown nu = 1/2 |grad u|^2 331c4762a1bSJed Brown 332c4762a1bSJed Brown so that 333c4762a1bSJed Brown 334c4762a1bSJed Brown -\div \nu \grad u + f = -16 (x^2 + y^2) + 16 (x^2 + y^2) = 0 335c4762a1bSJed Brown */ 336c4762a1bSJed Brown void f0_analytic_nonlinear_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 337c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 338c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 339c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 340c4762a1bSJed Brown { 341c4762a1bSJed Brown f0[0] = 16.0*(x[0]*x[0] + x[1]*x[1]); 342c4762a1bSJed Brown } 343c4762a1bSJed Brown 344c4762a1bSJed Brown /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */ 345c4762a1bSJed Brown void f1_analytic_nonlinear_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 346c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 347c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 348c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 349c4762a1bSJed Brown { 350c4762a1bSJed Brown PetscScalar nu = 0.0; 351c4762a1bSJed Brown PetscInt d; 352c4762a1bSJed Brown for (d = 0; d < dim; ++d) nu += u_x[d]*u_x[d]; 353c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = 0.5*nu*u_x[d]; 354c4762a1bSJed Brown } 355c4762a1bSJed Brown 356c4762a1bSJed Brown /* 357c4762a1bSJed Brown grad (u + eps w) - grad u = eps grad w 358c4762a1bSJed Brown 359c4762a1bSJed Brown 1/2 |grad (u + eps w)|^2 grad (u + eps w) - 1/2 |grad u|^2 grad u 360c4762a1bSJed Brown = 1/2 (|grad u|^2 + 2 eps <grad u,grad w>) (grad u + eps grad w) - 1/2 |grad u|^2 grad u 361c4762a1bSJed Brown = 1/2 (eps |grad u|^2 grad w + 2 eps <grad u,grad w> grad u) 362c4762a1bSJed Brown = eps (1/2 |grad u|^2 grad w + grad u <grad u,grad w>) 363c4762a1bSJed Brown */ 364c4762a1bSJed Brown void g3_analytic_nonlinear_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 365c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 366c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 367c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 368c4762a1bSJed Brown { 369c4762a1bSJed Brown PetscScalar nu = 0.0; 370c4762a1bSJed Brown PetscInt d, e; 371c4762a1bSJed Brown for (d = 0; d < dim; ++d) nu += u_x[d]*u_x[d]; 372c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 373c4762a1bSJed Brown g3[d*dim+d] = 0.5*nu; 374c4762a1bSJed Brown for (e = 0; e < dim; ++e) { 375c4762a1bSJed Brown g3[d*dim+e] += u_x[d]*u_x[e]; 376c4762a1bSJed Brown } 377c4762a1bSJed Brown } 378c4762a1bSJed Brown } 379c4762a1bSJed Brown 380c4762a1bSJed Brown /* 381c4762a1bSJed Brown In 3D for Dirichlet conditions we use exact solution: 382c4762a1bSJed Brown 383c4762a1bSJed Brown u = 2/3 (x^2 + y^2 + z^2) 384c4762a1bSJed Brown f = 4 385c4762a1bSJed Brown 386c4762a1bSJed Brown so that 387c4762a1bSJed Brown 388c4762a1bSJed Brown -\Delta u + f = -2/3 * 6 + 4 = 0 389c4762a1bSJed Brown 390c4762a1bSJed Brown For Neumann conditions, we have 391c4762a1bSJed Brown 392c4762a1bSJed Brown -\nabla u \cdot -\hat z |_{z=0} = (2z)|_{z=0} = 0 (bottom) 393c4762a1bSJed Brown -\nabla u \cdot \hat z |_{z=1} = -(2z)|_{z=1} = -2 (top) 394c4762a1bSJed Brown -\nabla u \cdot -\hat y |_{y=0} = (2y)|_{y=0} = 0 (front) 395c4762a1bSJed Brown -\nabla u \cdot \hat y |_{y=1} = -(2y)|_{y=1} = -2 (back) 396c4762a1bSJed Brown -\nabla u \cdot -\hat x |_{x=0} = (2x)|_{x=0} = 0 (left) 397c4762a1bSJed Brown -\nabla u \cdot \hat x |_{x=1} = -(2x)|_{x=1} = -2 (right) 398c4762a1bSJed Brown 399c4762a1bSJed Brown Which we can express as 400c4762a1bSJed Brown 401c4762a1bSJed Brown \nabla u \cdot \hat n|_\Gamma = {2 x, 2 y, 2z} \cdot \hat n = 2 (x + y + z) 402c4762a1bSJed Brown */ 403c4762a1bSJed Brown static PetscErrorCode quadratic_u_3d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 404c4762a1bSJed Brown { 405c4762a1bSJed Brown *u = 2.0*(x[0]*x[0] + x[1]*x[1] + x[2]*x[2])/3.0; 406c4762a1bSJed Brown return 0; 407c4762a1bSJed Brown } 408c4762a1bSJed Brown 4098d1b37daSJoe Wallwork static PetscErrorCode ball_u_3d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 4108d1b37daSJoe Wallwork { 4118d1b37daSJoe Wallwork const PetscReal alpha = 500.; 4128d1b37daSJoe Wallwork const PetscReal radius2 = PetscSqr(0.15); 4138d1b37daSJoe Wallwork const PetscReal r2 = PetscSqr(x[0] - 0.5) + PetscSqr(x[1] - 0.5) + PetscSqr(x[2] - 0.5); 4148d1b37daSJoe Wallwork const PetscReal xi = alpha*(radius2 - r2); 4158d1b37daSJoe Wallwork 4168d1b37daSJoe Wallwork *u = PetscTanhScalar(xi) + 1.0; 4178d1b37daSJoe Wallwork return 0; 4188d1b37daSJoe Wallwork } 4198d1b37daSJoe Wallwork 420c4762a1bSJed Brown static void quadratic_u_field_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 421c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 422c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 423c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar uexact[]) 424c4762a1bSJed Brown { 425c4762a1bSJed Brown uexact[0] = a[0]; 426c4762a1bSJed Brown } 427c4762a1bSJed Brown 4288d1b37daSJoe Wallwork static PetscErrorCode cross_u_3d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 4298d1b37daSJoe Wallwork { 4308d1b37daSJoe Wallwork const PetscReal alpha = 50*4; 4318d1b37daSJoe Wallwork const PetscReal xyz = (x[0]-0.5)*(x[1]-0.5)*(x[2]-0.5); 4328d1b37daSJoe Wallwork 4338d1b37daSJoe Wallwork *u = PetscSinReal(alpha*xyz) * (alpha*PetscAbsReal(xyz) < 2*PETSC_PI ? (alpha*PetscAbsReal(xyz) > -2*PETSC_PI ? 1.0 : 0.01) : 0.01); 4348d1b37daSJoe Wallwork return 0; 4358d1b37daSJoe Wallwork } 4368d1b37daSJoe Wallwork 4378d1b37daSJoe Wallwork static void f0_cross_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 4388d1b37daSJoe Wallwork const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 4398d1b37daSJoe Wallwork const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 4408d1b37daSJoe Wallwork PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 4418d1b37daSJoe Wallwork { 4428d1b37daSJoe Wallwork const PetscReal alpha = 50*4; 4438d1b37daSJoe Wallwork const PetscReal xyz = (x[0]-0.5)*(x[1]-0.5)*(x[2]-0.5); 4448d1b37daSJoe Wallwork 4458d1b37daSJoe Wallwork f0[0] = PetscSinReal(alpha*xyz) * (alpha*PetscAbsReal(xyz) < 2*PETSC_PI ? (alpha*PetscAbsReal(xyz) > -2*PETSC_PI ? 1.0 : 0.01) : 0.01); 4468d1b37daSJoe Wallwork } 4478d1b37daSJoe Wallwork 448c4762a1bSJed Brown static void bd_integral_2d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 449c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 450c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 451c4762a1bSJed Brown PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *uint) 452c4762a1bSJed Brown { 453c4762a1bSJed Brown uint[0] = u[0]; 454c4762a1bSJed Brown } 455c4762a1bSJed Brown 456c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 457c4762a1bSJed Brown { 458c4762a1bSJed Brown const char *bcTypes[3] = {"neumann", "dirichlet", "none"}; 459c4762a1bSJed Brown const char *runTypes[4] = {"full", "exact", "test", "perf"}; 4608d1b37daSJoe Wallwork const char *coeffTypes[8] = {"none", "analytic", "field", "nonlinear", "ball", "cross", "checkerboard_0", "checkerboard_1"}; 46130602db0SMatthew G. Knepley PetscInt bc, run, coeff; 462c4762a1bSJed Brown PetscErrorCode ierr; 463c4762a1bSJed Brown 464c4762a1bSJed Brown PetscFunctionBeginUser; 465c4762a1bSJed Brown options->runType = RUN_FULL; 466c4762a1bSJed Brown options->bcType = DIRICHLET; 467c4762a1bSJed Brown options->variableCoefficient = COEFF_NONE; 468c4762a1bSJed Brown options->fieldBC = PETSC_FALSE; 469c4762a1bSJed Brown options->jacobianMF = PETSC_FALSE; 470c4762a1bSJed Brown options->showInitial = PETSC_FALSE; 471c4762a1bSJed Brown options->showSolution = PETSC_FALSE; 472c4762a1bSJed Brown options->restart = PETSC_FALSE; 473c4762a1bSJed Brown options->quiet = PETSC_FALSE; 474c4762a1bSJed Brown options->nonzInit = PETSC_FALSE; 475c4762a1bSJed Brown options->bdIntegral = PETSC_FALSE; 476c4762a1bSJed Brown options->checkksp = PETSC_FALSE; 477d6837840SMatthew G. Knepley options->div = 4; 478d6837840SMatthew G. Knepley options->k = 1; 479d6837840SMatthew G. Knepley options->kgrid = NULL; 48030602db0SMatthew G. Knepley options->rand = PETSC_FALSE; 481c4762a1bSJed Brown 482c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");CHKERRQ(ierr); 483c4762a1bSJed Brown run = options->runType; 4845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEList("-run_type", "The run type", "ex12.c", runTypes, 4, runTypes[options->runType], &run, NULL)); 485c4762a1bSJed Brown options->runType = (RunType) run; 486c4762a1bSJed Brown bc = options->bcType; 4875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEList("-bc_type","Type of boundary condition","ex12.c",bcTypes,3,bcTypes[options->bcType],&bc,NULL)); 488c4762a1bSJed Brown options->bcType = (BCType) bc; 489c4762a1bSJed Brown coeff = options->variableCoefficient; 4905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEList("-variable_coefficient","Type of variable coefficent","ex12.c",coeffTypes,8,coeffTypes[options->variableCoefficient],&coeff,NULL)); 491c4762a1bSJed Brown options->variableCoefficient = (CoeffType) coeff; 492c4762a1bSJed Brown 4935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-field_bc", "Use a field representation for the BC", "ex12.c", options->fieldBC, &options->fieldBC, NULL)); 4945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-jacobian_mf", "Calculate the action of the Jacobian on the fly", "ex12.c", options->jacobianMF, &options->jacobianMF, NULL)); 4955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-show_initial", "Output the initial guess for verification", "ex12.c", options->showInitial, &options->showInitial, NULL)); 4965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-show_solution", "Output the solution for verification", "ex12.c", options->showSolution, &options->showSolution, NULL)); 4975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-restart", "Read in the mesh and solution from a file", "ex12.c", options->restart, &options->restart, NULL)); 4985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-quiet", "Don't print any vecs", "ex12.c", options->quiet, &options->quiet, NULL)); 4995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-nonzero_initial_guess", "nonzero initial guess", "ex12.c", options->nonzInit, &options->nonzInit, NULL)); 5005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-bd_integral", "Compute the integral of the solution on the boundary", "ex12.c", options->bdIntegral, &options->bdIntegral, NULL)); 501c4762a1bSJed Brown if (options->runType == RUN_TEST) { 5025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-run_test_check_ksp", "Check solution of KSP", "ex12.c", options->checkksp, &options->checkksp, NULL)); 503c4762a1bSJed Brown } 5045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-div", "The number of division for the checkerboard coefficient", "ex12.c", options->div, &options->div, NULL)); 5055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-k", "The exponent for the checkerboard coefficient", "ex12.c", options->k, &options->k, NULL)); 5065f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-k_random", "Assign random k values to checkerboard", "ex12.c", options->rand, &options->rand, NULL)); 5071e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 508c4762a1bSJed Brown PetscFunctionReturn(0); 509c4762a1bSJed Brown } 510c4762a1bSJed Brown 511c4762a1bSJed Brown static PetscErrorCode CreateBCLabel(DM dm, const char name[]) 512c4762a1bSJed Brown { 513408cafa0SMatthew G. Knepley DM plex; 514c4762a1bSJed Brown DMLabel label; 515c4762a1bSJed Brown 516c4762a1bSJed Brown PetscFunctionBeginUser; 5175f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLabel(dm, name)); 5185f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, name, &label)); 5195f80ce2aSJacob Faibussowitsch CHKERRQ(DMConvert(dm, DMPLEX, &plex)); 5205f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexMarkBoundaryFaces(plex, 1, label)); 5215f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&plex)); 522c4762a1bSJed Brown PetscFunctionReturn(0); 523c4762a1bSJed Brown } 524c4762a1bSJed Brown 525c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 526c4762a1bSJed Brown { 527c4762a1bSJed Brown PetscErrorCode ierr; 528c4762a1bSJed Brown 529c4762a1bSJed Brown PetscFunctionBeginUser; 5305f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 5315f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 5325f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 533c4762a1bSJed Brown { 534c4762a1bSJed Brown char convType[256]; 535c4762a1bSJed Brown PetscBool flg; 536c4762a1bSJed Brown 537c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");CHKERRQ(ierr); 5385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsFList("-dm_plex_convert_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg)); 5391e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 540c4762a1bSJed Brown if (flg) { 541c4762a1bSJed Brown DM dmConv; 542c4762a1bSJed Brown 5435f80ce2aSJacob Faibussowitsch CHKERRQ(DMConvert(*dm,convType,&dmConv)); 544c4762a1bSJed Brown if (dmConv) { 5455f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 546c4762a1bSJed Brown *dm = dmConv; 547c4762a1bSJed Brown } 5485f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 5495f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(*dm)); 55030602db0SMatthew G. Knepley } 55130602db0SMatthew G. Knepley } 5525f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 55330602db0SMatthew G. Knepley if (user->rand) { 55430602db0SMatthew G. Knepley PetscRandom r; 55530602db0SMatthew G. Knepley PetscReal val; 55630602db0SMatthew G. Knepley PetscInt dim, N, i; 557c4762a1bSJed Brown 5585f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(*dm, &dim)); 55930602db0SMatthew G. Knepley N = PetscPowInt(user->div, dim); 5605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(N, &user->kgrid)); 5615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF, &r)); 5625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(r)); 5635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetInterval(r, 0.0, user->k)); 5645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetSeed(r, 1973)); 5655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSeed(r)); 56630602db0SMatthew G. Knepley for (i = 0; i < N; ++i) { 5675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValueReal(r, &val)); 56830602db0SMatthew G. Knepley user->kgrid[i] = 1 + (PetscInt) val; 569c4762a1bSJed Brown } 5705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&r)); 571c4762a1bSJed Brown } 572c4762a1bSJed Brown PetscFunctionReturn(0); 573c4762a1bSJed Brown } 574c4762a1bSJed Brown 575c4762a1bSJed Brown static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 576c4762a1bSJed Brown { 57745480ffeSMatthew G. Knepley PetscDS ds; 57845480ffeSMatthew G. Knepley DMLabel label; 57945480ffeSMatthew G. Knepley PetscWeakForm wf; 58030602db0SMatthew G. Knepley const DMBoundaryType *periodicity; 581c4762a1bSJed Brown const PetscInt id = 1; 58230602db0SMatthew G. Knepley PetscInt bd, dim; 583c4762a1bSJed Brown 584c4762a1bSJed Brown PetscFunctionBeginUser; 5855f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm, &ds)); 5865f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 5875f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetPeriodicity(dm, NULL, NULL, NULL, &periodicity)); 588c4762a1bSJed Brown switch (user->variableCoefficient) { 589c4762a1bSJed Brown case COEFF_NONE: 59030602db0SMatthew G. Knepley if (periodicity && periodicity[0]) { 59130602db0SMatthew G. Knepley if (periodicity && periodicity[1]) { 5925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, f0_xytrig_u, f1_u)); 5935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 594c4762a1bSJed Brown } else { 5955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, f0_xtrig_u, f1_u)); 5965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 597c4762a1bSJed Brown } 598c4762a1bSJed Brown } else { 5995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, f0_u, f1_u)); 6005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 601c4762a1bSJed Brown } 602c4762a1bSJed Brown break; 603c4762a1bSJed Brown case COEFF_ANALYTIC: 6045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, f0_analytic_u, f1_analytic_u)); 6055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_analytic_uu)); 606c4762a1bSJed Brown break; 607c4762a1bSJed Brown case COEFF_FIELD: 6085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, f0_analytic_u, f1_field_u)); 6095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_field_uu)); 610c4762a1bSJed Brown break; 611c4762a1bSJed Brown case COEFF_NONLINEAR: 6125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, f0_analytic_nonlinear_u, f1_analytic_nonlinear_u)); 6135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_analytic_nonlinear_uu)); 614c4762a1bSJed Brown break; 6158d1b37daSJoe Wallwork case COEFF_BALL: 6165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, f0_ball_u, f1_u)); 6175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 618c4762a1bSJed Brown break; 619c4762a1bSJed Brown case COEFF_CROSS: 6208d1b37daSJoe Wallwork switch (dim) { 6218d1b37daSJoe Wallwork case 2: 6225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, f0_cross_u_2d, f1_u)); 6238d1b37daSJoe Wallwork break; 6248d1b37daSJoe Wallwork case 3: 6255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, f0_cross_u_3d, f1_u)); 6268d1b37daSJoe Wallwork break; 6278d1b37daSJoe Wallwork default: 62898921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %d", dim); 6298d1b37daSJoe Wallwork } 6305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 631c4762a1bSJed Brown break; 632d6837840SMatthew G. Knepley case COEFF_CHECKERBOARD_0: 6335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, f0_checkerboard_0_u, f1_field_u)); 6345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_field_uu)); 635d6837840SMatthew G. Knepley break; 63698921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid variable coefficient type %d", user->variableCoefficient); 637c4762a1bSJed Brown } 63830602db0SMatthew G. Knepley switch (dim) { 639c4762a1bSJed Brown case 2: 640c4762a1bSJed Brown switch (user->variableCoefficient) { 6418d1b37daSJoe Wallwork case COEFF_BALL: 6428d1b37daSJoe Wallwork user->exactFuncs[0] = ball_u_2d;break; 643c4762a1bSJed Brown case COEFF_CROSS: 644c4762a1bSJed Brown user->exactFuncs[0] = cross_u_2d;break; 645d6837840SMatthew G. Knepley case COEFF_CHECKERBOARD_0: 646d6837840SMatthew G. Knepley user->exactFuncs[0] = zero;break; 647c4762a1bSJed Brown default: 64830602db0SMatthew G. Knepley if (periodicity && periodicity[0]) { 64930602db0SMatthew G. Knepley if (periodicity && periodicity[1]) { 650c4762a1bSJed Brown user->exactFuncs[0] = xytrig_u_2d; 651c4762a1bSJed Brown } else { 652c4762a1bSJed Brown user->exactFuncs[0] = xtrig_u_2d; 653c4762a1bSJed Brown } 654c4762a1bSJed Brown } else { 655c4762a1bSJed Brown user->exactFuncs[0] = quadratic_u_2d; 656c4762a1bSJed Brown user->exactFields[0] = quadratic_u_field_2d; 657c4762a1bSJed Brown } 658c4762a1bSJed Brown } 65945480ffeSMatthew G. Knepley if (user->bcType == NEUMANN) { 6605f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "boundary", &label)); 6615f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddBoundary(dm, DM_BC_NATURAL, "wall", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 6625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 6635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_bd_u, 0, NULL)); 66445480ffeSMatthew G. Knepley } 665c4762a1bSJed Brown break; 666c4762a1bSJed Brown case 3: 6678d1b37daSJoe Wallwork switch (user->variableCoefficient) { 6688d1b37daSJoe Wallwork case COEFF_BALL: 6698d1b37daSJoe Wallwork user->exactFuncs[0] = ball_u_3d;break; 6708d1b37daSJoe Wallwork case COEFF_CROSS: 6718d1b37daSJoe Wallwork user->exactFuncs[0] = cross_u_3d;break; 6728d1b37daSJoe Wallwork default: 673c4762a1bSJed Brown user->exactFuncs[0] = quadratic_u_3d; 674c4762a1bSJed Brown user->exactFields[0] = quadratic_u_field_3d; 6758d1b37daSJoe Wallwork } 67645480ffeSMatthew G. Knepley if (user->bcType == NEUMANN) { 6775f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "boundary", &label)); 6785f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddBoundary(dm, DM_BC_NATURAL, "wall", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 6795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 6805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_bd_u, 0, NULL)); 68145480ffeSMatthew G. Knepley } 682c4762a1bSJed Brown break; 683c4762a1bSJed Brown default: 68498921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %d", dim); 685c4762a1bSJed Brown } 686d6837840SMatthew G. Knepley /* Setup constants */ 687d6837840SMatthew G. Knepley switch (user->variableCoefficient) { 688d6837840SMatthew G. Knepley case COEFF_CHECKERBOARD_0: 689d6837840SMatthew G. Knepley { 690d6837840SMatthew G. Knepley PetscScalar constants[2]; 691d6837840SMatthew G. Knepley 692d6837840SMatthew G. Knepley constants[0] = user->div; 693d6837840SMatthew G. Knepley constants[1] = user->k; 6945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetConstants(ds, 2, constants)); 695d6837840SMatthew G. Knepley } 696d6837840SMatthew G. Knepley break; 697d6837840SMatthew G. Knepley default: break; 698d6837840SMatthew G. Knepley } 6995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolution(ds, 0, user->exactFuncs[0], user)); 700d6837840SMatthew G. Knepley /* Setup Boundary Conditions */ 70145480ffeSMatthew G. Knepley if (user->bcType == DIRICHLET) { 7025f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "marker", &label)); 70345480ffeSMatthew G. Knepley if (!label) { 70445480ffeSMatthew G. Knepley /* Right now, p4est cannot create labels immediately */ 7055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSAddBoundaryByName(ds, user->fieldBC ? DM_BC_ESSENTIAL_FIELD : DM_BC_ESSENTIAL, "wall", "marker", 1, &id, 0, 0, NULL, user->fieldBC ? (void (*)(void)) user->exactFields[0] : (void (*)(void)) user->exactFuncs[0], NULL, user, NULL)); 70645480ffeSMatthew G. Knepley } else { 7075f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddBoundary(dm, user->fieldBC ? DM_BC_ESSENTIAL_FIELD : DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, user->fieldBC ? (void (*)(void)) user->exactFields[0] : (void (*)(void)) user->exactFuncs[0], NULL, user, NULL)); 70845480ffeSMatthew G. Knepley } 709c4762a1bSJed Brown } 710c4762a1bSJed Brown PetscFunctionReturn(0); 711c4762a1bSJed Brown } 712c4762a1bSJed Brown 713c4762a1bSJed Brown static PetscErrorCode SetupMaterial(DM dm, DM dmAux, AppCtx *user) 714c4762a1bSJed Brown { 715c4762a1bSJed Brown PetscErrorCode (*matFuncs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx) = {nu_2d}; 716d6837840SMatthew G. Knepley void *ctx[1]; 717c4762a1bSJed Brown Vec nu; 718c4762a1bSJed Brown 719c4762a1bSJed Brown PetscFunctionBegin; 720d6837840SMatthew G. Knepley ctx[0] = user; 721d6837840SMatthew G. Knepley if (user->variableCoefficient == COEFF_CHECKERBOARD_0) {matFuncs[0] = checkerboardCoeff;} 7225f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLocalVector(dmAux, &nu)); 7235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) nu, "Coefficient")); 7245f80ce2aSJacob Faibussowitsch CHKERRQ(DMProjectFunctionLocal(dmAux, 0.0, matFuncs, ctx, INSERT_ALL_VALUES, nu)); 7255f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetAuxiliaryVec(dm, NULL, 0, 0, nu)); 7265f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&nu)); 727c4762a1bSJed Brown PetscFunctionReturn(0); 728c4762a1bSJed Brown } 729c4762a1bSJed Brown 730c4762a1bSJed Brown static PetscErrorCode SetupBC(DM dm, DM dmAux, AppCtx *user) 731c4762a1bSJed Brown { 732c4762a1bSJed Brown PetscErrorCode (*bcFuncs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx); 733c4762a1bSJed Brown Vec uexact; 734c4762a1bSJed Brown PetscInt dim; 735c4762a1bSJed Brown 736c4762a1bSJed Brown PetscFunctionBegin; 7375f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 738c4762a1bSJed Brown if (dim == 2) bcFuncs[0] = quadratic_u_2d; 739c4762a1bSJed Brown else bcFuncs[0] = quadratic_u_3d; 7405f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLocalVector(dmAux, &uexact)); 7415f80ce2aSJacob Faibussowitsch CHKERRQ(DMProjectFunctionLocal(dmAux, 0.0, bcFuncs, NULL, INSERT_ALL_VALUES, uexact)); 7425f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetAuxiliaryVec(dm, NULL, 0, 0, uexact)); 7435f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&uexact)); 744c4762a1bSJed Brown PetscFunctionReturn(0); 745c4762a1bSJed Brown } 746c4762a1bSJed Brown 747c4762a1bSJed Brown static PetscErrorCode SetupAuxDM(DM dm, PetscFE feAux, AppCtx *user) 748c4762a1bSJed Brown { 749c4762a1bSJed Brown DM dmAux, coordDM; 750c4762a1bSJed Brown 751c4762a1bSJed Brown PetscFunctionBegin; 752c4762a1bSJed Brown /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */ 7535f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(dm, &coordDM)); 754c4762a1bSJed Brown if (!feAux) PetscFunctionReturn(0); 7555f80ce2aSJacob Faibussowitsch CHKERRQ(DMClone(dm, &dmAux)); 7565f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetCoordinateDM(dmAux, coordDM)); 7575f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dmAux, 0, NULL, (PetscObject) feAux)); 7585f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(dmAux)); 7595f80ce2aSJacob Faibussowitsch if (user->fieldBC) CHKERRQ(SetupBC(dm, dmAux, user)); 7605f80ce2aSJacob Faibussowitsch else CHKERRQ(SetupMaterial(dm, dmAux, user)); 7615f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dmAux)); 762c4762a1bSJed Brown PetscFunctionReturn(0); 763c4762a1bSJed Brown } 764c4762a1bSJed Brown 765c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 766c4762a1bSJed Brown { 76730602db0SMatthew G. Knepley DM plex, cdm = dm; 768c4762a1bSJed Brown PetscFE fe, feAux = NULL; 76930602db0SMatthew G. Knepley PetscBool simplex; 77030602db0SMatthew G. Knepley PetscInt dim; 771c4762a1bSJed Brown MPI_Comm comm; 772c4762a1bSJed Brown 773c4762a1bSJed Brown PetscFunctionBeginUser; 7745f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 7755f80ce2aSJacob Faibussowitsch CHKERRQ(DMConvert(dm, DMPLEX, &plex)); 7765f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexIsSimplex(plex, &simplex)); 7775f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&plex)); 7785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm)); 7795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, -1, &fe)); 7805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe, "potential")); 781d6837840SMatthew G. Knepley if (user->variableCoefficient == COEFF_FIELD || user->variableCoefficient == COEFF_CHECKERBOARD_0) { 7825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "mat_", -1, &feAux)); 7835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) feAux, "coefficient")); 7845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECopyQuadrature(fe, feAux)); 785c4762a1bSJed Brown } else if (user->fieldBC) { 7865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "bc_", -1, &feAux)); 7875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECopyQuadrature(fe, feAux)); 788c4762a1bSJed Brown } 789c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 7905f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe)); 7915f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(dm)); 7925f80ce2aSJacob Faibussowitsch CHKERRQ(SetupProblem(dm, user)); 793c4762a1bSJed Brown while (cdm) { 7945f80ce2aSJacob Faibussowitsch CHKERRQ(SetupAuxDM(cdm, feAux, user)); 79530602db0SMatthew G. Knepley if (user->bcType == DIRICHLET) { 796c4762a1bSJed Brown PetscBool hasLabel; 797c4762a1bSJed Brown 7985f80ce2aSJacob Faibussowitsch CHKERRQ(DMHasLabel(cdm, "marker", &hasLabel)); 7995f80ce2aSJacob Faibussowitsch if (!hasLabel) CHKERRQ(CreateBCLabel(cdm, "marker")); 800c4762a1bSJed Brown } 8015f80ce2aSJacob Faibussowitsch CHKERRQ(DMCopyDisc(dm, cdm)); 8025f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoarseDM(cdm, &cdm)); 803c4762a1bSJed Brown } 8045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe)); 8055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&feAux)); 806c4762a1bSJed Brown PetscFunctionReturn(0); 807c4762a1bSJed Brown } 808c4762a1bSJed Brown 809c4762a1bSJed Brown int main(int argc, char **argv) 810c4762a1bSJed Brown { 811c4762a1bSJed Brown DM dm; /* Problem specification */ 812c4762a1bSJed Brown SNES snes; /* nonlinear solver */ 813c4762a1bSJed Brown Vec u; /* solution vector */ 814c4762a1bSJed Brown Mat A,J; /* Jacobian matrix */ 815c4762a1bSJed Brown MatNullSpace nullSpace; /* May be necessary for Neumann conditions */ 816c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 817c4762a1bSJed Brown JacActionCtx userJ; /* context for Jacobian MF action */ 818c4762a1bSJed Brown PetscReal error = 0.0; /* L_2 error in the solution */ 819c4762a1bSJed Brown 820*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, NULL,help)); 8215f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user)); 8225f80ce2aSJacob Faibussowitsch CHKERRQ(SNESCreate(PETSC_COMM_WORLD, &snes)); 8235f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 8245f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetDM(snes, dm)); 8255f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetApplicationContext(dm, &user)); 826c4762a1bSJed Brown 8275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(1, &user.exactFuncs, 1, &user.exactFields)); 8285f80ce2aSJacob Faibussowitsch CHKERRQ(SetupDiscretization(dm, &user)); 829c4762a1bSJed Brown 8305f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dm, &u)); 8315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) u, "potential")); 832c4762a1bSJed Brown 8335f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(dm, &J)); 834c4762a1bSJed Brown if (user.jacobianMF) { 835c4762a1bSJed Brown PetscInt M, m, N, n; 836c4762a1bSJed Brown 8375f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(J, &M, &N)); 8385f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(J, &m, &n)); 8395f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD, &A)); 8405f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A, m, n, M, N)); 8415f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A, MATSHELL)); 8425f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 843c4762a1bSJed Brown #if 0 8445f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(A, MATOP_MULT, (void (*)(void))FormJacobianAction)); 845c4762a1bSJed Brown #endif 846c4762a1bSJed Brown 847c4762a1bSJed Brown userJ.dm = dm; 848c4762a1bSJed Brown userJ.J = J; 849c4762a1bSJed Brown userJ.user = &user; 850c4762a1bSJed Brown 8515f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLocalVector(dm, &userJ.u)); 8525f80ce2aSJacob Faibussowitsch if (user.fieldBC) CHKERRQ(DMProjectFieldLocal(dm, 0.0, userJ.u, user.exactFields, INSERT_BC_VALUES, userJ.u)); 8535f80ce2aSJacob Faibussowitsch else CHKERRQ(DMProjectFunctionLocal(dm, 0.0, user.exactFuncs, NULL, INSERT_BC_VALUES, userJ.u)); 8545f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetContext(A, &userJ)); 855c4762a1bSJed Brown } else { 856c4762a1bSJed Brown A = J; 857c4762a1bSJed Brown } 858c4762a1bSJed Brown 859c4762a1bSJed Brown nullSpace = NULL; 860c4762a1bSJed Brown if (user.bcType != DIRICHLET) { 8615f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_TRUE, 0, NULL, &nullSpace)); 8625f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetNullSpace(A, nullSpace)); 863c4762a1bSJed Brown } 864c4762a1bSJed Brown 8655f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetSNESLocalFEM(dm,&user,&user,&user)); 8665f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetJacobian(snes, A, J, NULL, NULL)); 867c4762a1bSJed Brown 8685f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFromOptions(snes)); 869c4762a1bSJed Brown 8705f80ce2aSJacob Faibussowitsch if (user.fieldBC) CHKERRQ(DMProjectField(dm, 0.0, u, user.exactFields, INSERT_ALL_VALUES, u)); 8715f80ce2aSJacob Faibussowitsch else CHKERRQ(DMProjectFunction(dm, 0.0, user.exactFuncs, NULL, INSERT_ALL_VALUES, u)); 872c4762a1bSJed Brown if (user.restart) { 873c4762a1bSJed Brown #if defined(PETSC_HAVE_HDF5) 874c4762a1bSJed Brown PetscViewer viewer; 87530602db0SMatthew G. Knepley char filename[PETSC_MAX_PATH_LEN]; 876c4762a1bSJed Brown 8775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL, NULL, "-dm_plex_filename", filename, sizeof(filename), NULL)); 8785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerCreate(PETSC_COMM_WORLD, &viewer)); 8795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerSetType(viewer, PETSCVIEWERHDF5)); 8805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFileSetMode(viewer, FILE_MODE_READ)); 8815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFileSetName(viewer, filename)); 8825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerHDF5PushGroup(viewer, "/fields")); 8835f80ce2aSJacob Faibussowitsch CHKERRQ(VecLoad(u, viewer)); 8845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerHDF5PopGroup(viewer)); 8855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 886c4762a1bSJed Brown #endif 887c4762a1bSJed Brown } 888c4762a1bSJed Brown if (user.showInitial) { 889c4762a1bSJed Brown Vec lv; 8905f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dm, &lv)); 8915f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(dm, u, INSERT_VALUES, lv)); 8925f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(dm, u, INSERT_VALUES, lv)); 8935f80ce2aSJacob Faibussowitsch CHKERRQ(DMPrintLocalVec(dm, "Local function", 1.0e-10, lv)); 8945f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dm, &lv)); 895c4762a1bSJed Brown } 896c4762a1bSJed Brown if (user.runType == RUN_FULL || user.runType == RUN_EXACT) { 897c4762a1bSJed Brown PetscErrorCode (*initialGuess[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx) = {zero}; 898c4762a1bSJed Brown 899c4762a1bSJed Brown if (user.nonzInit) initialGuess[0] = ecks; 900c4762a1bSJed Brown if (user.runType == RUN_FULL) { 9015f80ce2aSJacob Faibussowitsch CHKERRQ(DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u)); 902c4762a1bSJed Brown } 9035f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(u, NULL, "-guess_vec_view")); 9045f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSolve(snes, NULL, u)); 9055f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetSolution(snes, &u)); 9065f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetDM(snes, &dm)); 907c4762a1bSJed Brown 908c4762a1bSJed Brown if (user.showSolution) { 9095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Solution\n")); 9105f80ce2aSJacob Faibussowitsch CHKERRQ(VecChop(u, 3.0e-9)); 9115f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(u, PETSC_VIEWER_STDOUT_WORLD)); 912c4762a1bSJed Brown } 913c4762a1bSJed Brown } else if (user.runType == RUN_PERF) { 914c4762a1bSJed Brown Vec r; 915c4762a1bSJed Brown PetscReal res = 0.0; 916c4762a1bSJed Brown 9175f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetFunction(snes, &r, NULL, NULL)); 9185f80ce2aSJacob Faibussowitsch CHKERRQ(SNESComputeFunction(snes, u, r)); 9195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n")); 9205f80ce2aSJacob Faibussowitsch CHKERRQ(VecChop(r, 1.0e-10)); 9215f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(r, NORM_2, &res)); 9225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", (double)res)); 923c4762a1bSJed Brown } else { 924c4762a1bSJed Brown Vec r; 925c4762a1bSJed Brown PetscReal res = 0.0, tol = 1.0e-11; 926c4762a1bSJed Brown 927c4762a1bSJed Brown /* Check discretization error */ 9285f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetFunction(snes, &r, NULL, NULL)); 9295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n")); 9305f80ce2aSJacob Faibussowitsch if (!user.quiet) CHKERRQ(VecView(u, PETSC_VIEWER_STDOUT_WORLD)); 9315f80ce2aSJacob Faibussowitsch CHKERRQ(DMComputeL2Diff(dm, 0.0, user.exactFuncs, NULL, u, &error)); 9325f80ce2aSJacob Faibussowitsch if (error < tol) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: < %2.1e\n", (double)tol)); 9335f80ce2aSJacob Faibussowitsch else CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g\n", (double)error)); 934c4762a1bSJed Brown /* Check residual */ 9355f80ce2aSJacob Faibussowitsch CHKERRQ(SNESComputeFunction(snes, u, r)); 9365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n")); 9375f80ce2aSJacob Faibussowitsch CHKERRQ(VecChop(r, 1.0e-10)); 9385f80ce2aSJacob Faibussowitsch if (!user.quiet) CHKERRQ(VecView(r, PETSC_VIEWER_STDOUT_WORLD)); 9395f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(r, NORM_2, &res)); 9405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", (double)res)); 941c4762a1bSJed Brown /* Check Jacobian */ 942c4762a1bSJed Brown { 943c4762a1bSJed Brown Vec b; 944c4762a1bSJed Brown 9455f80ce2aSJacob Faibussowitsch CHKERRQ(SNESComputeJacobian(snes, u, A, A)); 9465f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u, &b)); 9475f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(r, 0.0)); 9485f80ce2aSJacob Faibussowitsch CHKERRQ(SNESComputeFunction(snes, r, b)); 9495f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A, u, r)); 9505f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(r, 1.0, b)); 9515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Au - b = Au + F(0)\n")); 9525f80ce2aSJacob Faibussowitsch CHKERRQ(VecChop(r, 1.0e-10)); 9535f80ce2aSJacob Faibussowitsch if (!user.quiet) CHKERRQ(VecView(r, PETSC_VIEWER_STDOUT_WORLD)); 9545f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(r, NORM_2, &res)); 9555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", (double)res)); 956c4762a1bSJed Brown /* check solver */ 957c4762a1bSJed Brown if (user.checkksp) { 958c4762a1bSJed Brown KSP ksp; 959c4762a1bSJed Brown 960c4762a1bSJed Brown if (nullSpace) { 9615f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceRemove(nullSpace, u)); 962c4762a1bSJed Brown } 9635f80ce2aSJacob Faibussowitsch CHKERRQ(SNESComputeJacobian(snes, u, A, J)); 9645f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A, u, b)); 9655f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetKSP(snes, &ksp)); 9665f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetOperators(ksp, A, J)); 9675f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSolve(ksp, b, r)); 9685f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(r, -1.0, u)); 9695f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(r, NORM_2, &res)); 9705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "KSP Error: %g\n", (double)res)); 971c4762a1bSJed Brown } 9725f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b)); 973c4762a1bSJed Brown } 974c4762a1bSJed Brown } 9755f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(u, NULL, "-vec_view")); 976d6837840SMatthew G. Knepley { 977d6837840SMatthew G. Knepley Vec nu; 978d6837840SMatthew G. Knepley 9795f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetAuxiliaryVec(dm, NULL, 0, 0, &nu)); 9805f80ce2aSJacob Faibussowitsch if (nu) CHKERRQ(VecViewFromOptions(nu, NULL, "-coeff_view")); 981d6837840SMatthew G. Knepley } 982c4762a1bSJed Brown 983c4762a1bSJed Brown if (user.bdIntegral) { 984c4762a1bSJed Brown DMLabel label; 985c4762a1bSJed Brown PetscInt id = 1; 986c4762a1bSJed Brown PetscScalar bdInt = 0.0; 987c4762a1bSJed Brown PetscReal exact = 3.3333333333; 988c4762a1bSJed Brown 9895f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "marker", &label)); 9905f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexComputeBdIntegral(dm, u, label, 1, &id, bd_integral_2d, &bdInt, NULL)); 9915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Solution boundary integral: %.4g\n", (double) PetscAbsScalar(bdInt))); 9922c71b3e2SJacob Faibussowitsch PetscCheckFalse(PetscAbsReal(PetscAbsScalar(bdInt) - exact) > PETSC_SQRT_MACHINE_EPSILON,PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Invalid boundary integral %g != %g", (double) PetscAbsScalar(bdInt), (double)exact); 993c4762a1bSJed Brown } 994c4762a1bSJed Brown 9955f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceDestroy(&nullSpace)); 9965f80ce2aSJacob Faibussowitsch if (user.jacobianMF) CHKERRQ(VecDestroy(&userJ.u)); 9975f80ce2aSJacob Faibussowitsch if (A != J) CHKERRQ(MatDestroy(&A)); 9985f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&J)); 9995f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 10005f80ce2aSJacob Faibussowitsch CHKERRQ(SNESDestroy(&snes)); 10015f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 10025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(user.exactFuncs, user.exactFields)); 10035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(user.kgrid)); 1004*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 1005*b122ec5aSJacob Faibussowitsch return 0; 1006c4762a1bSJed Brown } 1007c4762a1bSJed Brown 1008c4762a1bSJed Brown /*TEST 1009c4762a1bSJed Brown # 2D serial P1 test 0-4 1010c4762a1bSJed Brown test: 1011c4762a1bSJed Brown suffix: 2d_p1_0 1012c4762a1bSJed Brown requires: triangle 101330602db0SMatthew G. Knepley args: -run_type test -bc_type dirichlet -dm_plex_interpolate 0 -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1014c4762a1bSJed Brown 1015c4762a1bSJed Brown test: 1016c4762a1bSJed Brown suffix: 2d_p1_1 1017c4762a1bSJed Brown requires: triangle 101830602db0SMatthew G. Knepley args: -run_type test -bc_type dirichlet -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1019c4762a1bSJed Brown 1020c4762a1bSJed Brown test: 1021c4762a1bSJed Brown suffix: 2d_p1_2 1022c4762a1bSJed Brown requires: triangle 102330602db0SMatthew G. Knepley args: -run_type test -dm_refine_volume_limit_pre 0.0625 -bc_type dirichlet -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1024c4762a1bSJed Brown 1025c4762a1bSJed Brown test: 1026c4762a1bSJed Brown suffix: 2d_p1_neumann_0 1027c4762a1bSJed Brown requires: triangle 102830602db0SMatthew G. Knepley args: -dm_coord_space 0 -run_type test -bc_type neumann -dm_plex_boundary_label boundary -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 -dm_view ascii::ascii_info_detail 1029c4762a1bSJed Brown 1030c4762a1bSJed Brown test: 1031c4762a1bSJed Brown suffix: 2d_p1_neumann_1 1032c4762a1bSJed Brown requires: triangle 103330602db0SMatthew G. Knepley args: -run_type test -dm_refine_volume_limit_pre 0.0625 -bc_type neumann -dm_plex_boundary_label boundary -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1034c4762a1bSJed Brown 1035c4762a1bSJed Brown # 2D serial P2 test 5-8 1036c4762a1bSJed Brown test: 1037c4762a1bSJed Brown suffix: 2d_p2_0 1038c4762a1bSJed Brown requires: triangle 103930602db0SMatthew G. Knepley args: -run_type test -bc_type dirichlet -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1040c4762a1bSJed Brown 1041c4762a1bSJed Brown test: 1042c4762a1bSJed Brown suffix: 2d_p2_1 1043c4762a1bSJed Brown requires: triangle 104430602db0SMatthew G. Knepley args: -run_type test -dm_refine_volume_limit_pre 0.0625 -bc_type dirichlet -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1045c4762a1bSJed Brown 1046c4762a1bSJed Brown test: 1047c4762a1bSJed Brown suffix: 2d_p2_neumann_0 1048c4762a1bSJed Brown requires: triangle 104930602db0SMatthew G. Knepley args: -dm_coord_space 0 -run_type test -bc_type neumann -dm_plex_boundary_label boundary -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 -dm_view ascii::ascii_info_detail 1050c4762a1bSJed Brown 1051c4762a1bSJed Brown test: 1052c4762a1bSJed Brown suffix: 2d_p2_neumann_1 1053c4762a1bSJed Brown requires: triangle 105430602db0SMatthew G. Knepley args: -dm_coord_space 0 -run_type test -dm_refine_volume_limit_pre 0.0625 -bc_type neumann -dm_plex_boundary_label boundary -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 -dm_view ascii::ascii_info_detail 1055c4762a1bSJed Brown 1056c4762a1bSJed Brown test: 1057c4762a1bSJed Brown suffix: bd_int_0 1058c4762a1bSJed Brown requires: triangle 105930602db0SMatthew G. Knepley args: -run_type test -bc_type dirichlet -petscspace_degree 2 -bd_integral -dm_view -quiet 1060c4762a1bSJed Brown 1061c4762a1bSJed Brown test: 1062c4762a1bSJed Brown suffix: bd_int_1 1063c4762a1bSJed Brown requires: triangle 106430602db0SMatthew G. Knepley args: -run_type test -dm_refine 2 -bc_type dirichlet -petscspace_degree 2 -bd_integral -dm_view -quiet 1065c4762a1bSJed Brown 1066c4762a1bSJed Brown # 3D serial P1 test 9-12 1067c4762a1bSJed Brown test: 1068c4762a1bSJed Brown suffix: 3d_p1_0 1069c4762a1bSJed Brown requires: ctetgen 107030602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -bc_type dirichlet -dm_plex_interpolate 0 -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 -dm_view 1071c4762a1bSJed Brown 1072c4762a1bSJed Brown test: 1073c4762a1bSJed Brown suffix: 3d_p1_1 1074c4762a1bSJed Brown requires: ctetgen 107530602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -bc_type dirichlet -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 -dm_view 1076c4762a1bSJed Brown 1077c4762a1bSJed Brown test: 1078c4762a1bSJed Brown suffix: 3d_p1_2 1079c4762a1bSJed Brown requires: ctetgen 108030602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -dm_refine_volume_limit_pre 0.0125 -bc_type dirichlet -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 -dm_view 1081c4762a1bSJed Brown 1082c4762a1bSJed Brown test: 1083c4762a1bSJed Brown suffix: 3d_p1_neumann_0 1084c4762a1bSJed Brown requires: ctetgen 108530602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -bc_type neumann -dm_plex_boundary_label boundary -petscspace_degree 1 -snes_fd -show_initial -dm_plex_print_fem 1 -dm_view 1086c4762a1bSJed Brown 1087c4762a1bSJed Brown # Analytic variable coefficient 13-20 1088c4762a1bSJed Brown test: 1089c4762a1bSJed Brown suffix: 13 1090c4762a1bSJed Brown requires: triangle 109130602db0SMatthew G. Knepley args: -run_type test -variable_coefficient analytic -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1092c4762a1bSJed Brown test: 1093c4762a1bSJed Brown suffix: 14 1094c4762a1bSJed Brown requires: triangle 109530602db0SMatthew G. Knepley args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient analytic -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1096c4762a1bSJed Brown test: 1097c4762a1bSJed Brown suffix: 15 1098c4762a1bSJed Brown requires: triangle 109930602db0SMatthew G. Knepley args: -run_type test -variable_coefficient analytic -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1100c4762a1bSJed Brown test: 1101c4762a1bSJed Brown suffix: 16 1102c4762a1bSJed Brown requires: triangle 110330602db0SMatthew G. Knepley args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient analytic -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1104c4762a1bSJed Brown test: 1105c4762a1bSJed Brown suffix: 17 1106c4762a1bSJed Brown requires: ctetgen 110730602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -variable_coefficient analytic -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1108c4762a1bSJed Brown 1109c4762a1bSJed Brown test: 1110c4762a1bSJed Brown suffix: 18 1111c4762a1bSJed Brown requires: ctetgen 111230602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -dm_refine_volume_limit_pre 0.0125 -variable_coefficient analytic -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1113c4762a1bSJed Brown 1114c4762a1bSJed Brown test: 1115c4762a1bSJed Brown suffix: 19 1116c4762a1bSJed Brown requires: ctetgen 111730602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -variable_coefficient analytic -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1118c4762a1bSJed Brown 1119c4762a1bSJed Brown test: 1120c4762a1bSJed Brown suffix: 20 1121c4762a1bSJed Brown requires: ctetgen 112230602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -dm_refine_volume_limit_pre 0.0125 -variable_coefficient analytic -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1123c4762a1bSJed Brown 1124c4762a1bSJed Brown # P1 variable coefficient 21-28 1125c4762a1bSJed Brown test: 1126c4762a1bSJed Brown suffix: 21 1127c4762a1bSJed Brown requires: triangle 112830602db0SMatthew G. Knepley args: -run_type test -variable_coefficient field -petscspace_degree 1 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1129c4762a1bSJed Brown 1130c4762a1bSJed Brown test: 1131c4762a1bSJed Brown suffix: 22 1132c4762a1bSJed Brown requires: triangle 113330602db0SMatthew G. Knepley args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient field -petscspace_degree 1 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1134c4762a1bSJed Brown 1135c4762a1bSJed Brown test: 1136c4762a1bSJed Brown suffix: 23 1137c4762a1bSJed Brown requires: triangle 113830602db0SMatthew G. Knepley args: -run_type test -variable_coefficient field -petscspace_degree 2 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1139c4762a1bSJed Brown 1140c4762a1bSJed Brown test: 1141c4762a1bSJed Brown suffix: 24 1142c4762a1bSJed Brown requires: triangle 114330602db0SMatthew G. Knepley args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient field -petscspace_degree 2 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1144c4762a1bSJed Brown 1145c4762a1bSJed Brown test: 1146c4762a1bSJed Brown suffix: 25 1147c4762a1bSJed Brown requires: ctetgen 114830602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -variable_coefficient field -petscspace_degree 1 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1149c4762a1bSJed Brown 1150c4762a1bSJed Brown test: 1151c4762a1bSJed Brown suffix: 26 1152c4762a1bSJed Brown requires: ctetgen 115330602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -dm_refine_volume_limit_pre 0.0125 -variable_coefficient field -petscspace_degree 1 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1154c4762a1bSJed Brown 1155c4762a1bSJed Brown test: 1156c4762a1bSJed Brown suffix: 27 1157c4762a1bSJed Brown requires: ctetgen 115830602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -variable_coefficient field -petscspace_degree 2 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1159c4762a1bSJed Brown 1160c4762a1bSJed Brown test: 1161c4762a1bSJed Brown suffix: 28 1162c4762a1bSJed Brown requires: ctetgen 116330602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -dm_refine_volume_limit_pre 0.0125 -variable_coefficient field -petscspace_degree 2 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1164c4762a1bSJed Brown 1165c4762a1bSJed Brown # P0 variable coefficient 29-36 1166c4762a1bSJed Brown test: 1167c4762a1bSJed Brown suffix: 29 1168c4762a1bSJed Brown requires: triangle 116930602db0SMatthew G. Knepley args: -run_type test -variable_coefficient field -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1170c4762a1bSJed Brown 1171c4762a1bSJed Brown test: 1172c4762a1bSJed Brown suffix: 30 1173c4762a1bSJed Brown requires: triangle 117430602db0SMatthew G. Knepley args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient field -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1175c4762a1bSJed Brown 1176c4762a1bSJed Brown test: 1177c4762a1bSJed Brown suffix: 31 1178c4762a1bSJed Brown requires: triangle 117930602db0SMatthew G. Knepley args: -run_type test -variable_coefficient field -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1180c4762a1bSJed Brown 1181c4762a1bSJed Brown test: 1182c4762a1bSJed Brown requires: triangle 1183c4762a1bSJed Brown suffix: 32 118430602db0SMatthew G. Knepley args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient field -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1185c4762a1bSJed Brown 1186c4762a1bSJed Brown test: 1187c4762a1bSJed Brown requires: ctetgen 1188c4762a1bSJed Brown suffix: 33 118930602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -variable_coefficient field -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1190c4762a1bSJed Brown 1191c4762a1bSJed Brown test: 1192c4762a1bSJed Brown suffix: 34 1193c4762a1bSJed Brown requires: ctetgen 119430602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -dm_refine_volume_limit_pre 0.0125 -variable_coefficient field -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1195c4762a1bSJed Brown 1196c4762a1bSJed Brown test: 1197c4762a1bSJed Brown suffix: 35 1198c4762a1bSJed Brown requires: ctetgen 119930602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -variable_coefficient field -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1200c4762a1bSJed Brown 1201c4762a1bSJed Brown test: 1202c4762a1bSJed Brown suffix: 36 1203c4762a1bSJed Brown requires: ctetgen 120430602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -dm_refine_volume_limit_pre 0.0125 -variable_coefficient field -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1205c4762a1bSJed Brown 1206c4762a1bSJed Brown # Full solve 39-44 1207c4762a1bSJed Brown test: 1208c4762a1bSJed Brown suffix: 39 1209c4762a1bSJed Brown requires: triangle !single 121073f7197eSJed Brown args: -run_type full -dm_refine_volume_limit_pre 0.015625 -petscspace_degree 2 -pc_type gamg -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 -ksp_rtol 1.0e-10 -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason ::ascii_info_detail 1211c4762a1bSJed Brown test: 1212c4762a1bSJed Brown suffix: 40 1213c4762a1bSJed Brown requires: triangle !single 121430602db0SMatthew G. Knepley args: -run_type full -dm_refine_volume_limit_pre 0.015625 -variable_coefficient nonlinear -petscspace_degree 2 -pc_type svd -ksp_rtol 1.0e-10 -snes_monitor_short -snes_converged_reason ::ascii_info_detail 1215c4762a1bSJed Brown test: 1216c4762a1bSJed Brown suffix: 41 1217c4762a1bSJed Brown requires: triangle !single 121830602db0SMatthew G. Knepley args: -run_type full -dm_refine_volume_limit_pre 0.03125 -variable_coefficient nonlinear -petscspace_degree 1 -snes_type fas -snes_fas_levels 2 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -dm_refine_hierarchy 1 -snes_view -fas_levels_1_snes_type newtonls -fas_levels_1_pc_type svd -fas_levels_1_ksp_rtol 1.0e-10 -fas_levels_1_snes_monitor_short 1219c4762a1bSJed Brown test: 1220c4762a1bSJed Brown suffix: 42 1221c4762a1bSJed Brown requires: triangle !single 122230602db0SMatthew G. Knepley args: -run_type full -dm_refine_volume_limit_pre 0.0625 -variable_coefficient nonlinear -petscspace_degree 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -dm_refine_hierarchy 2 -snes_view -fas_levels_1_snes_type newtonls -fas_levels_1_pc_type svd -fas_levels_1_ksp_rtol 1.0e-10 -fas_levels_1_snes_monitor_short -fas_levels_2_snes_type newtonls -fas_levels_2_pc_type svd -fas_levels_2_ksp_rtol 1.0e-10 -fas_levels_2_snes_atol 1.0e-11 -fas_levels_2_snes_monitor_short 1223c4762a1bSJed Brown test: 1224c4762a1bSJed Brown suffix: 43 1225c4762a1bSJed Brown requires: triangle !single 1226c4762a1bSJed Brown nsize: 2 1227e600fa54SMatthew G. Knepley args: -run_type full -dm_refine_volume_limit_pre 0.03125 -variable_coefficient nonlinear -petscspace_degree 1 -snes_type fas -snes_fas_levels 2 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -dm_refine_hierarchy 1 -snes_view -fas_levels_1_snes_type newtonls -fas_levels_1_pc_type svd -fas_levels_1_ksp_rtol 1.0e-10 -fas_levels_1_snes_monitor_short 1228c4762a1bSJed Brown 1229c4762a1bSJed Brown test: 1230c4762a1bSJed Brown suffix: 44 1231c4762a1bSJed Brown requires: triangle !single 1232c4762a1bSJed Brown nsize: 2 1233e600fa54SMatthew G. Knepley args: -run_type full -dm_refine_volume_limit_pre 0.0625 -variable_coefficient nonlinear -petscspace_degree 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -dm_refine_hierarchy 2 -dm_plex_print_fem 0 -snes_view -fas_levels_1_snes_type newtonls -fas_levels_1_pc_type svd -fas_levels_1_ksp_rtol 1.0e-10 -fas_levels_1_snes_monitor_short -fas_levels_2_snes_type newtonls -fas_levels_2_pc_type svd -fas_levels_2_ksp_rtol 1.0e-10 -fas_levels_2_snes_atol 1.0e-11 -fas_levels_2_snes_monitor_short 1234c4762a1bSJed Brown 1235c4762a1bSJed Brown # These tests use a loose tolerance just to exercise the PtAP operations for MATIS and multiple PCBDDC setup calls inside PCMG 1236c4762a1bSJed Brown testset: 1237c4762a1bSJed Brown requires: triangle !single 1238c4762a1bSJed Brown nsize: 3 1239e600fa54SMatthew G. Knepley args: -run_type full -petscspace_degree 1 -dm_mat_type is -pc_type mg -pc_mg_levels 2 -mg_coarse_pc_type bddc -pc_mg_galerkin pmat -ksp_rtol 1.0e-2 -snes_converged_reason -dm_refine_hierarchy 2 -snes_max_it 4 1240c4762a1bSJed Brown test: 1241c4762a1bSJed Brown suffix: gmg_bddc 1242c4762a1bSJed Brown filter: sed -e "s/CONVERGED_FNORM_RELATIVE iterations 3/CONVERGED_FNORM_RELATIVE iterations 4/g" 1243c4762a1bSJed Brown args: -mg_levels_pc_type jacobi 1244c4762a1bSJed Brown test: 1245c4762a1bSJed Brown filter: sed -e "s/iterations [0-4]/iterations 4/g" 1246c4762a1bSJed Brown suffix: gmg_bddc_lev 1247c4762a1bSJed Brown args: -mg_levels_pc_type bddc 1248c4762a1bSJed Brown 1249c4762a1bSJed Brown # Restarting 1250c4762a1bSJed Brown testset: 1251c4762a1bSJed Brown suffix: restart 1252c4762a1bSJed Brown requires: hdf5 triangle !complex 125330602db0SMatthew G. Knepley args: -run_type test -bc_type dirichlet -petscspace_degree 1 1254c4762a1bSJed Brown test: 1255c4762a1bSJed Brown args: -dm_view hdf5:sol.h5 -vec_view hdf5:sol.h5::append 1256c4762a1bSJed Brown test: 1257cd7e8a5eSksagiyam args: -dm_plex_filename sol.h5 -dm_plex_name box -restart 1258c4762a1bSJed Brown 1259c4762a1bSJed Brown # Periodicity 1260c4762a1bSJed Brown test: 1261c4762a1bSJed Brown suffix: periodic_0 1262c4762a1bSJed Brown requires: triangle 126330602db0SMatthew G. Knepley args: -run_type full -bc_type dirichlet -petscspace_degree 1 -snes_converged_reason ::ascii_info_detail 1264c4762a1bSJed Brown 1265c4762a1bSJed Brown test: 1266c4762a1bSJed Brown requires: !complex 1267c4762a1bSJed Brown suffix: periodic_1 126830602db0SMatthew G. Knepley args: -quiet -run_type test -dm_plex_simplex 0 -dm_plex_box_faces 3,3 -dm_plex_box_bd periodic,periodic -vec_view vtk:test.vtu:vtk_vtu -petscspace_degree 1 -dm_refine 1 1269c4762a1bSJed Brown 1270c4762a1bSJed Brown # 2D serial P1 test with field bc 1271c4762a1bSJed Brown test: 1272c4762a1bSJed Brown suffix: field_bc_2d_p1_0 1273c4762a1bSJed Brown requires: triangle 127430602db0SMatthew G. Knepley args: -run_type test -bc_type dirichlet -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1275c4762a1bSJed Brown 1276c4762a1bSJed Brown test: 1277c4762a1bSJed Brown suffix: field_bc_2d_p1_1 1278c4762a1bSJed Brown requires: triangle 127930602db0SMatthew G. Knepley args: -run_type test -dm_refine 1 -bc_type dirichlet -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1280c4762a1bSJed Brown 1281c4762a1bSJed Brown test: 1282c4762a1bSJed Brown suffix: field_bc_2d_p1_neumann_0 1283c4762a1bSJed Brown requires: triangle 128430602db0SMatthew G. Knepley args: -run_type test -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1285c4762a1bSJed Brown 1286c4762a1bSJed Brown test: 1287c4762a1bSJed Brown suffix: field_bc_2d_p1_neumann_1 1288c4762a1bSJed Brown requires: triangle 128930602db0SMatthew G. Knepley args: -run_type test -dm_refine 1 -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1290c4762a1bSJed Brown 1291c4762a1bSJed Brown # 3D serial P1 test with field bc 1292c4762a1bSJed Brown test: 1293c4762a1bSJed Brown suffix: field_bc_3d_p1_0 1294c4762a1bSJed Brown requires: ctetgen 129530602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -bc_type dirichlet -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1296c4762a1bSJed Brown 1297c4762a1bSJed Brown test: 1298c4762a1bSJed Brown suffix: field_bc_3d_p1_1 1299c4762a1bSJed Brown requires: ctetgen 130030602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -dm_refine 1 -bc_type dirichlet -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1301c4762a1bSJed Brown 1302c4762a1bSJed Brown test: 1303c4762a1bSJed Brown suffix: field_bc_3d_p1_neumann_0 1304c4762a1bSJed Brown requires: ctetgen 130530602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1306c4762a1bSJed Brown 1307c4762a1bSJed Brown test: 1308c4762a1bSJed Brown suffix: field_bc_3d_p1_neumann_1 1309c4762a1bSJed Brown requires: ctetgen 131030602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -dm_refine 1 -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1311c4762a1bSJed Brown 1312c4762a1bSJed Brown # 2D serial P2 test with field bc 1313c4762a1bSJed Brown test: 1314c4762a1bSJed Brown suffix: field_bc_2d_p2_0 1315c4762a1bSJed Brown requires: triangle 131630602db0SMatthew G. Knepley args: -run_type test -bc_type dirichlet -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1317c4762a1bSJed Brown 1318c4762a1bSJed Brown test: 1319c4762a1bSJed Brown suffix: field_bc_2d_p2_1 1320c4762a1bSJed Brown requires: triangle 132130602db0SMatthew G. Knepley args: -run_type test -dm_refine 1 -bc_type dirichlet -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1322c4762a1bSJed Brown 1323c4762a1bSJed Brown test: 1324c4762a1bSJed Brown suffix: field_bc_2d_p2_neumann_0 1325c4762a1bSJed Brown requires: triangle 132630602db0SMatthew G. Knepley args: -run_type test -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1327c4762a1bSJed Brown 1328c4762a1bSJed Brown test: 1329c4762a1bSJed Brown suffix: field_bc_2d_p2_neumann_1 1330c4762a1bSJed Brown requires: triangle 133130602db0SMatthew G. Knepley args: -run_type test -dm_refine 1 -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1332c4762a1bSJed Brown 1333c4762a1bSJed Brown # 3D serial P2 test with field bc 1334c4762a1bSJed Brown test: 1335c4762a1bSJed Brown suffix: field_bc_3d_p2_0 1336c4762a1bSJed Brown requires: ctetgen 133730602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -bc_type dirichlet -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1338c4762a1bSJed Brown 1339c4762a1bSJed Brown test: 1340c4762a1bSJed Brown suffix: field_bc_3d_p2_1 1341c4762a1bSJed Brown requires: ctetgen 134230602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -dm_refine 1 -bc_type dirichlet -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1343c4762a1bSJed Brown 1344c4762a1bSJed Brown test: 1345c4762a1bSJed Brown suffix: field_bc_3d_p2_neumann_0 1346c4762a1bSJed Brown requires: ctetgen 134730602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1348c4762a1bSJed Brown 1349c4762a1bSJed Brown test: 1350c4762a1bSJed Brown suffix: field_bc_3d_p2_neumann_1 1351c4762a1bSJed Brown requires: ctetgen 135230602db0SMatthew G. Knepley args: -run_type test -dm_plex_dim 3 -dm_refine 1 -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1353c4762a1bSJed Brown 1354c4762a1bSJed Brown # Full solve simplex: Convergence 1355c4762a1bSJed Brown test: 13560fdc7489SMatthew Knepley suffix: 3d_p1_conv 1357c4762a1bSJed Brown requires: ctetgen 135830602db0SMatthew G. Knepley args: -run_type full -dm_plex_dim 3 -dm_refine 1 -bc_type dirichlet -petscspace_degree 1 \ 13590fdc7489SMatthew Knepley -snes_convergence_estimate -convest_num_refine 1 -pc_type lu 1360c4762a1bSJed Brown 1361c4762a1bSJed Brown # Full solve simplex: PCBDDC 1362c4762a1bSJed Brown test: 1363c4762a1bSJed Brown suffix: tri_bddc 1364c4762a1bSJed Brown requires: triangle !single 1365c4762a1bSJed Brown nsize: 5 1366e600fa54SMatthew G. Knepley args: -run_type full -petscpartitioner_type simple -dm_refine 2 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -dm_mat_type is -pc_type bddc -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0 1367c4762a1bSJed Brown 1368c4762a1bSJed Brown # Full solve simplex: PCBDDC 1369c4762a1bSJed Brown test: 1370c4762a1bSJed Brown suffix: tri_parmetis_bddc 1371c4762a1bSJed Brown requires: triangle !single parmetis 1372c4762a1bSJed Brown nsize: 4 1373e600fa54SMatthew G. Knepley args: -run_type full -petscpartitioner_type parmetis -dm_refine 2 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -dm_mat_type is -pc_type bddc -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0 1374c4762a1bSJed Brown 1375c4762a1bSJed Brown testset: 1376e600fa54SMatthew G. Knepley args: -run_type full -dm_plex_simplex 0 -dm_plex_box_faces 3,3 -petscpartitioner_type simple -dm_refine 2 -bc_type dirichlet -petscspace_degree 2 -dm_mat_type is -pc_type bddc -ksp_type gmres -snes_monitor_short -ksp_monitor_short -snes_view -petscspace_poly_tensor -pc_bddc_corner_selection -ksp_rtol 1.e-9 -pc_bddc_use_edges 0 1377c4762a1bSJed Brown nsize: 5 1378c4762a1bSJed Brown output_file: output/ex12_quad_bddc.out 1379c4762a1bSJed Brown filter: sed -e "s/aijcusparse/aij/g" -e "s/aijviennacl/aij/g" -e "s/factorization: cusparse/factorization: petsc/g" 1380c4762a1bSJed Brown test: 1381c4762a1bSJed Brown requires: !single 1382c4762a1bSJed Brown suffix: quad_bddc 1383c4762a1bSJed Brown test: 1384c4762a1bSJed Brown requires: !single cuda 1385c4762a1bSJed Brown suffix: quad_bddc_cuda 1386c4762a1bSJed Brown args: -matis_localmat_type aijcusparse -pc_bddc_dirichlet_pc_factor_mat_solver_type cusparse -pc_bddc_neumann_pc_factor_mat_solver_type cusparse 1387c4762a1bSJed Brown test: 1388c4762a1bSJed Brown requires: !single viennacl 1389c4762a1bSJed Brown suffix: quad_bddc_viennacl 1390c4762a1bSJed Brown args: -matis_localmat_type aijviennacl 1391c4762a1bSJed Brown 1392c4762a1bSJed Brown # Full solve simplex: ASM 1393c4762a1bSJed Brown test: 1394c4762a1bSJed Brown suffix: tri_q2q1_asm_lu 1395c4762a1bSJed Brown requires: triangle !single 139630602db0SMatthew G. Knepley args: -run_type full -dm_refine 3 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type asm -pc_asm_type restrict -pc_asm_blocks 4 -sub_pc_type lu -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0 1397c4762a1bSJed Brown 1398c4762a1bSJed Brown test: 1399c4762a1bSJed Brown suffix: tri_q2q1_msm_lu 1400c4762a1bSJed Brown requires: triangle !single 140130602db0SMatthew G. Knepley args: -run_type full -dm_refine 3 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type asm -pc_asm_type restrict -pc_asm_local_type multiplicative -pc_asm_blocks 4 -sub_pc_type lu -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0 1402c4762a1bSJed Brown 1403c4762a1bSJed Brown test: 1404c4762a1bSJed Brown suffix: tri_q2q1_asm_sor 1405c4762a1bSJed Brown requires: triangle !single 140630602db0SMatthew G. Knepley args: -run_type full -dm_refine 3 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type asm -pc_asm_type restrict -pc_asm_blocks 4 -sub_pc_type sor -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0 1407c4762a1bSJed Brown 1408c4762a1bSJed Brown test: 1409c4762a1bSJed Brown suffix: tri_q2q1_msm_sor 1410c4762a1bSJed Brown requires: triangle !single 141130602db0SMatthew G. Knepley args: -run_type full -dm_refine 3 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type asm -pc_asm_type restrict -pc_asm_local_type multiplicative -pc_asm_blocks 4 -sub_pc_type sor -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0 1412c4762a1bSJed Brown 1413c4762a1bSJed Brown # Full solve simplex: FAS 1414c4762a1bSJed Brown test: 1415c4762a1bSJed Brown suffix: fas_newton_0 1416c4762a1bSJed Brown requires: triangle !single 141730602db0SMatthew G. Knepley args: -run_type full -variable_coefficient nonlinear -petscspace_degree 1 -snes_type fas -snes_fas_levels 2 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -dm_refine_hierarchy 1 -snes_view -fas_levels_1_snes_type newtonls -fas_levels_1_pc_type svd -fas_levels_1_ksp_rtol 1.0e-10 -fas_levels_1_snes_monitor_short 1418c4762a1bSJed Brown 1419c4762a1bSJed Brown test: 1420c4762a1bSJed Brown suffix: fas_newton_1 1421c4762a1bSJed Brown requires: triangle !single 142230602db0SMatthew G. Knepley args: -run_type full -dm_refine_hierarchy 3 -petscspace_degree 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type lu -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_snes_linesearch_type basic -fas_levels_ksp_rtol 1.0e-10 -fas_levels_snes_monitor_short 1423c4ef839dSSatish Balay filter: sed -e "s/total number of linear solver iterations=14/total number of linear solver iterations=15/g" 1424c4762a1bSJed Brown 1425c4762a1bSJed Brown test: 1426c4762a1bSJed Brown suffix: fas_ngs_0 1427c4762a1bSJed Brown requires: triangle !single 142830602db0SMatthew G. Knepley args: -run_type full -variable_coefficient nonlinear -petscspace_degree 1 -snes_type fas -snes_fas_levels 2 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -dm_refine_hierarchy 1 -snes_view -fas_levels_1_snes_type ngs -fas_levels_1_snes_monitor_short 1429c4762a1bSJed Brown 1430071b71afSMatthew G. Knepley # These two tests are broken because DMPlexComputeInjectorFEM() only works for regularly refined meshes 1431c4762a1bSJed Brown test: 1432c4762a1bSJed Brown suffix: fas_newton_coarse_0 1433c4762a1bSJed Brown requires: pragmatic triangle 1434c4762a1bSJed Brown TODO: broken 1435071b71afSMatthew G. Knepley args: -run_type full -variable_coefficient nonlinear -petscspace_degree 1 \ 143634b6e994SJoe Wallwork -dm_refine 2 -dm_coarsen_hierarchy 1 -dm_plex_hash_location -dm_adaptor pragmatic \ 1437071b71afSMatthew G. Knepley -snes_type fas -snes_fas_levels 2 -snes_converged_reason ::ascii_info_detail -snes_monitor_short -snes_view \ 1438071b71afSMatthew G. Knepley -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -fas_coarse_snes_linesearch_type basic \ 1439071b71afSMatthew G. Knepley -fas_levels_1_snes_type newtonls -fas_levels_1_pc_type svd -fas_levels_1_ksp_rtol 1.0e-10 -fas_levels_1_snes_monitor_short 1440c4762a1bSJed Brown 1441c4762a1bSJed Brown test: 1442c4762a1bSJed Brown suffix: mg_newton_coarse_0 1443c4762a1bSJed Brown requires: triangle pragmatic 1444c4762a1bSJed Brown TODO: broken 1445071b71afSMatthew G. Knepley args: -run_type full -petscspace_degree 1 \ 144634b6e994SJoe Wallwork -dm_refine 3 -dm_coarsen_hierarchy 3 -dm_plex_hash_location -dm_adaptor pragmatic \ 1447071b71afSMatthew G. Knepley -snes_atol 1.0e-8 -snes_rtol 0.0 -snes_monitor_short -snes_converged_reason ::ascii_info_detail -snes_view \ 1448071b71afSMatthew G. Knepley -ksp_type richardson -ksp_atol 1.0e-8 -ksp_rtol 0.0 -ksp_norm_type unpreconditioned -ksp_monitor_true_residual \ 1449071b71afSMatthew G. Knepley -pc_type mg -pc_mg_levels 4 \ 1450071b71afSMatthew G. Knepley -mg_levels_ksp_type gmres -mg_levels_pc_type ilu -mg_levels_ksp_max_it 10 1451c4762a1bSJed Brown 1452c4762a1bSJed Brown # Full solve tensor 1453c4762a1bSJed Brown test: 1454c4762a1bSJed Brown suffix: tensor_plex_2d 145530602db0SMatthew G. Knepley args: -run_type test -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -dm_refine_hierarchy 2 1456c4762a1bSJed Brown 1457c4762a1bSJed Brown test: 1458c4762a1bSJed Brown suffix: tensor_p4est_2d 1459c4762a1bSJed Brown requires: p4est 146030602db0SMatthew G. Knepley args: -run_type test -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -dm_forest_initial_refinement 2 -dm_forest_minimum_refinement 0 -dm_plex_convert_type p4est 1461c4762a1bSJed Brown 1462c4762a1bSJed Brown test: 1463c4762a1bSJed Brown suffix: tensor_plex_3d 146430602db0SMatthew G. Knepley args: -run_type test -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -dm_plex_dim 3 -dm_refine_hierarchy 1 -dm_plex_box_faces 2,2,2 1465c4762a1bSJed Brown 1466c4762a1bSJed Brown test: 1467c4762a1bSJed Brown suffix: tensor_p4est_3d 1468c4762a1bSJed Brown requires: p4est 146930602db0SMatthew G. Knepley args: -run_type test -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -dm_forest_initial_refinement 1 -dm_forest_minimum_refinement 0 -dm_plex_dim 3 -dm_plex_convert_type p8est -dm_plex_box_faces 2,2,2 1470c4762a1bSJed Brown 1471c4762a1bSJed Brown test: 1472c4762a1bSJed Brown suffix: p4est_test_q2_conformal_serial 1473c4762a1bSJed Brown requires: p4est 147430602db0SMatthew G. Knepley args: -run_type test -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 1475c4762a1bSJed Brown 1476c4762a1bSJed Brown test: 1477c4762a1bSJed Brown suffix: p4est_test_q2_conformal_parallel 1478c4762a1bSJed Brown requires: p4est 1479c4762a1bSJed Brown nsize: 7 1480e600fa54SMatthew G. Knepley args: -run_type test -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -petscpartitioner_type simple 1481c4762a1bSJed Brown 1482c4762a1bSJed Brown test: 1483c4762a1bSJed Brown suffix: p4est_test_q2_conformal_parallel_parmetis 1484c4762a1bSJed Brown requires: parmetis p4est 1485c4762a1bSJed Brown nsize: 4 1486e600fa54SMatthew G. Knepley args: -run_type test -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -petscpartitioner_type parmetis 1487c4762a1bSJed Brown 1488c4762a1bSJed Brown test: 1489c4762a1bSJed Brown suffix: p4est_test_q2_nonconformal_serial 1490c4762a1bSJed Brown requires: p4est 1491c4762a1bSJed Brown filter: grep -v "CG or CGNE: variant" 149230602db0SMatthew G. Knepley args: -run_type test -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash 1493c4762a1bSJed Brown 1494c4762a1bSJed Brown test: 1495c4762a1bSJed Brown suffix: p4est_test_q2_nonconformal_parallel 1496c4762a1bSJed Brown requires: p4est 1497c4762a1bSJed Brown filter: grep -v "CG or CGNE: variant" 1498c4762a1bSJed Brown nsize: 7 1499e600fa54SMatthew G. Knepley args: -run_type test -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple 1500c4762a1bSJed Brown 1501c4762a1bSJed Brown test: 1502c4762a1bSJed Brown suffix: p4est_test_q2_nonconformal_parallel_parmetis 1503c4762a1bSJed Brown requires: parmetis p4est 1504c4762a1bSJed Brown nsize: 4 1505e600fa54SMatthew G. Knepley args: -run_type test -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type parmetis 1506c4762a1bSJed Brown 1507c4762a1bSJed Brown test: 1508c4762a1bSJed Brown suffix: p4est_exact_q2_conformal_serial 1509c4762a1bSJed Brown requires: p4est !single !complex !__float128 151030602db0SMatthew G. Knepley args: -run_type exact -petscspace_degree 2 -fas_levels_snes_atol 1.e-10 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 1511c4762a1bSJed Brown 1512c4762a1bSJed Brown test: 1513c4762a1bSJed Brown suffix: p4est_exact_q2_conformal_parallel 1514c4762a1bSJed Brown requires: p4est !single !complex !__float128 1515c4762a1bSJed Brown nsize: 4 1516e600fa54SMatthew G. Knepley args: -run_type exact -petscspace_degree 2 -fas_levels_snes_atol 1.e-10 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 1517c4762a1bSJed Brown 1518c4762a1bSJed Brown test: 1519c4762a1bSJed Brown suffix: p4est_exact_q2_conformal_parallel_parmetis 1520c4762a1bSJed Brown requires: parmetis p4est !single 1521c4762a1bSJed Brown nsize: 4 1522e600fa54SMatthew G. Knepley args: -run_type exact -petscspace_degree 2 -fas_levels_snes_atol 1.e-10 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -petscpartitioner_type parmetis 1523c4762a1bSJed Brown 1524c4762a1bSJed Brown test: 1525c4762a1bSJed Brown suffix: p4est_exact_q2_nonconformal_serial 1526c4762a1bSJed Brown requires: p4est 152730602db0SMatthew G. Knepley args: -run_type exact -petscspace_degree 2 -fas_levels_snes_atol 1.e-10 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash 1528c4762a1bSJed Brown 1529c4762a1bSJed Brown test: 1530c4762a1bSJed Brown suffix: p4est_exact_q2_nonconformal_parallel 1531c4762a1bSJed Brown requires: p4est 1532c4762a1bSJed Brown nsize: 7 1533e600fa54SMatthew G. Knepley args: -run_type exact -petscspace_degree 2 -fas_levels_snes_atol 1.e-10 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple 1534c4762a1bSJed Brown 1535c4762a1bSJed Brown test: 1536c4762a1bSJed Brown suffix: p4est_exact_q2_nonconformal_parallel_parmetis 1537c4762a1bSJed Brown requires: parmetis p4est 1538c4762a1bSJed Brown nsize: 4 1539e600fa54SMatthew G. Knepley args: -run_type exact -petscspace_degree 2 -fas_levels_snes_atol 1.e-10 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type parmetis 1540c4762a1bSJed Brown 1541c4762a1bSJed Brown test: 1542c4762a1bSJed Brown suffix: p4est_full_q2_nonconformal_serial 1543c4762a1bSJed Brown requires: p4est !single 1544c4762a1bSJed Brown filter: grep -v "variant HERMITIAN" 154530602db0SMatthew G. Knepley args: -run_type full -petscspace_degree 2 -snes_max_it 20 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type jacobi -fas_coarse_ksp_type cg -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type jacobi -fas_levels_ksp_type cg -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash 1546c4762a1bSJed Brown 1547c4762a1bSJed Brown test: 1548c4762a1bSJed Brown suffix: p4est_full_q2_nonconformal_parallel 1549c4762a1bSJed Brown requires: p4est !single 1550c4762a1bSJed Brown filter: grep -v "variant HERMITIAN" 1551c4762a1bSJed Brown nsize: 7 1552e600fa54SMatthew G. Knepley args: -run_type full -petscspace_degree 2 -snes_max_it 20 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type jacobi -fas_coarse_ksp_type cg -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type jacobi -fas_levels_ksp_type cg -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple 1553c4762a1bSJed Brown 1554c4762a1bSJed Brown test: 1555c4762a1bSJed Brown suffix: p4est_full_q2_nonconformal_parallel_bddcfas 1556c4762a1bSJed Brown requires: p4est !single 1557c4762a1bSJed Brown filter: grep -v "variant HERMITIAN" 1558c4762a1bSJed Brown nsize: 7 1559e600fa54SMatthew G. Knepley args: -run_type full -petscspace_degree 2 -snes_max_it 20 -snes_type fas -snes_fas_levels 3 -dm_mat_type is -fas_coarse_pc_type bddc -fas_coarse_ksp_type cg -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type bddc -fas_levels_ksp_type cg -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple 1560c4762a1bSJed Brown 1561c4762a1bSJed Brown test: 1562c4762a1bSJed Brown suffix: p4est_full_q2_nonconformal_parallel_bddc 1563c4762a1bSJed Brown requires: p4est !single 1564c4762a1bSJed Brown filter: grep -v "variant HERMITIAN" 1565c4762a1bSJed Brown nsize: 7 1566e600fa54SMatthew G. Knepley args: -run_type full -petscspace_degree 2 -snes_max_it 20 -snes_type newtonls -dm_mat_type is -pc_type bddc -ksp_type cg -snes_monitor_short -snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple 1567c4762a1bSJed Brown 1568c4762a1bSJed Brown test: 1569c4762a1bSJed Brown TODO: broken 1570c4762a1bSJed Brown suffix: p4est_fas_q2_conformal_serial 1571c4762a1bSJed Brown requires: p4est !complex !__float128 157230602db0SMatthew G. Knepley args: -run_type full -variable_coefficient nonlinear -petscspace_degree 2 -snes_max_it 20 -snes_type fas -snes_fas_levels 3 -pc_type jacobi -ksp_type gmres -fas_coarse_pc_type svd -fas_coarse_ksp_type gmres -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type svd -fas_levels_ksp_type gmres -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_refine_hierarchy 3 1573c4762a1bSJed Brown 1574c4762a1bSJed Brown test: 1575c4762a1bSJed Brown TODO: broken 1576c4762a1bSJed Brown suffix: p4est_fas_q2_nonconformal_serial 1577c4762a1bSJed Brown requires: p4est 157830602db0SMatthew G. Knepley args: -run_type full -variable_coefficient nonlinear -petscspace_degree 2 -snes_max_it 20 -snes_type fas -snes_fas_levels 3 -pc_type jacobi -ksp_type gmres -fas_coarse_pc_type jacobi -fas_coarse_ksp_type gmres -fas_coarse_ksp_monitor_true_residual -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type jacobi -fas_levels_ksp_type gmres -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash 1579c4762a1bSJed Brown 1580c4762a1bSJed Brown test: 1581c4762a1bSJed Brown suffix: fas_newton_0_p4est 1582c4762a1bSJed Brown requires: p4est !single !__float128 158330602db0SMatthew G. Knepley args: -run_type full -variable_coefficient nonlinear -petscspace_degree 1 -snes_type fas -snes_fas_levels 2 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_1_snes_type newtonls -fas_levels_1_pc_type svd -fas_levels_1_ksp_rtol 1.0e-10 -fas_levels_1_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash 1584c4762a1bSJed Brown 1585c4762a1bSJed Brown # Full solve simplicial AMR 1586c4762a1bSJed Brown test: 1587ab5a7ff4SJoe Wallwork suffix: tri_p1_adapt_init_pragmatic 1588c4762a1bSJed Brown requires: pragmatic 15898d1b37daSJoe Wallwork args: -run_type exact -dm_refine 5 -bc_type dirichlet -petscspace_degree 1 -variable_coefficient ball -snes_converged_reason ::ascii_info_detail -pc_type lu -snes_adapt_initial 1 -adaptor_target_num 4000 -dm_plex_metric_h_max 0.5 -dm_adaptor pragmatic 1590c4762a1bSJed Brown 1591c4762a1bSJed Brown test: 15920383c1e7SJoe Wallwork suffix: tri_p2_adapt_init_pragmatic 15930383c1e7SJoe Wallwork requires: pragmatic 15940383c1e7SJoe Wallwork args: -run_type exact -dm_refine 5 -bc_type dirichlet -petscspace_degree 2 -variable_coefficient ball -snes_converged_reason ::ascii_info_detail -pc_type lu -snes_adapt_initial 1 -adaptor_target_num 4000 -dm_plex_metric_h_max 0.5 -dm_adaptor pragmatic 15950383c1e7SJoe Wallwork 15960383c1e7SJoe Wallwork test: 1597ab5a7ff4SJoe Wallwork suffix: tri_p1_adapt_init_mmg 1598ab5a7ff4SJoe Wallwork requires: mmg 15998d1b37daSJoe Wallwork args: -run_type exact -dm_refine 5 -bc_type dirichlet -petscspace_degree 1 -variable_coefficient ball -snes_converged_reason ::ascii_info_detail -pc_type lu -snes_adapt_initial 1 -adaptor_target_num 4000 -dm_plex_metric_h_max 0.5 -dm_adaptor mmg 1600c4762a1bSJed Brown 1601c4762a1bSJed Brown test: 16020383c1e7SJoe Wallwork suffix: tri_p2_adapt_init_mmg 16030383c1e7SJoe Wallwork requires: mmg 16040383c1e7SJoe Wallwork args: -run_type exact -dm_refine 5 -bc_type dirichlet -petscspace_degree 2 -variable_coefficient ball -snes_converged_reason ::ascii_info_detail -pc_type lu -snes_adapt_initial 1 -adaptor_target_num 4000 -dm_plex_metric_h_max 0.5 -dm_adaptor mmg 16050383c1e7SJoe Wallwork 16060383c1e7SJoe Wallwork test: 1607ab5a7ff4SJoe Wallwork suffix: tri_p1_adapt_seq_pragmatic 1608c4762a1bSJed Brown requires: pragmatic 16098d1b37daSJoe Wallwork args: -run_type exact -dm_refine 5 -bc_type dirichlet -petscspace_degree 1 -variable_coefficient ball -snes_converged_reason ::ascii_info_detail -pc_type lu -snes_adapt_sequence 2 -adaptor_target_num 4000 -dm_plex_metric_h_max 0.5 -dm_adaptor pragmatic 1610ab5a7ff4SJoe Wallwork 1611ab5a7ff4SJoe Wallwork test: 16120383c1e7SJoe Wallwork suffix: tri_p2_adapt_seq_pragmatic 16130383c1e7SJoe Wallwork requires: pragmatic 16140383c1e7SJoe Wallwork args: -run_type exact -dm_refine 5 -bc_type dirichlet -petscspace_degree 2 -variable_coefficient ball -snes_converged_reason ::ascii_info_detail -pc_type lu -snes_adapt_sequence 2 -adaptor_target_num 4000 -dm_plex_metric_h_max 0.5 -dm_adaptor pragmatic 16150383c1e7SJoe Wallwork 16160383c1e7SJoe Wallwork test: 1617ab5a7ff4SJoe Wallwork suffix: tri_p1_adapt_seq_mmg 1618ab5a7ff4SJoe Wallwork requires: mmg 16198d1b37daSJoe Wallwork args: -run_type exact -dm_refine 5 -bc_type dirichlet -petscspace_degree 1 -variable_coefficient ball -snes_converged_reason ::ascii_info_detail -pc_type lu -snes_adapt_sequence 2 -adaptor_target_num 4000 -dm_plex_metric_h_max 0.5 -dm_adaptor mmg 1620ab5a7ff4SJoe Wallwork 1621ab5a7ff4SJoe Wallwork test: 16220383c1e7SJoe Wallwork suffix: tri_p2_adapt_seq_mmg 16230383c1e7SJoe Wallwork requires: mmg 16240383c1e7SJoe Wallwork args: -run_type exact -dm_refine 5 -bc_type dirichlet -petscspace_degree 2 -variable_coefficient ball -snes_converged_reason ::ascii_info_detail -pc_type lu -snes_adapt_sequence 2 -adaptor_target_num 4000 -dm_plex_metric_h_max 0.5 -dm_adaptor mmg 16250383c1e7SJoe Wallwork 16260383c1e7SJoe Wallwork test: 1627ab5a7ff4SJoe Wallwork suffix: tri_p1_adapt_analytic_pragmatic 1628ab5a7ff4SJoe Wallwork requires: pragmatic 1629ab5a7ff4SJoe Wallwork args: -run_type exact -dm_refine 3 -bc_type dirichlet -petscspace_degree 1 -variable_coefficient cross -snes_adapt_initial 4 -adaptor_target_num 500 -adaptor_monitor -dm_plex_metric_h_min 0.0001 -dm_plex_metric_h_max 0.05 -dm_adaptor pragmatic 1630ab5a7ff4SJoe Wallwork 1631ab5a7ff4SJoe Wallwork test: 16320383c1e7SJoe Wallwork suffix: tri_p2_adapt_analytic_pragmatic 16330383c1e7SJoe Wallwork requires: pragmatic 16340383c1e7SJoe Wallwork args: -run_type exact -dm_refine 3 -bc_type dirichlet -petscspace_degree 2 -variable_coefficient cross -snes_adapt_initial 4 -adaptor_target_num 500 -adaptor_monitor -dm_plex_metric_h_min 0.0001 -dm_plex_metric_h_max 0.05 -dm_adaptor pragmatic 16350383c1e7SJoe Wallwork 16360383c1e7SJoe Wallwork test: 1637ab5a7ff4SJoe Wallwork suffix: tri_p1_adapt_analytic_mmg 1638ab5a7ff4SJoe Wallwork requires: mmg 16398d1b37daSJoe Wallwork args: -run_type exact -dm_refine 3 -bc_type dirichlet -petscspace_degree 1 -variable_coefficient cross -snes_adapt_initial 4 -adaptor_target_num 500 -adaptor_monitor -dm_plex_metric_h_max 0.5 -dm_adaptor mmg 1640c4762a1bSJed Brown 1641b8d0c900SJoe Wallwork test: 16420383c1e7SJoe Wallwork suffix: tri_p2_adapt_analytic_mmg 16430383c1e7SJoe Wallwork requires: mmg 16440383c1e7SJoe Wallwork args: -run_type exact -dm_refine 3 -bc_type dirichlet -petscspace_degree 2 -variable_coefficient cross -snes_adapt_initial 4 -adaptor_target_num 500 -adaptor_monitor -dm_plex_metric_h_max 0.5 -dm_adaptor mmg 16450383c1e7SJoe Wallwork 16460383c1e7SJoe Wallwork test: 1647b8d0c900SJoe Wallwork suffix: tri_p1_adapt_uniform_pragmatic 1648b8d0c900SJoe Wallwork requires: pragmatic tetgen 1649dc13bed2SJoe Wallwork nsize: 2 1650e600fa54SMatthew G. Knepley args: -run_type full -dm_plex_box_faces 8,8,8 -bc_type dirichlet -petscspace_degree 1 -variable_coefficient none -snes_converged_reason ::ascii_info_detail -ksp_type cg -pc_type sor -snes_adapt_sequence 3 -adaptor_target_num 400 -dm_plex_metric_h_max 0.5 -dm_plex_dim 3 -dm_adaptor pragmatic 1651b8d0c900SJoe Wallwork timeoutfactor: 2 1652b8d0c900SJoe Wallwork 1653b8d0c900SJoe Wallwork test: 16540383c1e7SJoe Wallwork suffix: tri_p2_adapt_uniform_pragmatic 16550383c1e7SJoe Wallwork requires: pragmatic tetgen 1656dc13bed2SJoe Wallwork nsize: 2 1657e600fa54SMatthew G. Knepley args: -run_type full -dm_plex_box_faces 8,8,8 -bc_type dirichlet -petscspace_degree 2 -variable_coefficient none -snes_converged_reason ::ascii_info_detail -ksp_type cg -pc_type sor -snes_adapt_sequence 1 -adaptor_target_num 400 -dm_plex_metric_h_max 0.5 -dm_plex_dim 3 -dm_adaptor pragmatic 16580383c1e7SJoe Wallwork timeoutfactor: 1 16590383c1e7SJoe Wallwork 16600383c1e7SJoe Wallwork test: 1661b8d0c900SJoe Wallwork suffix: tri_p1_adapt_uniform_mmg 1662b8d0c900SJoe Wallwork requires: mmg tetgen 16638d1b37daSJoe Wallwork args: -run_type full -dm_plex_box_faces 4,4,4 -bc_type dirichlet -petscspace_degree 1 -variable_coefficient none -snes_converged_reason ::ascii_info_detail -ksp_type cg -pc_type sor -snes_adapt_sequence 3 -adaptor_target_num 400 -dm_plex_metric_h_max 0.5 -dm_plex_dim 3 -dm_adaptor mmg 1664b8d0c900SJoe Wallwork timeoutfactor: 2 1665b8d0c900SJoe Wallwork 1666b8d0c900SJoe Wallwork test: 16670383c1e7SJoe Wallwork suffix: tri_p2_adapt_uniform_mmg 16680383c1e7SJoe Wallwork requires: mmg tetgen 16690383c1e7SJoe Wallwork args: -run_type full -dm_plex_box_faces 4,4,4 -bc_type dirichlet -petscspace_degree 2 -variable_coefficient none -snes_converged_reason ::ascii_info_detail -ksp_type cg -pc_type sor -snes_adapt_sequence 1 -adaptor_target_num 400 -dm_plex_metric_h_max 0.5 -dm_plex_dim 3 -dm_adaptor mmg 16700383c1e7SJoe Wallwork timeoutfactor: 1 16710383c1e7SJoe Wallwork 16720383c1e7SJoe Wallwork test: 1673b8d0c900SJoe Wallwork suffix: tri_p1_adapt_uniform_parmmg 1674b8d0c900SJoe Wallwork requires: parmmg tetgen 1675dc13bed2SJoe Wallwork nsize: 2 1676e600fa54SMatthew G. Knepley args: -run_type full -dm_plex_box_faces 8,8,8 -bc_type dirichlet -petscspace_degree 1 -variable_coefficient none -snes_converged_reason ::ascii_info_detail -ksp_type cg -pc_type sor -snes_adapt_sequence 3 -adaptor_target_num 400 -dm_plex_metric_h_max 0.5 -dm_plex_dim 3 -dm_adaptor parmmg 1677b8d0c900SJoe Wallwork timeoutfactor: 2 1678b8d0c900SJoe Wallwork 16790383c1e7SJoe Wallwork test: 16800383c1e7SJoe Wallwork suffix: tri_p2_adapt_uniform_parmmg 16810383c1e7SJoe Wallwork requires: parmmg tetgen 1682dc13bed2SJoe Wallwork nsize: 2 1683e600fa54SMatthew G. Knepley args: -run_type full -dm_plex_box_faces 8,8,8 -bc_type dirichlet -petscspace_degree 2 -variable_coefficient none -snes_converged_reason ::ascii_info_detail -ksp_type cg -pc_type sor -snes_adapt_sequence 1 -adaptor_target_num 400 -dm_plex_metric_h_max 0.5 -dm_plex_dim 3 -dm_adaptor parmmg 16840383c1e7SJoe Wallwork timeoutfactor: 1 16850383c1e7SJoe Wallwork 1686c4762a1bSJed Brown # Full solve tensor AMR 1687c4762a1bSJed Brown test: 1688c4762a1bSJed Brown suffix: quad_q1_adapt_0 1689c4762a1bSJed Brown requires: p4est 16908d1b37daSJoe Wallwork args: -run_type exact -dm_plex_simplex 0 -dm_plex_convert_type p4est -bc_type dirichlet -petscspace_degree 1 -variable_coefficient ball -snes_converged_reason ::ascii_info_detail -pc_type lu -dm_forest_initial_refinement 4 -snes_adapt_initial 1 -dm_view 1691c4762a1bSJed Brown filter: grep -v DM_ 1692c4762a1bSJed Brown 1693c4762a1bSJed Brown test: 1694c4762a1bSJed Brown suffix: amr_0 1695c4762a1bSJed Brown nsize: 5 1696e600fa54SMatthew G. Knepley args: -run_type test -petscpartitioner_type simple -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -dm_refine 1 1697c4762a1bSJed Brown 1698c4762a1bSJed Brown test: 1699c4762a1bSJed Brown suffix: amr_1 1700c4762a1bSJed Brown requires: p4est !complex 170130602db0SMatthew G. Knepley args: -run_type test -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -dm_plex_convert_type p4est -dm_p4est_refine_pattern center -dm_forest_maximum_refinement 5 -dm_view vtk:amr.vtu:vtk_vtu -vec_view vtk:amr.vtu:vtk_vtu:append 1702c4762a1bSJed Brown 1703c4762a1bSJed Brown test: 1704c4762a1bSJed Brown suffix: p4est_solve_bddc 1705c4762a1bSJed Brown requires: p4est !complex 1706e600fa54SMatthew G. Knepley args: -run_type full -variable_coefficient nonlinear -nonzero_initial_guess 1 -petscspace_degree 2 -snes_max_it 20 -snes_type newtonls -dm_mat_type is -pc_type bddc -ksp_type cg -snes_monitor_short -ksp_monitor -snes_linesearch_type bt -snes_converged_reason -snes_view -dm_plex_simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple -pc_bddc_detect_disconnected 1707c4762a1bSJed Brown nsize: 4 1708c4762a1bSJed Brown 1709c4762a1bSJed Brown test: 1710c4762a1bSJed Brown suffix: p4est_solve_fas 1711c4762a1bSJed Brown requires: p4est 1712e600fa54SMatthew G. Knepley args: -run_type full -variable_coefficient nonlinear -nonzero_initial_guess 1 -petscspace_degree 2 -snes_max_it 10 -snes_type fas -snes_linesearch_type bt -snes_fas_levels 3 -fas_coarse_snes_type newtonls -fas_coarse_snes_linesearch_type basic -fas_coarse_ksp_type cg -fas_coarse_pc_type jacobi -fas_coarse_snes_monitor_short -fas_levels_snes_max_it 4 -fas_levels_snes_type newtonls -fas_levels_snes_linesearch_type bt -fas_levels_ksp_type cg -fas_levels_pc_type jacobi -fas_levels_snes_monitor_short -fas_levels_cycle_snes_linesearch_type bt -snes_monitor_short -snes_converged_reason -snes_view -dm_plex_simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash 1713c4762a1bSJed Brown nsize: 4 1714c4762a1bSJed Brown TODO: identical machine two runs produce slightly different solver trackers 1715c4762a1bSJed Brown 1716c4762a1bSJed Brown test: 1717c4762a1bSJed Brown suffix: p4est_convergence_test_1 1718c4762a1bSJed Brown requires: p4est 1719e600fa54SMatthew G. Knepley args: -quiet -run_type test -petscspace_degree 1 -dm_plex_simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 2 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash 1720c4762a1bSJed Brown nsize: 4 1721c4762a1bSJed Brown 1722c4762a1bSJed Brown test: 1723c4762a1bSJed Brown suffix: p4est_convergence_test_2 1724c4762a1bSJed Brown requires: p4est 172530602db0SMatthew G. Knepley args: -quiet -run_type test -petscspace_degree 1 -dm_plex_simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 3 -dm_forest_initial_refinement 3 -dm_forest_maximum_refinement 5 -dm_p4est_refine_pattern hash 1726c4762a1bSJed Brown 1727c4762a1bSJed Brown test: 1728c4762a1bSJed Brown suffix: p4est_convergence_test_3 1729c4762a1bSJed Brown requires: p4est 173030602db0SMatthew G. Knepley args: -quiet -run_type test -petscspace_degree 1 -dm_plex_simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 4 -dm_forest_initial_refinement 4 -dm_forest_maximum_refinement 6 -dm_p4est_refine_pattern hash 1731c4762a1bSJed Brown 1732c4762a1bSJed Brown test: 1733c4762a1bSJed Brown suffix: p4est_convergence_test_4 1734c4762a1bSJed Brown requires: p4est 173530602db0SMatthew G. Knepley args: -quiet -run_type test -petscspace_degree 1 -dm_plex_simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 5 -dm_forest_initial_refinement 5 -dm_forest_maximum_refinement 7 -dm_p4est_refine_pattern hash 1736c4762a1bSJed Brown timeoutfactor: 5 1737c4762a1bSJed Brown 1738c4762a1bSJed Brown # Serial tests with GLVis visualization 1739c4762a1bSJed Brown test: 1740c4762a1bSJed Brown suffix: glvis_2d_tet_p1 174130602db0SMatthew G. Knepley args: -quiet -run_type test -bc_type dirichlet -petscspace_degree 1 -vec_view glvis: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_boundary_label marker -dm_plex_gmsh_periodic 0 -dm_coord_space 0 1742c4762a1bSJed Brown test: 1743c4762a1bSJed Brown suffix: glvis_2d_tet_p2 174430602db0SMatthew G. Knepley args: -quiet -run_type test -bc_type dirichlet -petscspace_degree 2 -vec_view glvis: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_boundary_label marker -dm_plex_gmsh_periodic 0 -dm_coord_space 0 1745c4762a1bSJed Brown test: 1746c4762a1bSJed Brown suffix: glvis_2d_hex_p1 174730602db0SMatthew G. Knepley args: -quiet -run_type test -bc_type dirichlet -petscspace_degree 1 -vec_view glvis: -dm_plex_simplex 0 -dm_refine 1 -dm_coord_space 0 1748c4762a1bSJed Brown test: 1749c4762a1bSJed Brown suffix: glvis_2d_hex_p2 175030602db0SMatthew G. Knepley args: -quiet -run_type test -bc_type dirichlet -petscspace_degree 2 -vec_view glvis: -dm_plex_simplex 0 -dm_refine 1 -dm_coord_space 0 1751c4762a1bSJed Brown test: 1752c4762a1bSJed Brown suffix: glvis_2d_hex_p2_p4est 1753c4762a1bSJed Brown requires: p4est 175430602db0SMatthew G. Knepley args: -quiet -run_type test -bc_type dirichlet -petscspace_degree 2 -vec_view glvis: -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 1 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -viewer_glvis_dm_plex_enable_ncmesh 1755c4762a1bSJed Brown test: 1756c4762a1bSJed Brown suffix: glvis_2d_tet_p0 175730602db0SMatthew G. Knepley args: -run_type exact -guess_vec_view glvis: -nonzero_initial_guess 1 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_boundary_label marker -petscspace_degree 0 -dm_coord_space 0 1758c4762a1bSJed Brown test: 1759c4762a1bSJed Brown suffix: glvis_2d_hex_p0 176030602db0SMatthew G. Knepley args: -run_type exact -guess_vec_view glvis: -nonzero_initial_guess 1 -dm_plex_box_faces 5,7 -dm_plex_simplex 0 -petscspace_degree 0 -dm_coord_space 0 1761c4762a1bSJed Brown 1762c4762a1bSJed Brown # PCHPDDM tests 1763c4762a1bSJed Brown testset: 1764c4762a1bSJed Brown nsize: 4 1765dfd57a17SPierre Jolivet requires: hpddm slepc !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES) 1766e600fa54SMatthew G. Knepley args: -run_type test -run_test_check_ksp -quiet -petscspace_degree 1 -petscpartitioner_type simple -bc_type none -dm_plex_simplex 0 -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_nev 2 -pc_hpddm_coarse_p 1 -pc_hpddm_coarse_pc_type svd -ksp_rtol 1.e-10 -pc_hpddm_levels_1_st_pc_factor_shift_type INBLOCKS -ksp_converged_reason 1767c4762a1bSJed Brown test: 1768c4762a1bSJed Brown suffix: quad_singular_hpddm 176930602db0SMatthew G. Knepley args: -dm_plex_box_faces 6,7 1770c4762a1bSJed Brown test: 1771c4762a1bSJed Brown requires: p4est 1772c4762a1bSJed Brown suffix: p4est_singular_2d_hpddm 1773c4762a1bSJed Brown args: -dm_plex_convert_type p4est -dm_forest_minimum_refinement 1 -dm_forest_initial_refinement 3 -dm_forest_maximum_refinement 3 1774c4762a1bSJed Brown test: 1775c4762a1bSJed Brown requires: p4est 1776c4762a1bSJed Brown suffix: p4est_nc_singular_2d_hpddm 1777c4762a1bSJed Brown args: -dm_plex_convert_type p4est -dm_forest_minimum_refinement 1 -dm_forest_initial_refinement 1 -dm_forest_maximum_refinement 3 -dm_p4est_refine_pattern hash 1778c4762a1bSJed Brown testset: 1779c4762a1bSJed Brown nsize: 4 1780dfd57a17SPierre Jolivet requires: hpddm slepc triangle !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES) 1781e600fa54SMatthew G. Knepley args: -run_type full -petscpartitioner_type simple -dm_refine 2 -bc_type dirichlet -petscspace_degree 2 -ksp_type gmres -ksp_gmres_restart 100 -pc_type hpddm -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0 -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_nev 4 -pc_hpddm_coarse_p 2 -pc_hpddm_coarse_pc_type redundant -ksp_rtol 1.e-1 1782c4762a1bSJed Brown test: 1783c4762a1bSJed Brown args: -pc_hpddm_coarse_mat_type baij -options_left no 1784c4762a1bSJed Brown suffix: tri_hpddm_reuse_baij 1785c4762a1bSJed Brown test: 1786c4762a1bSJed Brown requires: !complex 1787c4762a1bSJed Brown suffix: tri_hpddm_reuse 1788c4762a1bSJed Brown testset: 1789c4762a1bSJed Brown nsize: 4 1790dfd57a17SPierre Jolivet requires: hpddm slepc !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES) 1791e600fa54SMatthew G. Knepley args: -run_type full -petscpartitioner_type simple -dm_plex_box_faces 7,5 -dm_refine 2 -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 2 -ksp_type gmres -ksp_gmres_restart 100 -pc_type hpddm -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0 -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_nev 4 -pc_hpddm_coarse_p 2 -pc_hpddm_coarse_pc_type redundant -ksp_rtol 1.e-1 1792c4762a1bSJed Brown test: 1793c4762a1bSJed Brown args: -pc_hpddm_coarse_mat_type baij -options_left no 1794c4762a1bSJed Brown suffix: quad_hpddm_reuse_baij 1795c4762a1bSJed Brown test: 1796c4762a1bSJed Brown requires: !complex 1797c4762a1bSJed Brown suffix: quad_hpddm_reuse 1798c4762a1bSJed Brown testset: 1799c4762a1bSJed Brown nsize: 4 1800dfd57a17SPierre Jolivet requires: hpddm slepc !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES) 1801e600fa54SMatthew G. Knepley args: -run_type full -petscpartitioner_type simple -dm_plex_box_faces 7,5 -dm_refine 2 -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -pc_type hpddm -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0 -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_threshold 0.1 -pc_hpddm_coarse_p 2 -pc_hpddm_coarse_pc_type redundant -ksp_rtol 1.e-1 1802c4762a1bSJed Brown test: 1803c4762a1bSJed Brown args: -pc_hpddm_coarse_mat_type baij -options_left no 1804c4762a1bSJed Brown suffix: quad_hpddm_reuse_threshold_baij 1805c4762a1bSJed Brown test: 1806c4762a1bSJed Brown requires: !complex 1807c4762a1bSJed Brown suffix: quad_hpddm_reuse_threshold 1808c4762a1bSJed Brown testset: 1809c4762a1bSJed Brown nsize: 4 1810dfd57a17SPierre Jolivet requires: hpddm slepc parmetis !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES) 1811117ef88eSStefano Zampini filter: sed -e "s/linear solver iterations=17/linear solver iterations=16/g" 1812e600fa54SMatthew G. Knepley args: -run_type full -petscpartitioner_type parmetis -dm_refine 3 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -pc_type hpddm -snes_monitor_short -snes_converged_reason ::ascii_info_detail -snes_view -show_solution 0 -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type icc -pc_hpddm_levels_1_eps_nev 20 -pc_hpddm_coarse_p 2 -pc_hpddm_coarse_pc_type redundant -ksp_rtol 1.e-10 -dm_plex_filename ${PETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_boundary_label marker -pc_hpddm_levels_1_sub_pc_factor_levels 3 -variable_coefficient ball -dm_plex_gmsh_periodic 0 1813c4762a1bSJed Brown test: 1814c4762a1bSJed Brown args: -pc_hpddm_coarse_mat_type baij -options_left no 18156ba0327bSPierre Jolivet filter: grep -v " total: nonzeros=" | grep -v " rows=" | sed -e "s/total number of linear solver iterations=1[5-7]/total number of linear solver iterations=16/g" 1816c4762a1bSJed Brown suffix: tri_parmetis_hpddm_baij 1817c4762a1bSJed Brown test: 18186ba0327bSPierre Jolivet filter: grep -v " total: nonzeros=" | grep -v " rows=" | sed -e "s/total number of linear solver iterations=1[5-7]/total number of linear solver iterations=16/g" 1819c4762a1bSJed Brown requires: !complex 1820c4762a1bSJed Brown suffix: tri_parmetis_hpddm 1821d6837840SMatthew G. Knepley 1822d6837840SMatthew G. Knepley # 2D serial P1 tests for adaptive MG 1823d6837840SMatthew G. Knepley test: 1824d6837840SMatthew G. Knepley suffix: 2d_p1_adaptmg_0 1825d6837840SMatthew G. Knepley requires: triangle bamg 182630602db0SMatthew G. Knepley args: -dm_refine_hierarchy 3 -dm_plex_box_faces 4,4 -bc_type dirichlet -petscspace_degree 1 \ 1827d6837840SMatthew G. Knepley -variable_coefficient checkerboard_0 -mat_petscspace_degree 0 -div 16 -k 3 \ 1828d6837840SMatthew G. Knepley -snes_max_it 1 -ksp_converged_reason \ 1829d6837840SMatthew G. Knepley -ksp_rtol 1e-8 -pc_type mg 1830d6837840SMatthew G. Knepley # -ksp_monitor_true_residual -ksp_converged_reason -mg_levels_ksp_monitor_true_residual -pc_mg_mesp_monitor -dm_adapt_interp_view_fine draw -dm_adapt_interp_view_coarse draw -draw_pause 1 1831d6837840SMatthew G. Knepley test: 1832d6837840SMatthew G. Knepley suffix: 2d_p1_adaptmg_1 1833d6837840SMatthew G. Knepley requires: triangle bamg 183430602db0SMatthew G. Knepley args: -dm_refine_hierarchy 3 -dm_plex_box_faces 4,4 -bc_type dirichlet -petscspace_degree 1 \ 1835d6837840SMatthew G. Knepley -variable_coefficient checkerboard_0 -mat_petscspace_degree 0 -div 16 -k 3 \ 1836d6837840SMatthew G. Knepley -snes_max_it 1 -ksp_converged_reason \ 1837d6837840SMatthew G. Knepley -ksp_rtol 1e-8 -pc_type mg -pc_mg_galerkin -pc_mg_adapt_interp -pc_mg_adapt_interp_coarse_space eigenvector -pc_mg_adapt_interp_n 1 \ 1838d6837840SMatthew G. Knepley -pc_mg_mesp_ksp_type richardson -pc_mg_mesp_ksp_richardson_self_scale -pc_mg_mesp_ksp_max_it 100 -pc_mg_mesp_pc_type none 1839d6837840SMatthew G. Knepley 1840c4762a1bSJed Brown TEST*/ 1841