1 static char help[] = "Poisson Problem in 2d and 3d with simplicial finite elements.\n\ 2 We solve the Poisson problem in a rectangular\n\ 3 domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\ 4 This example supports discretized auxiliary fields (conductivity) as well as\n\ 5 multilevel nonlinear solvers.\n\n\n"; 6 7 /* 8 A visualization of the adaptation can be accomplished using: 9 10 -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 11 12 Information on refinement: 13 14 -info :~sys,vec,is,mat,ksp,snes,ts 15 */ 16 17 #include <petscdmplex.h> 18 #include <petscdmadaptor.h> 19 #include <petscsnes.h> 20 #include <petscds.h> 21 #include <petscviewerhdf5.h> 22 23 typedef enum {NEUMANN, DIRICHLET, NONE} BCType; 24 typedef enum {RUN_FULL, RUN_EXACT, RUN_TEST, RUN_PERF} RunType; 25 typedef enum {COEFF_NONE, COEFF_ANALYTIC, COEFF_FIELD, COEFF_NONLINEAR, COEFF_BALL, COEFF_CROSS, COEFF_CHECKERBOARD_0, COEFF_CHECKERBOARD_1} CoeffType; 26 27 typedef struct { 28 RunType runType; /* Whether to run tests, or solve the full problem */ 29 PetscBool jacobianMF; /* Whether to calculate the Jacobian action on the fly */ 30 PetscBool showInitial, showSolution, restart, quiet, nonzInit; 31 /* Problem definition */ 32 BCType bcType; 33 CoeffType variableCoefficient; 34 PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx); 35 PetscBool fieldBC; 36 void (**exactFields)(PetscInt, PetscInt, PetscInt, 37 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 38 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 39 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]); 40 PetscBool bdIntegral; /* Compute the integral of the solution on the boundary */ 41 /* Reproducing tests from SISC 40(3), pp. A1473-A1493, 2018 */ 42 PetscInt div; /* Number of divisions */ 43 PetscInt k; /* Parameter for checkerboard coefficient */ 44 PetscInt *kgrid; /* Random parameter grid */ 45 PetscBool rand; /* Make random assignments */ 46 /* Solver */ 47 PC pcmg; /* This is needed for error monitoring */ 48 PetscBool checkksp; /* Whether to check the KSPSolve for runType == RUN_TEST */ 49 } AppCtx; 50 51 static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 52 { 53 u[0] = 0.0; 54 return 0; 55 } 56 57 static PetscErrorCode ecks(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 58 { 59 u[0] = x[0]; 60 return 0; 61 } 62 63 /* 64 In 2D for Dirichlet conditions, we use exact solution: 65 66 u = x^2 + y^2 67 f = 4 68 69 so that 70 71 -\Delta u + f = -4 + 4 = 0 72 73 For Neumann conditions, we have 74 75 -\nabla u \cdot -\hat y |_{y=0} = (2y)|_{y=0} = 0 (bottom) 76 -\nabla u \cdot \hat y |_{y=1} = -(2y)|_{y=1} = -2 (top) 77 -\nabla u \cdot -\hat x |_{x=0} = (2x)|_{x=0} = 0 (left) 78 -\nabla u \cdot \hat x |_{x=1} = -(2x)|_{x=1} = -2 (right) 79 80 Which we can express as 81 82 \nabla u \cdot \hat n|_\Gamma = {2 x, 2 y} \cdot \hat n = 2 (x + y) 83 84 The boundary integral of this solution is (assuming we are not orienting the edges) 85 86 \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 87 */ 88 static PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 89 { 90 *u = x[0]*x[0] + x[1]*x[1]; 91 return 0; 92 } 93 94 static void quadratic_u_field_2d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 95 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 96 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 97 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar uexact[]) 98 { 99 uexact[0] = a[0]; 100 } 101 102 static PetscErrorCode ball_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 103 { 104 const PetscReal alpha = 500.; 105 const PetscReal radius2 = PetscSqr(0.15); 106 const PetscReal r2 = PetscSqr(x[0] - 0.5) + PetscSqr(x[1] - 0.5); 107 const PetscReal xi = alpha*(radius2 - r2); 108 109 *u = PetscTanhScalar(xi) + 1.0; 110 return 0; 111 } 112 113 static PetscErrorCode cross_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 114 { 115 const PetscReal alpha = 50*4; 116 const PetscReal xy = (x[0]-0.5)*(x[1]-0.5); 117 118 *u = PetscSinReal(alpha*xy) * (alpha*PetscAbsReal(xy) < 2*PETSC_PI ? 1 : 0.01); 119 return 0; 120 } 121 122 static void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 123 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 124 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 125 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 126 { 127 f0[0] = 4.0; 128 } 129 130 static void f0_ball_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 131 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 132 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 133 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 134 { 135 PetscInt d; 136 const PetscReal alpha = 500., radius2 = PetscSqr(0.15); 137 PetscReal r2, xi; 138 139 for (d = 0, r2 = 0.0; d < dim; ++d) r2 += PetscSqr(x[d] - 0.5); 140 xi = alpha*(radius2 - r2); 141 f0[0] = (-2.0*dim*alpha - 8.0*PetscSqr(alpha)*r2*PetscTanhReal(xi)) * PetscSqr(1.0/PetscCoshReal(xi)); 142 } 143 144 static void f0_cross_u_2d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 145 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 146 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 147 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 148 { 149 const PetscReal alpha = 50*4; 150 const PetscReal xy = (x[0]-0.5)*(x[1]-0.5); 151 152 f0[0] = PetscSinReal(alpha*xy) * (alpha*PetscAbsReal(xy) < 2*PETSC_PI ? 1 : 0.01); 153 } 154 155 static void f0_checkerboard_0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 156 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 157 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 158 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 159 { 160 f0[0] = -20.0*PetscExpReal(-(PetscSqr(x[0] - 0.5) + PetscSqr(x[1] - 0.5))); 161 } 162 163 static void f0_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 164 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 165 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 166 PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 167 { 168 PetscInt d; 169 for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += -n[d]*2.0*x[d]; 170 } 171 172 /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */ 173 static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 174 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 175 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 176 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 177 { 178 PetscInt d; 179 for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 180 } 181 182 /* < \nabla v, \nabla u + {\nabla u}^T > 183 This just gives \nabla u, give the perdiagonal for the transpose */ 184 static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 185 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 186 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 187 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 188 { 189 PetscInt d; 190 for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 191 } 192 193 /* 194 In 2D for x periodicity and y Dirichlet conditions, we use exact solution: 195 196 u = sin(2 pi x) 197 f = -4 pi^2 sin(2 pi x) 198 199 so that 200 201 -\Delta u + f = 4 pi^2 sin(2 pi x) - 4 pi^2 sin(2 pi x) = 0 202 */ 203 static PetscErrorCode xtrig_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 204 { 205 *u = PetscSinReal(2.0*PETSC_PI*x[0]); 206 return 0; 207 } 208 209 static void f0_xtrig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 210 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 211 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 212 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 213 { 214 f0[0] = -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[0]); 215 } 216 217 /* 218 In 2D for x-y periodicity, we use exact solution: 219 220 u = sin(2 pi x) sin(2 pi y) 221 f = -8 pi^2 sin(2 pi x) 222 223 so that 224 225 -\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 226 */ 227 static PetscErrorCode xytrig_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 228 { 229 *u = PetscSinReal(2.0*PETSC_PI*x[0])*PetscSinReal(2.0*PETSC_PI*x[1]); 230 return 0; 231 } 232 233 static void f0_xytrig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 234 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 235 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 236 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 237 { 238 f0[0] = -8.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[0]); 239 } 240 241 /* 242 In 2D for Dirichlet conditions with a variable coefficient, we use exact solution: 243 244 u = x^2 + y^2 245 f = 6 (x + y) 246 nu = (x + y) 247 248 so that 249 250 -\div \nu \grad u + f = -6 (x + y) + 6 (x + y) = 0 251 */ 252 static PetscErrorCode nu_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 253 { 254 *u = x[0] + x[1]; 255 return 0; 256 } 257 258 static PetscErrorCode checkerboardCoeff(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 259 { 260 AppCtx *user = (AppCtx *) ctx; 261 PetscInt div = user->div; 262 PetscInt k = user->k; 263 PetscInt mask = 0, ind = 0, d; 264 265 PetscFunctionBeginUser; 266 for (d = 0; d < dim; ++d) mask = (mask + (PetscInt) (x[d]*div)) % 2; 267 if (user->kgrid) { 268 for (d = 0; d < dim; ++d) { 269 if (d > 0) ind *= dim; 270 ind += (PetscInt) (x[d]*div); 271 } 272 k = user->kgrid[ind]; 273 } 274 u[0] = mask ? 1.0 : PetscPowRealInt(10.0, -k); 275 PetscFunctionReturn(0); 276 } 277 278 void f0_analytic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 279 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 280 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 281 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 282 { 283 f0[0] = 6.0*(x[0] + x[1]); 284 } 285 286 /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */ 287 void f1_analytic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 288 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 289 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 290 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 291 { 292 PetscInt d; 293 for (d = 0; d < dim; ++d) f1[d] = (x[0] + x[1])*u_x[d]; 294 } 295 296 void f1_field_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 297 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 298 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 299 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 300 { 301 PetscInt d; 302 for (d = 0; d < dim; ++d) f1[d] = a[0]*u_x[d]; 303 } 304 305 /* < \nabla v, \nabla u + {\nabla u}^T > 306 This just gives \nabla u, give the perdiagonal for the transpose */ 307 void g3_analytic_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 308 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 309 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 310 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 311 { 312 PetscInt d; 313 for (d = 0; d < dim; ++d) g3[d*dim+d] = x[0] + x[1]; 314 } 315 316 void g3_field_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 317 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 318 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 319 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 320 { 321 PetscInt d; 322 for (d = 0; d < dim; ++d) g3[d*dim+d] = a[0]; 323 } 324 325 /* 326 In 2D for Dirichlet conditions with a nonlinear coefficient (p-Laplacian with p = 4), we use exact solution: 327 328 u = x^2 + y^2 329 f = 16 (x^2 + y^2) 330 nu = 1/2 |grad u|^2 331 332 so that 333 334 -\div \nu \grad u + f = -16 (x^2 + y^2) + 16 (x^2 + y^2) = 0 335 */ 336 void f0_analytic_nonlinear_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 337 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 338 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 339 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 340 { 341 f0[0] = 16.0*(x[0]*x[0] + x[1]*x[1]); 342 } 343 344 /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */ 345 void f1_analytic_nonlinear_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 346 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 347 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 348 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 349 { 350 PetscScalar nu = 0.0; 351 PetscInt d; 352 for (d = 0; d < dim; ++d) nu += u_x[d]*u_x[d]; 353 for (d = 0; d < dim; ++d) f1[d] = 0.5*nu*u_x[d]; 354 } 355 356 /* 357 grad (u + eps w) - grad u = eps grad w 358 359 1/2 |grad (u + eps w)|^2 grad (u + eps w) - 1/2 |grad u|^2 grad u 360 = 1/2 (|grad u|^2 + 2 eps <grad u,grad w>) (grad u + eps grad w) - 1/2 |grad u|^2 grad u 361 = 1/2 (eps |grad u|^2 grad w + 2 eps <grad u,grad w> grad u) 362 = eps (1/2 |grad u|^2 grad w + grad u <grad u,grad w>) 363 */ 364 void g3_analytic_nonlinear_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 365 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 366 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 367 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 368 { 369 PetscScalar nu = 0.0; 370 PetscInt d, e; 371 for (d = 0; d < dim; ++d) nu += u_x[d]*u_x[d]; 372 for (d = 0; d < dim; ++d) { 373 g3[d*dim+d] = 0.5*nu; 374 for (e = 0; e < dim; ++e) { 375 g3[d*dim+e] += u_x[d]*u_x[e]; 376 } 377 } 378 } 379 380 /* 381 In 3D for Dirichlet conditions we use exact solution: 382 383 u = 2/3 (x^2 + y^2 + z^2) 384 f = 4 385 386 so that 387 388 -\Delta u + f = -2/3 * 6 + 4 = 0 389 390 For Neumann conditions, we have 391 392 -\nabla u \cdot -\hat z |_{z=0} = (2z)|_{z=0} = 0 (bottom) 393 -\nabla u \cdot \hat z |_{z=1} = -(2z)|_{z=1} = -2 (top) 394 -\nabla u \cdot -\hat y |_{y=0} = (2y)|_{y=0} = 0 (front) 395 -\nabla u \cdot \hat y |_{y=1} = -(2y)|_{y=1} = -2 (back) 396 -\nabla u \cdot -\hat x |_{x=0} = (2x)|_{x=0} = 0 (left) 397 -\nabla u \cdot \hat x |_{x=1} = -(2x)|_{x=1} = -2 (right) 398 399 Which we can express as 400 401 \nabla u \cdot \hat n|_\Gamma = {2 x, 2 y, 2z} \cdot \hat n = 2 (x + y + z) 402 */ 403 static PetscErrorCode quadratic_u_3d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 404 { 405 *u = 2.0*(x[0]*x[0] + x[1]*x[1] + x[2]*x[2])/3.0; 406 return 0; 407 } 408 409 static PetscErrorCode ball_u_3d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 410 { 411 const PetscReal alpha = 500.; 412 const PetscReal radius2 = PetscSqr(0.15); 413 const PetscReal r2 = PetscSqr(x[0] - 0.5) + PetscSqr(x[1] - 0.5) + PetscSqr(x[2] - 0.5); 414 const PetscReal xi = alpha*(radius2 - r2); 415 416 *u = PetscTanhScalar(xi) + 1.0; 417 return 0; 418 } 419 420 static void quadratic_u_field_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 421 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 422 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 423 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar uexact[]) 424 { 425 uexact[0] = a[0]; 426 } 427 428 static PetscErrorCode cross_u_3d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 429 { 430 const PetscReal alpha = 50*4; 431 const PetscReal xyz = (x[0]-0.5)*(x[1]-0.5)*(x[2]-0.5); 432 433 *u = PetscSinReal(alpha*xyz) * (alpha*PetscAbsReal(xyz) < 2*PETSC_PI ? (alpha*PetscAbsReal(xyz) > -2*PETSC_PI ? 1.0 : 0.01) : 0.01); 434 return 0; 435 } 436 437 static void f0_cross_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 438 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 439 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 440 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 441 { 442 const PetscReal alpha = 50*4; 443 const PetscReal xyz = (x[0]-0.5)*(x[1]-0.5)*(x[2]-0.5); 444 445 f0[0] = PetscSinReal(alpha*xyz) * (alpha*PetscAbsReal(xyz) < 2*PETSC_PI ? (alpha*PetscAbsReal(xyz) > -2*PETSC_PI ? 1.0 : 0.01) : 0.01); 446 } 447 448 static void bd_integral_2d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 449 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 450 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 451 PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *uint) 452 { 453 uint[0] = u[0]; 454 } 455 456 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 457 { 458 const char *bcTypes[3] = {"neumann", "dirichlet", "none"}; 459 const char *runTypes[4] = {"full", "exact", "test", "perf"}; 460 const char *coeffTypes[8] = {"none", "analytic", "field", "nonlinear", "ball", "cross", "checkerboard_0", "checkerboard_1"}; 461 PetscInt bc, run, coeff; 462 PetscErrorCode ierr; 463 464 PetscFunctionBeginUser; 465 options->runType = RUN_FULL; 466 options->bcType = DIRICHLET; 467 options->variableCoefficient = COEFF_NONE; 468 options->fieldBC = PETSC_FALSE; 469 options->jacobianMF = PETSC_FALSE; 470 options->showInitial = PETSC_FALSE; 471 options->showSolution = PETSC_FALSE; 472 options->restart = PETSC_FALSE; 473 options->quiet = PETSC_FALSE; 474 options->nonzInit = PETSC_FALSE; 475 options->bdIntegral = PETSC_FALSE; 476 options->checkksp = PETSC_FALSE; 477 options->div = 4; 478 options->k = 1; 479 options->kgrid = NULL; 480 options->rand = PETSC_FALSE; 481 482 ierr = PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");CHKERRQ(ierr); 483 run = options->runType; 484 ierr = PetscOptionsEList("-run_type", "The run type", "ex12.c", runTypes, 4, runTypes[options->runType], &run, NULL);CHKERRQ(ierr); 485 options->runType = (RunType) run; 486 bc = options->bcType; 487 ierr = PetscOptionsEList("-bc_type","Type of boundary condition","ex12.c",bcTypes,3,bcTypes[options->bcType],&bc,NULL);CHKERRQ(ierr); 488 options->bcType = (BCType) bc; 489 coeff = options->variableCoefficient; 490 ierr = PetscOptionsEList("-variable_coefficient","Type of variable coefficent","ex12.c",coeffTypes,8,coeffTypes[options->variableCoefficient],&coeff,NULL);CHKERRQ(ierr); 491 options->variableCoefficient = (CoeffType) coeff; 492 493 ierr = PetscOptionsBool("-field_bc", "Use a field representation for the BC", "ex12.c", options->fieldBC, &options->fieldBC, NULL);CHKERRQ(ierr); 494 ierr = PetscOptionsBool("-jacobian_mf", "Calculate the action of the Jacobian on the fly", "ex12.c", options->jacobianMF, &options->jacobianMF, NULL);CHKERRQ(ierr); 495 ierr = PetscOptionsBool("-show_initial", "Output the initial guess for verification", "ex12.c", options->showInitial, &options->showInitial, NULL);CHKERRQ(ierr); 496 ierr = PetscOptionsBool("-show_solution", "Output the solution for verification", "ex12.c", options->showSolution, &options->showSolution, NULL);CHKERRQ(ierr); 497 ierr = PetscOptionsBool("-restart", "Read in the mesh and solution from a file", "ex12.c", options->restart, &options->restart, NULL);CHKERRQ(ierr); 498 ierr = PetscOptionsBool("-quiet", "Don't print any vecs", "ex12.c", options->quiet, &options->quiet, NULL);CHKERRQ(ierr); 499 ierr = PetscOptionsBool("-nonzero_initial_guess", "nonzero initial guess", "ex12.c", options->nonzInit, &options->nonzInit, NULL);CHKERRQ(ierr); 500 ierr = PetscOptionsBool("-bd_integral", "Compute the integral of the solution on the boundary", "ex12.c", options->bdIntegral, &options->bdIntegral, NULL);CHKERRQ(ierr); 501 if (options->runType == RUN_TEST) { 502 ierr = PetscOptionsBool("-run_test_check_ksp", "Check solution of KSP", "ex12.c", options->checkksp, &options->checkksp, NULL);CHKERRQ(ierr); 503 } 504 ierr = PetscOptionsInt("-div", "The number of division for the checkerboard coefficient", "ex12.c", options->div, &options->div, NULL);CHKERRQ(ierr); 505 ierr = PetscOptionsInt("-k", "The exponent for the checkerboard coefficient", "ex12.c", options->k, &options->k, NULL);CHKERRQ(ierr); 506 ierr = PetscOptionsBool("-k_random", "Assign random k values to checkerboard", "ex12.c", options->rand, &options->rand, NULL);CHKERRQ(ierr); 507 ierr = PetscOptionsEnd();CHKERRQ(ierr); 508 PetscFunctionReturn(0); 509 } 510 511 static PetscErrorCode CreateBCLabel(DM dm, const char name[]) 512 { 513 DM plex; 514 DMLabel label; 515 PetscErrorCode ierr; 516 517 PetscFunctionBeginUser; 518 ierr = DMCreateLabel(dm, name);CHKERRQ(ierr); 519 ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr); 520 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 521 ierr = DMPlexMarkBoundaryFaces(plex, 1, label);CHKERRQ(ierr); 522 ierr = DMDestroy(&plex);CHKERRQ(ierr); 523 PetscFunctionReturn(0); 524 } 525 526 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 527 { 528 PetscErrorCode ierr; 529 530 PetscFunctionBeginUser; 531 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 532 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 533 ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 534 { 535 char convType[256]; 536 PetscBool flg; 537 538 ierr = PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");CHKERRQ(ierr); 539 ierr = PetscOptionsFList("-dm_plex_convert_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);CHKERRQ(ierr); 540 ierr = PetscOptionsEnd();CHKERRQ(ierr); 541 if (flg) { 542 DM dmConv; 543 544 ierr = DMConvert(*dm,convType,&dmConv);CHKERRQ(ierr); 545 if (dmConv) { 546 ierr = DMDestroy(dm);CHKERRQ(ierr); 547 *dm = dmConv; 548 } 549 ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 550 ierr = DMSetUp(*dm);CHKERRQ(ierr); 551 } 552 } 553 ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 554 if (user->rand) { 555 PetscRandom r; 556 PetscReal val; 557 PetscInt dim, N, i; 558 559 ierr = DMGetDimension(*dm, &dim);CHKERRQ(ierr); 560 N = PetscPowInt(user->div, dim); 561 ierr = PetscMalloc1(N, &user->kgrid);CHKERRQ(ierr); 562 ierr = PetscRandomCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr); 563 ierr = PetscRandomSetFromOptions(r);CHKERRQ(ierr); 564 ierr = PetscRandomSetInterval(r, 0.0, user->k);CHKERRQ(ierr); 565 ierr = PetscRandomSetSeed(r, 1973);CHKERRQ(ierr); 566 ierr = PetscRandomSeed(r);CHKERRQ(ierr); 567 for (i = 0; i < N; ++i) { 568 ierr = PetscRandomGetValueReal(r, &val);CHKERRQ(ierr); 569 user->kgrid[i] = 1 + (PetscInt) val; 570 } 571 ierr = PetscRandomDestroy(&r);CHKERRQ(ierr); 572 } 573 PetscFunctionReturn(0); 574 } 575 576 static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 577 { 578 PetscDS ds; 579 DMLabel label; 580 PetscWeakForm wf; 581 const DMBoundaryType *periodicity; 582 const PetscInt id = 1; 583 PetscInt bd, dim; 584 PetscErrorCode ierr; 585 586 PetscFunctionBeginUser; 587 ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 588 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 589 ierr = DMGetPeriodicity(dm, NULL, NULL, NULL, &periodicity);CHKERRQ(ierr); 590 switch (user->variableCoefficient) { 591 case COEFF_NONE: 592 if (periodicity && periodicity[0]) { 593 if (periodicity && periodicity[1]) { 594 ierr = PetscDSSetResidual(ds, 0, f0_xytrig_u, f1_u);CHKERRQ(ierr); 595 ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr); 596 } else { 597 ierr = PetscDSSetResidual(ds, 0, f0_xtrig_u, f1_u);CHKERRQ(ierr); 598 ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr); 599 } 600 } else { 601 ierr = PetscDSSetResidual(ds, 0, f0_u, f1_u);CHKERRQ(ierr); 602 ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr); 603 } 604 break; 605 case COEFF_ANALYTIC: 606 ierr = PetscDSSetResidual(ds, 0, f0_analytic_u, f1_analytic_u);CHKERRQ(ierr); 607 ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_analytic_uu);CHKERRQ(ierr); 608 break; 609 case COEFF_FIELD: 610 ierr = PetscDSSetResidual(ds, 0, f0_analytic_u, f1_field_u);CHKERRQ(ierr); 611 ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_field_uu);CHKERRQ(ierr); 612 break; 613 case COEFF_NONLINEAR: 614 ierr = PetscDSSetResidual(ds, 0, f0_analytic_nonlinear_u, f1_analytic_nonlinear_u);CHKERRQ(ierr); 615 ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_analytic_nonlinear_uu);CHKERRQ(ierr); 616 break; 617 case COEFF_BALL: 618 ierr = PetscDSSetResidual(ds, 0, f0_ball_u, f1_u);CHKERRQ(ierr); 619 ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr); 620 break; 621 case COEFF_CROSS: 622 switch (dim) { 623 case 2: 624 ierr = PetscDSSetResidual(ds, 0, f0_cross_u_2d, f1_u);CHKERRQ(ierr); 625 break; 626 case 3: 627 ierr = PetscDSSetResidual(ds, 0, f0_cross_u_3d, f1_u);CHKERRQ(ierr); 628 break; 629 default: 630 SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %d", dim); 631 } 632 ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr); 633 break; 634 case COEFF_CHECKERBOARD_0: 635 ierr = PetscDSSetResidual(ds, 0, f0_checkerboard_0_u, f1_field_u);CHKERRQ(ierr); 636 ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_field_uu);CHKERRQ(ierr); 637 break; 638 default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid variable coefficient type %d", user->variableCoefficient); 639 } 640 switch (dim) { 641 case 2: 642 switch (user->variableCoefficient) { 643 case COEFF_BALL: 644 user->exactFuncs[0] = ball_u_2d;break; 645 case COEFF_CROSS: 646 user->exactFuncs[0] = cross_u_2d;break; 647 case COEFF_CHECKERBOARD_0: 648 user->exactFuncs[0] = zero;break; 649 default: 650 if (periodicity && periodicity[0]) { 651 if (periodicity && periodicity[1]) { 652 user->exactFuncs[0] = xytrig_u_2d; 653 } else { 654 user->exactFuncs[0] = xtrig_u_2d; 655 } 656 } else { 657 user->exactFuncs[0] = quadratic_u_2d; 658 user->exactFields[0] = quadratic_u_field_2d; 659 } 660 } 661 if (user->bcType == NEUMANN) { 662 ierr = DMGetLabel(dm, "boundary", &label);CHKERRQ(ierr); 663 ierr = DMAddBoundary(dm, DM_BC_NATURAL, "wall", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd);CHKERRQ(ierr); 664 ierr = PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 665 ierr = PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_bd_u, 0, NULL);CHKERRQ(ierr); 666 } 667 break; 668 case 3: 669 switch (user->variableCoefficient) { 670 case COEFF_BALL: 671 user->exactFuncs[0] = ball_u_3d;break; 672 case COEFF_CROSS: 673 user->exactFuncs[0] = cross_u_3d;break; 674 default: 675 user->exactFuncs[0] = quadratic_u_3d; 676 user->exactFields[0] = quadratic_u_field_3d; 677 } 678 if (user->bcType == NEUMANN) { 679 ierr = DMGetLabel(dm, "boundary", &label);CHKERRQ(ierr); 680 ierr = DMAddBoundary(dm, DM_BC_NATURAL, "wall", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd);CHKERRQ(ierr); 681 ierr = PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 682 ierr = PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_bd_u, 0, NULL);CHKERRQ(ierr); 683 } 684 break; 685 default: 686 SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %d", dim); 687 } 688 /* Setup constants */ 689 switch (user->variableCoefficient) { 690 case COEFF_CHECKERBOARD_0: 691 { 692 PetscScalar constants[2]; 693 694 constants[0] = user->div; 695 constants[1] = user->k; 696 ierr = PetscDSSetConstants(ds, 2, constants);CHKERRQ(ierr); 697 } 698 break; 699 default: break; 700 } 701 ierr = PetscDSSetExactSolution(ds, 0, user->exactFuncs[0], user);CHKERRQ(ierr); 702 /* Setup Boundary Conditions */ 703 if (user->bcType == DIRICHLET) { 704 ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr); 705 if (!label) { 706 /* Right now, p4est cannot create labels immediately */ 707 ierr = 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);CHKERRQ(ierr); 708 } else { 709 ierr = 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);CHKERRQ(ierr); 710 } 711 } 712 PetscFunctionReturn(0); 713 } 714 715 static PetscErrorCode SetupMaterial(DM dm, DM dmAux, AppCtx *user) 716 { 717 PetscErrorCode (*matFuncs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx) = {nu_2d}; 718 void *ctx[1]; 719 Vec nu; 720 PetscErrorCode ierr; 721 722 PetscFunctionBegin; 723 ctx[0] = user; 724 if (user->variableCoefficient == COEFF_CHECKERBOARD_0) {matFuncs[0] = checkerboardCoeff;} 725 ierr = DMCreateLocalVector(dmAux, &nu);CHKERRQ(ierr); 726 ierr = PetscObjectSetName((PetscObject) nu, "Coefficient");CHKERRQ(ierr); 727 ierr = DMProjectFunctionLocal(dmAux, 0.0, matFuncs, ctx, INSERT_ALL_VALUES, nu);CHKERRQ(ierr); 728 ierr = DMSetAuxiliaryVec(dm, NULL, 0, nu);CHKERRQ(ierr); 729 ierr = VecDestroy(&nu);CHKERRQ(ierr); 730 PetscFunctionReturn(0); 731 } 732 733 static PetscErrorCode SetupBC(DM dm, DM dmAux, AppCtx *user) 734 { 735 PetscErrorCode (*bcFuncs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx); 736 Vec uexact; 737 PetscInt dim; 738 PetscErrorCode ierr; 739 740 PetscFunctionBegin; 741 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 742 if (dim == 2) bcFuncs[0] = quadratic_u_2d; 743 else bcFuncs[0] = quadratic_u_3d; 744 ierr = DMCreateLocalVector(dmAux, &uexact);CHKERRQ(ierr); 745 ierr = DMProjectFunctionLocal(dmAux, 0.0, bcFuncs, NULL, INSERT_ALL_VALUES, uexact);CHKERRQ(ierr); 746 ierr = DMSetAuxiliaryVec(dm, NULL, 0, uexact);CHKERRQ(ierr); 747 ierr = VecDestroy(&uexact);CHKERRQ(ierr); 748 PetscFunctionReturn(0); 749 } 750 751 static PetscErrorCode SetupAuxDM(DM dm, PetscFE feAux, AppCtx *user) 752 { 753 DM dmAux, coordDM; 754 PetscErrorCode ierr; 755 756 PetscFunctionBegin; 757 /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */ 758 ierr = DMGetCoordinateDM(dm, &coordDM);CHKERRQ(ierr); 759 if (!feAux) PetscFunctionReturn(0); 760 ierr = DMClone(dm, &dmAux);CHKERRQ(ierr); 761 ierr = DMSetCoordinateDM(dmAux, coordDM);CHKERRQ(ierr); 762 ierr = DMSetField(dmAux, 0, NULL, (PetscObject) feAux);CHKERRQ(ierr); 763 ierr = DMCreateDS(dmAux);CHKERRQ(ierr); 764 if (user->fieldBC) {ierr = SetupBC(dm, dmAux, user);CHKERRQ(ierr);} 765 else {ierr = SetupMaterial(dm, dmAux, user);CHKERRQ(ierr);} 766 ierr = DMDestroy(&dmAux);CHKERRQ(ierr); 767 PetscFunctionReturn(0); 768 } 769 770 static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 771 { 772 DM plex, cdm = dm; 773 PetscFE fe, feAux = NULL; 774 PetscBool simplex; 775 PetscInt dim; 776 MPI_Comm comm; 777 PetscErrorCode ierr; 778 779 PetscFunctionBeginUser; 780 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 781 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 782 ierr = DMPlexIsSimplex(plex, &simplex);CHKERRQ(ierr); 783 ierr = DMDestroy(&plex);CHKERRQ(ierr); 784 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 785 ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, -1, &fe);CHKERRQ(ierr); 786 ierr = PetscObjectSetName((PetscObject) fe, "potential");CHKERRQ(ierr); 787 if (user->variableCoefficient == COEFF_FIELD || user->variableCoefficient == COEFF_CHECKERBOARD_0) { 788 ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "mat_", -1, &feAux);CHKERRQ(ierr); 789 ierr = PetscObjectSetName((PetscObject) feAux, "coefficient");CHKERRQ(ierr); 790 ierr = PetscFECopyQuadrature(fe, feAux);CHKERRQ(ierr); 791 } else if (user->fieldBC) { 792 ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "bc_", -1, &feAux);CHKERRQ(ierr); 793 ierr = PetscFECopyQuadrature(fe, feAux);CHKERRQ(ierr); 794 } 795 /* Set discretization and boundary conditions for each mesh */ 796 ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 797 ierr = DMCreateDS(dm);CHKERRQ(ierr); 798 ierr = SetupProblem(dm, user);CHKERRQ(ierr); 799 while (cdm) { 800 ierr = SetupAuxDM(cdm, feAux, user);CHKERRQ(ierr); 801 if (user->bcType == DIRICHLET) { 802 PetscBool hasLabel; 803 804 ierr = DMHasLabel(cdm, "marker", &hasLabel);CHKERRQ(ierr); 805 if (!hasLabel) {ierr = CreateBCLabel(cdm, "marker");CHKERRQ(ierr);} 806 } 807 ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr); 808 ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 809 } 810 ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 811 ierr = PetscFEDestroy(&feAux);CHKERRQ(ierr); 812 PetscFunctionReturn(0); 813 } 814 815 int main(int argc, char **argv) 816 { 817 DM dm; /* Problem specification */ 818 SNES snes; /* nonlinear solver */ 819 Vec u; /* solution vector */ 820 Mat A,J; /* Jacobian matrix */ 821 MatNullSpace nullSpace; /* May be necessary for Neumann conditions */ 822 AppCtx user; /* user-defined work context */ 823 JacActionCtx userJ; /* context for Jacobian MF action */ 824 PetscReal error = 0.0; /* L_2 error in the solution */ 825 PetscErrorCode ierr; 826 827 ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 828 ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 829 ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr); 830 ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 831 ierr = SNESSetDM(snes, dm);CHKERRQ(ierr); 832 ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr); 833 834 ierr = PetscMalloc2(1, &user.exactFuncs, 1, &user.exactFields);CHKERRQ(ierr); 835 ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr); 836 837 ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 838 ierr = PetscObjectSetName((PetscObject) u, "potential");CHKERRQ(ierr); 839 840 ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr); 841 if (user.jacobianMF) { 842 PetscInt M, m, N, n; 843 844 ierr = MatGetSize(J, &M, &N);CHKERRQ(ierr); 845 ierr = MatGetLocalSize(J, &m, &n);CHKERRQ(ierr); 846 ierr = MatCreate(PETSC_COMM_WORLD, &A);CHKERRQ(ierr); 847 ierr = MatSetSizes(A, m, n, M, N);CHKERRQ(ierr); 848 ierr = MatSetType(A, MATSHELL);CHKERRQ(ierr); 849 ierr = MatSetUp(A);CHKERRQ(ierr); 850 #if 0 851 ierr = MatShellSetOperation(A, MATOP_MULT, (void (*)(void))FormJacobianAction);CHKERRQ(ierr); 852 #endif 853 854 userJ.dm = dm; 855 userJ.J = J; 856 userJ.user = &user; 857 858 ierr = DMCreateLocalVector(dm, &userJ.u);CHKERRQ(ierr); 859 if (user.fieldBC) {ierr = DMProjectFieldLocal(dm, 0.0, userJ.u, user.exactFields, INSERT_BC_VALUES, userJ.u);CHKERRQ(ierr);} 860 else {ierr = DMProjectFunctionLocal(dm, 0.0, user.exactFuncs, NULL, INSERT_BC_VALUES, userJ.u);CHKERRQ(ierr);} 861 ierr = MatShellSetContext(A, &userJ);CHKERRQ(ierr); 862 } else { 863 A = J; 864 } 865 866 nullSpace = NULL; 867 if (user.bcType != DIRICHLET) { 868 ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_TRUE, 0, NULL, &nullSpace);CHKERRQ(ierr); 869 ierr = MatSetNullSpace(A, nullSpace);CHKERRQ(ierr); 870 } 871 872 ierr = DMPlexSetSNESLocalFEM(dm,&user,&user,&user);CHKERRQ(ierr); 873 ierr = SNESSetJacobian(snes, A, J, NULL, NULL);CHKERRQ(ierr); 874 875 ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 876 877 if (user.fieldBC) {ierr = DMProjectField(dm, 0.0, u, user.exactFields, INSERT_ALL_VALUES, u);CHKERRQ(ierr);} 878 else {ierr = DMProjectFunction(dm, 0.0, user.exactFuncs, NULL, INSERT_ALL_VALUES, u);CHKERRQ(ierr);} 879 if (user.restart) { 880 #if defined(PETSC_HAVE_HDF5) 881 PetscViewer viewer; 882 char filename[PETSC_MAX_PATH_LEN]; 883 884 ierr = PetscOptionsGetString(NULL, NULL, "-dm_plex_filename", filename, sizeof(filename), NULL);CHKERRQ(ierr); 885 ierr = PetscViewerCreate(PETSC_COMM_WORLD, &viewer);CHKERRQ(ierr); 886 ierr = PetscViewerSetType(viewer, PETSCVIEWERHDF5);CHKERRQ(ierr); 887 ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr); 888 ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr); 889 ierr = PetscViewerHDF5PushGroup(viewer, "/fields");CHKERRQ(ierr); 890 ierr = VecLoad(u, viewer);CHKERRQ(ierr); 891 ierr = PetscViewerHDF5PopGroup(viewer);CHKERRQ(ierr); 892 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 893 #endif 894 } 895 if (user.showInitial) { 896 Vec lv; 897 ierr = DMGetLocalVector(dm, &lv);CHKERRQ(ierr); 898 ierr = DMGlobalToLocalBegin(dm, u, INSERT_VALUES, lv);CHKERRQ(ierr); 899 ierr = DMGlobalToLocalEnd(dm, u, INSERT_VALUES, lv);CHKERRQ(ierr); 900 ierr = DMPrintLocalVec(dm, "Local function", 1.0e-10, lv);CHKERRQ(ierr); 901 ierr = DMRestoreLocalVector(dm, &lv);CHKERRQ(ierr); 902 } 903 if (user.runType == RUN_FULL || user.runType == RUN_EXACT) { 904 PetscErrorCode (*initialGuess[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx) = {zero}; 905 906 if (user.nonzInit) initialGuess[0] = ecks; 907 if (user.runType == RUN_FULL) { 908 ierr = DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u);CHKERRQ(ierr); 909 } 910 ierr = VecViewFromOptions(u, NULL, "-guess_vec_view");CHKERRQ(ierr); 911 ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr); 912 ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr); 913 ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr); 914 915 if (user.showSolution) { 916 ierr = PetscPrintf(PETSC_COMM_WORLD, "Solution\n");CHKERRQ(ierr); 917 ierr = VecChop(u, 3.0e-9);CHKERRQ(ierr); 918 ierr = VecView(u, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 919 } 920 } else if (user.runType == RUN_PERF) { 921 Vec r; 922 PetscReal res = 0.0; 923 924 ierr = SNESGetFunction(snes, &r, NULL, NULL);CHKERRQ(ierr); 925 ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr); 926 ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n");CHKERRQ(ierr); 927 ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr); 928 ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr); 929 ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", (double)res);CHKERRQ(ierr); 930 } else { 931 Vec r; 932 PetscReal res = 0.0, tol = 1.0e-11; 933 934 /* Check discretization error */ 935 ierr = SNESGetFunction(snes, &r, NULL, NULL);CHKERRQ(ierr); 936 ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n");CHKERRQ(ierr); 937 if (!user.quiet) {ierr = VecView(u, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);} 938 ierr = DMComputeL2Diff(dm, 0.0, user.exactFuncs, NULL, u, &error);CHKERRQ(ierr); 939 if (error < tol) {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: < %2.1e\n", (double)tol);CHKERRQ(ierr);} 940 else {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g\n", (double)error);CHKERRQ(ierr);} 941 /* Check residual */ 942 ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr); 943 ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n");CHKERRQ(ierr); 944 ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr); 945 if (!user.quiet) {ierr = VecView(r, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);} 946 ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr); 947 ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", (double)res);CHKERRQ(ierr); 948 /* Check Jacobian */ 949 { 950 Vec b; 951 952 ierr = SNESComputeJacobian(snes, u, A, A);CHKERRQ(ierr); 953 ierr = VecDuplicate(u, &b);CHKERRQ(ierr); 954 ierr = VecSet(r, 0.0);CHKERRQ(ierr); 955 ierr = SNESComputeFunction(snes, r, b);CHKERRQ(ierr); 956 ierr = MatMult(A, u, r);CHKERRQ(ierr); 957 ierr = VecAXPY(r, 1.0, b);CHKERRQ(ierr); 958 ierr = PetscPrintf(PETSC_COMM_WORLD, "Au - b = Au + F(0)\n");CHKERRQ(ierr); 959 ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr); 960 if (!user.quiet) {ierr = VecView(r, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);} 961 ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr); 962 ierr = PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", (double)res);CHKERRQ(ierr); 963 /* check solver */ 964 if (user.checkksp) { 965 KSP ksp; 966 967 if (nullSpace) { 968 ierr = MatNullSpaceRemove(nullSpace, u);CHKERRQ(ierr); 969 } 970 ierr = SNESComputeJacobian(snes, u, A, J);CHKERRQ(ierr); 971 ierr = MatMult(A, u, b);CHKERRQ(ierr); 972 ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr); 973 ierr = KSPSetOperators(ksp, A, J);CHKERRQ(ierr); 974 ierr = KSPSolve(ksp, b, r);CHKERRQ(ierr); 975 ierr = VecAXPY(r, -1.0, u);CHKERRQ(ierr); 976 ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr); 977 ierr = PetscPrintf(PETSC_COMM_WORLD, "KSP Error: %g\n", (double)res);CHKERRQ(ierr); 978 } 979 ierr = VecDestroy(&b);CHKERRQ(ierr); 980 } 981 } 982 ierr = VecViewFromOptions(u, NULL, "-vec_view");CHKERRQ(ierr); 983 { 984 Vec nu; 985 986 ierr = DMGetAuxiliaryVec(dm, NULL, 0, &nu);CHKERRQ(ierr); 987 if (nu) {ierr = VecViewFromOptions(nu, NULL, "-coeff_view");CHKERRQ(ierr);} 988 } 989 990 if (user.bdIntegral) { 991 DMLabel label; 992 PetscInt id = 1; 993 PetscScalar bdInt = 0.0; 994 PetscReal exact = 3.3333333333; 995 996 ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr); 997 ierr = DMPlexComputeBdIntegral(dm, u, label, 1, &id, bd_integral_2d, &bdInt, NULL);CHKERRQ(ierr); 998 ierr = PetscPrintf(PETSC_COMM_WORLD, "Solution boundary integral: %.4g\n", (double) PetscAbsScalar(bdInt));CHKERRQ(ierr); 999 if (PetscAbsReal(PetscAbsScalar(bdInt) - exact) > PETSC_SQRT_MACHINE_EPSILON) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Invalid boundary integral %g != %g", (double) PetscAbsScalar(bdInt), (double)exact); 1000 } 1001 1002 ierr = MatNullSpaceDestroy(&nullSpace);CHKERRQ(ierr); 1003 if (user.jacobianMF) {ierr = VecDestroy(&userJ.u);CHKERRQ(ierr);} 1004 if (A != J) {ierr = MatDestroy(&A);CHKERRQ(ierr);} 1005 ierr = MatDestroy(&J);CHKERRQ(ierr); 1006 ierr = VecDestroy(&u);CHKERRQ(ierr); 1007 ierr = SNESDestroy(&snes);CHKERRQ(ierr); 1008 ierr = DMDestroy(&dm);CHKERRQ(ierr); 1009 ierr = PetscFree2(user.exactFuncs, user.exactFields);CHKERRQ(ierr); 1010 ierr = PetscFree(user.kgrid);CHKERRQ(ierr); 1011 ierr = PetscFinalize(); 1012 return ierr; 1013 } 1014 1015 /*TEST 1016 # 2D serial P1 test 0-4 1017 test: 1018 suffix: 2d_p1_0 1019 requires: triangle 1020 args: -run_type test -bc_type dirichlet -dm_plex_interpolate 0 -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1021 1022 test: 1023 suffix: 2d_p1_1 1024 requires: triangle 1025 args: -run_type test -bc_type dirichlet -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1026 1027 test: 1028 suffix: 2d_p1_2 1029 requires: triangle 1030 args: -run_type test -dm_refine_volume_limit_pre 0.0625 -bc_type dirichlet -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1031 1032 test: 1033 suffix: 2d_p1_neumann_0 1034 requires: triangle 1035 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 1036 1037 test: 1038 suffix: 2d_p1_neumann_1 1039 requires: triangle 1040 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 1041 1042 # 2D serial P2 test 5-8 1043 test: 1044 suffix: 2d_p2_0 1045 requires: triangle 1046 args: -run_type test -bc_type dirichlet -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1047 1048 test: 1049 suffix: 2d_p2_1 1050 requires: triangle 1051 args: -run_type test -dm_refine_volume_limit_pre 0.0625 -bc_type dirichlet -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1052 1053 test: 1054 suffix: 2d_p2_neumann_0 1055 requires: triangle 1056 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 1057 1058 test: 1059 suffix: 2d_p2_neumann_1 1060 requires: triangle 1061 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 1062 1063 test: 1064 suffix: bd_int_0 1065 requires: triangle 1066 args: -run_type test -bc_type dirichlet -petscspace_degree 2 -bd_integral -dm_view -quiet 1067 1068 test: 1069 suffix: bd_int_1 1070 requires: triangle 1071 args: -run_type test -dm_refine 2 -bc_type dirichlet -petscspace_degree 2 -bd_integral -dm_view -quiet 1072 1073 # 3D serial P1 test 9-12 1074 test: 1075 suffix: 3d_p1_0 1076 requires: ctetgen 1077 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 1078 1079 test: 1080 suffix: 3d_p1_1 1081 requires: ctetgen 1082 args: -run_type test -dm_plex_dim 3 -bc_type dirichlet -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 -dm_view 1083 1084 test: 1085 suffix: 3d_p1_2 1086 requires: ctetgen 1087 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 1088 1089 test: 1090 suffix: 3d_p1_neumann_0 1091 requires: ctetgen 1092 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 1093 1094 # Analytic variable coefficient 13-20 1095 test: 1096 suffix: 13 1097 requires: triangle 1098 args: -run_type test -variable_coefficient analytic -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1099 test: 1100 suffix: 14 1101 requires: triangle 1102 args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient analytic -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1103 test: 1104 suffix: 15 1105 requires: triangle 1106 args: -run_type test -variable_coefficient analytic -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1107 test: 1108 suffix: 16 1109 requires: triangle 1110 args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient analytic -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1111 test: 1112 suffix: 17 1113 requires: ctetgen 1114 args: -run_type test -dm_plex_dim 3 -variable_coefficient analytic -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1115 1116 test: 1117 suffix: 18 1118 requires: ctetgen 1119 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 1120 1121 test: 1122 suffix: 19 1123 requires: ctetgen 1124 args: -run_type test -dm_plex_dim 3 -variable_coefficient analytic -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1125 1126 test: 1127 suffix: 20 1128 requires: ctetgen 1129 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 1130 1131 # P1 variable coefficient 21-28 1132 test: 1133 suffix: 21 1134 requires: triangle 1135 args: -run_type test -variable_coefficient field -petscspace_degree 1 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1136 1137 test: 1138 suffix: 22 1139 requires: triangle 1140 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 1141 1142 test: 1143 suffix: 23 1144 requires: triangle 1145 args: -run_type test -variable_coefficient field -petscspace_degree 2 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1146 1147 test: 1148 suffix: 24 1149 requires: triangle 1150 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 1151 1152 test: 1153 suffix: 25 1154 requires: ctetgen 1155 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 1156 1157 test: 1158 suffix: 26 1159 requires: ctetgen 1160 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 1161 1162 test: 1163 suffix: 27 1164 requires: ctetgen 1165 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 1166 1167 test: 1168 suffix: 28 1169 requires: ctetgen 1170 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 1171 1172 # P0 variable coefficient 29-36 1173 test: 1174 suffix: 29 1175 requires: triangle 1176 args: -run_type test -variable_coefficient field -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1177 1178 test: 1179 suffix: 30 1180 requires: triangle 1181 args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient field -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1182 1183 test: 1184 suffix: 31 1185 requires: triangle 1186 args: -run_type test -variable_coefficient field -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1187 1188 test: 1189 requires: triangle 1190 suffix: 32 1191 args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient field -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1192 1193 test: 1194 requires: ctetgen 1195 suffix: 33 1196 args: -run_type test -dm_plex_dim 3 -variable_coefficient field -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1197 1198 test: 1199 suffix: 34 1200 requires: ctetgen 1201 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 1202 1203 test: 1204 suffix: 35 1205 requires: ctetgen 1206 args: -run_type test -dm_plex_dim 3 -variable_coefficient field -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1207 1208 test: 1209 suffix: 36 1210 requires: ctetgen 1211 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 1212 1213 # Full solve 39-44 1214 test: 1215 suffix: 39 1216 requires: triangle !single 1217 args: -run_type full -dm_refine_volume_limit_pre 0.015625 -petscspace_degree 2 -pc_type gamg -ksp_rtol 1.0e-10 -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason ::ascii_info_detail 1218 test: 1219 suffix: 40 1220 requires: triangle !single 1221 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 1222 test: 1223 suffix: 41 1224 requires: triangle !single 1225 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 1226 test: 1227 suffix: 42 1228 requires: triangle !single 1229 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 1230 test: 1231 suffix: 43 1232 requires: triangle !single 1233 nsize: 2 1234 args: -run_type full -dm_distribute -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 1235 1236 test: 1237 suffix: 44 1238 requires: triangle !single 1239 nsize: 2 1240 args: -run_type full -dm_distribute -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 1241 1242 # These tests use a loose tolerance just to exercise the PtAP operations for MATIS and multiple PCBDDC setup calls inside PCMG 1243 testset: 1244 requires: triangle !single 1245 nsize: 3 1246 args: -run_type full -dm_distribute -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 1247 test: 1248 suffix: gmg_bddc 1249 filter: sed -e "s/CONVERGED_FNORM_RELATIVE iterations 3/CONVERGED_FNORM_RELATIVE iterations 4/g" 1250 args: -mg_levels_pc_type jacobi 1251 test: 1252 filter: sed -e "s/iterations [0-4]/iterations 4/g" 1253 suffix: gmg_bddc_lev 1254 args: -mg_levels_pc_type bddc 1255 1256 # Restarting 1257 testset: 1258 suffix: restart 1259 requires: hdf5 triangle !complex 1260 args: -run_type test -bc_type dirichlet -petscspace_degree 1 1261 test: 1262 args: -dm_view hdf5:sol.h5 -vec_view hdf5:sol.h5::append 1263 test: 1264 args: -dm_plex_filename sol.h5 -dm_plex_name box -restart 1265 1266 # Periodicity 1267 test: 1268 suffix: periodic_0 1269 requires: triangle 1270 args: -run_type full -bc_type dirichlet -petscspace_degree 1 -snes_converged_reason ::ascii_info_detail 1271 1272 test: 1273 requires: !complex 1274 suffix: periodic_1 1275 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 1276 1277 # 2D serial P1 test with field bc 1278 test: 1279 suffix: field_bc_2d_p1_0 1280 requires: triangle 1281 args: -run_type test -bc_type dirichlet -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1282 1283 test: 1284 suffix: field_bc_2d_p1_1 1285 requires: triangle 1286 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 1287 1288 test: 1289 suffix: field_bc_2d_p1_neumann_0 1290 requires: triangle 1291 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 1292 1293 test: 1294 suffix: field_bc_2d_p1_neumann_1 1295 requires: triangle 1296 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 1297 1298 # 3D serial P1 test with field bc 1299 test: 1300 suffix: field_bc_3d_p1_0 1301 requires: ctetgen 1302 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 1303 1304 test: 1305 suffix: field_bc_3d_p1_1 1306 requires: ctetgen 1307 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 1308 1309 test: 1310 suffix: field_bc_3d_p1_neumann_0 1311 requires: ctetgen 1312 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 1313 1314 test: 1315 suffix: field_bc_3d_p1_neumann_1 1316 requires: ctetgen 1317 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 1318 1319 # 2D serial P2 test with field bc 1320 test: 1321 suffix: field_bc_2d_p2_0 1322 requires: triangle 1323 args: -run_type test -bc_type dirichlet -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1324 1325 test: 1326 suffix: field_bc_2d_p2_1 1327 requires: triangle 1328 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 1329 1330 test: 1331 suffix: field_bc_2d_p2_neumann_0 1332 requires: triangle 1333 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 1334 1335 test: 1336 suffix: field_bc_2d_p2_neumann_1 1337 requires: triangle 1338 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 1339 1340 # 3D serial P2 test with field bc 1341 test: 1342 suffix: field_bc_3d_p2_0 1343 requires: ctetgen 1344 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 1345 1346 test: 1347 suffix: field_bc_3d_p2_1 1348 requires: ctetgen 1349 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 1350 1351 test: 1352 suffix: field_bc_3d_p2_neumann_0 1353 requires: ctetgen 1354 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 1355 1356 test: 1357 suffix: field_bc_3d_p2_neumann_1 1358 requires: ctetgen 1359 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 1360 1361 # Full solve simplex: Convergence 1362 test: 1363 suffix: 3d_p1_conv 1364 requires: ctetgen 1365 args: -run_type full -dm_plex_dim 3 -dm_refine 1 -bc_type dirichlet -petscspace_degree 1 \ 1366 -snes_convergence_estimate -convest_num_refine 1 -pc_type lu 1367 1368 # Full solve simplex: PCBDDC 1369 test: 1370 suffix: tri_bddc 1371 requires: triangle !single 1372 nsize: 5 1373 args: -run_type full -dm_distribute -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 1374 1375 # Full solve simplex: PCBDDC 1376 test: 1377 suffix: tri_parmetis_bddc 1378 requires: triangle !single parmetis 1379 nsize: 4 1380 args: -run_type full -dm_distribute -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 1381 1382 testset: 1383 args: -run_type full -dm_distribute -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 1384 nsize: 5 1385 output_file: output/ex12_quad_bddc.out 1386 filter: sed -e "s/aijcusparse/aij/g" -e "s/aijviennacl/aij/g" -e "s/factorization: cusparse/factorization: petsc/g" 1387 test: 1388 requires: !single 1389 suffix: quad_bddc 1390 test: 1391 requires: !single cuda 1392 suffix: quad_bddc_cuda 1393 args: -matis_localmat_type aijcusparse -pc_bddc_dirichlet_pc_factor_mat_solver_type cusparse -pc_bddc_neumann_pc_factor_mat_solver_type cusparse 1394 test: 1395 requires: !single viennacl 1396 suffix: quad_bddc_viennacl 1397 args: -matis_localmat_type aijviennacl 1398 1399 # Full solve simplex: ASM 1400 test: 1401 suffix: tri_q2q1_asm_lu 1402 requires: triangle !single 1403 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 1404 1405 test: 1406 suffix: tri_q2q1_msm_lu 1407 requires: triangle !single 1408 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 1409 1410 test: 1411 suffix: tri_q2q1_asm_sor 1412 requires: triangle !single 1413 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 1414 1415 test: 1416 suffix: tri_q2q1_msm_sor 1417 requires: triangle !single 1418 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 1419 1420 # Full solve simplex: FAS 1421 test: 1422 suffix: fas_newton_0 1423 requires: triangle !single 1424 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 1425 1426 test: 1427 suffix: fas_newton_1 1428 requires: triangle !single 1429 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 1430 filter: sed -e "s/total number of linear solver iterations=14/total number of linear solver iterations=15/g" 1431 1432 test: 1433 suffix: fas_ngs_0 1434 requires: triangle !single 1435 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 1436 1437 # These two tests are broken because DMPlexComputeInjectorFEM() only works for regularly refined meshes 1438 test: 1439 suffix: fas_newton_coarse_0 1440 requires: pragmatic triangle 1441 TODO: broken 1442 args: -run_type full -variable_coefficient nonlinear -petscspace_degree 1 \ 1443 -dm_refine 2 -dm_coarsen_hierarchy 1 -dm_plex_hash_location -dm_adaptor pragmatic \ 1444 -snes_type fas -snes_fas_levels 2 -snes_converged_reason ::ascii_info_detail -snes_monitor_short -snes_view \ 1445 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -fas_coarse_snes_linesearch_type basic \ 1446 -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 1447 1448 test: 1449 suffix: mg_newton_coarse_0 1450 requires: triangle pragmatic 1451 TODO: broken 1452 args: -run_type full -petscspace_degree 1 \ 1453 -dm_refine 3 -dm_coarsen_hierarchy 3 -dm_plex_hash_location -dm_adaptor pragmatic \ 1454 -snes_atol 1.0e-8 -snes_rtol 0.0 -snes_monitor_short -snes_converged_reason ::ascii_info_detail -snes_view \ 1455 -ksp_type richardson -ksp_atol 1.0e-8 -ksp_rtol 0.0 -ksp_norm_type unpreconditioned -ksp_monitor_true_residual \ 1456 -pc_type mg -pc_mg_levels 4 \ 1457 -mg_levels_ksp_type gmres -mg_levels_pc_type ilu -mg_levels_ksp_max_it 10 1458 1459 # Full solve tensor 1460 test: 1461 suffix: tensor_plex_2d 1462 args: -run_type test -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -dm_refine_hierarchy 2 1463 1464 test: 1465 suffix: tensor_p4est_2d 1466 requires: p4est 1467 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 1468 1469 test: 1470 suffix: tensor_plex_3d 1471 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 1472 1473 test: 1474 suffix: tensor_p4est_3d 1475 requires: p4est 1476 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 1477 1478 test: 1479 suffix: p4est_test_q2_conformal_serial 1480 requires: p4est 1481 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 1482 1483 test: 1484 suffix: p4est_test_q2_conformal_parallel 1485 requires: p4est 1486 nsize: 7 1487 args: -run_type test -dm_distribute -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 1488 1489 test: 1490 suffix: p4est_test_q2_conformal_parallel_parmetis 1491 requires: parmetis p4est 1492 nsize: 4 1493 args: -run_type test -dm_distribute -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 1494 1495 test: 1496 suffix: p4est_test_q2_nonconformal_serial 1497 requires: p4est 1498 filter: grep -v "CG or CGNE: variant" 1499 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 1500 1501 test: 1502 suffix: p4est_test_q2_nonconformal_parallel 1503 requires: p4est 1504 filter: grep -v "CG or CGNE: variant" 1505 nsize: 7 1506 args: -run_type test -dm_distribute -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 1507 1508 test: 1509 suffix: p4est_test_q2_nonconformal_parallel_parmetis 1510 requires: parmetis p4est 1511 nsize: 4 1512 args: -run_type test -dm_distribute -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 1513 1514 test: 1515 suffix: p4est_exact_q2_conformal_serial 1516 requires: p4est !single !complex !__float128 1517 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 1518 1519 test: 1520 suffix: p4est_exact_q2_conformal_parallel 1521 requires: p4est !single !complex !__float128 1522 nsize: 4 1523 args: -run_type exact -dm_distribute -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 1524 1525 test: 1526 suffix: p4est_exact_q2_conformal_parallel_parmetis 1527 requires: parmetis p4est !single 1528 nsize: 4 1529 args: -run_type exact -dm_distribute -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 1530 1531 test: 1532 suffix: p4est_exact_q2_nonconformal_serial 1533 requires: p4est 1534 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 1535 1536 test: 1537 suffix: p4est_exact_q2_nonconformal_parallel 1538 requires: p4est 1539 nsize: 7 1540 args: -run_type exact -dm_distribute -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 1541 1542 test: 1543 suffix: p4est_exact_q2_nonconformal_parallel_parmetis 1544 requires: parmetis p4est 1545 nsize: 4 1546 args: -run_type exact -dm_distribute -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 1547 1548 test: 1549 suffix: p4est_full_q2_nonconformal_serial 1550 requires: p4est !single 1551 filter: grep -v "variant HERMITIAN" 1552 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 1553 1554 test: 1555 suffix: p4est_full_q2_nonconformal_parallel 1556 requires: p4est !single 1557 filter: grep -v "variant HERMITIAN" 1558 nsize: 7 1559 args: -run_type full -dm_distribute -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 1560 1561 test: 1562 suffix: p4est_full_q2_nonconformal_parallel_bddcfas 1563 requires: p4est !single 1564 filter: grep -v "variant HERMITIAN" 1565 nsize: 7 1566 args: -run_type full -dm_distribute -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 1567 1568 test: 1569 suffix: p4est_full_q2_nonconformal_parallel_bddc 1570 requires: p4est !single 1571 filter: grep -v "variant HERMITIAN" 1572 nsize: 7 1573 args: -run_type full -dm_distribute -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 1574 1575 test: 1576 TODO: broken 1577 suffix: p4est_fas_q2_conformal_serial 1578 requires: p4est !complex !__float128 1579 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 1580 1581 test: 1582 TODO: broken 1583 suffix: p4est_fas_q2_nonconformal_serial 1584 requires: p4est 1585 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 1586 1587 test: 1588 suffix: fas_newton_0_p4est 1589 requires: p4est !single !__float128 1590 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 1591 1592 # Full solve simplicial AMR 1593 test: 1594 suffix: tri_p1_adapt_init_pragmatic 1595 requires: pragmatic 1596 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 1597 1598 test: 1599 suffix: tri_p2_adapt_init_pragmatic 1600 requires: pragmatic 1601 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 1602 1603 test: 1604 suffix: tri_p1_adapt_init_mmg 1605 requires: mmg 1606 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 1607 1608 test: 1609 suffix: tri_p2_adapt_init_mmg 1610 requires: mmg 1611 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 1612 1613 test: 1614 suffix: tri_p1_adapt_seq_pragmatic 1615 requires: pragmatic 1616 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 1617 1618 test: 1619 suffix: tri_p2_adapt_seq_pragmatic 1620 requires: pragmatic 1621 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 1622 1623 test: 1624 suffix: tri_p1_adapt_seq_mmg 1625 requires: mmg 1626 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 1627 1628 test: 1629 suffix: tri_p2_adapt_seq_mmg 1630 requires: mmg 1631 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 1632 1633 test: 1634 suffix: tri_p1_adapt_analytic_pragmatic 1635 requires: pragmatic 1636 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 1637 1638 test: 1639 suffix: tri_p2_adapt_analytic_pragmatic 1640 requires: pragmatic 1641 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 1642 1643 test: 1644 suffix: tri_p1_adapt_analytic_mmg 1645 requires: mmg 1646 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 1647 1648 test: 1649 suffix: tri_p2_adapt_analytic_mmg 1650 requires: mmg 1651 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 1652 1653 test: 1654 suffix: tri_p1_adapt_uniform_pragmatic 1655 requires: pragmatic tetgen 1656 nsize: 4 1657 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_distribute -dm_adaptor pragmatic 1658 timeoutfactor: 2 1659 1660 test: 1661 suffix: tri_p2_adapt_uniform_pragmatic 1662 requires: pragmatic tetgen 1663 nsize: 4 1664 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_distribute -dm_adaptor pragmatic 1665 timeoutfactor: 1 1666 1667 test: 1668 suffix: tri_p1_adapt_uniform_mmg 1669 requires: mmg tetgen 1670 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 1671 timeoutfactor: 2 1672 1673 test: 1674 suffix: tri_p2_adapt_uniform_mmg 1675 requires: mmg tetgen 1676 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 1677 timeoutfactor: 1 1678 1679 test: 1680 suffix: tri_p1_adapt_uniform_parmmg 1681 requires: parmmg tetgen 1682 nsize: 4 1683 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_distribute -dm_adaptor parmmg 1684 timeoutfactor: 2 1685 1686 test: 1687 suffix: tri_p2_adapt_uniform_parmmg 1688 requires: parmmg tetgen 1689 nsize: 4 1690 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_distribute -dm_adaptor parmmg 1691 timeoutfactor: 1 1692 1693 # Full solve tensor AMR 1694 test: 1695 suffix: quad_q1_adapt_0 1696 requires: p4est 1697 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 1698 filter: grep -v DM_ 1699 1700 test: 1701 suffix: amr_0 1702 nsize: 5 1703 args: -run_type test -dm_distribute -petscpartitioner_type simple -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -dm_refine 1 1704 1705 test: 1706 suffix: amr_1 1707 requires: p4est !complex 1708 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 1709 1710 test: 1711 suffix: p4est_solve_bddc 1712 requires: p4est !complex 1713 args: -run_type full -dm_distribute -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 1714 nsize: 4 1715 1716 test: 1717 suffix: p4est_solve_fas 1718 requires: p4est 1719 args: -run_type full -dm_distribute -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 1720 nsize: 4 1721 TODO: identical machine two runs produce slightly different solver trackers 1722 1723 test: 1724 suffix: p4est_convergence_test_1 1725 requires: p4est 1726 args: -quiet -run_type test -dm_distribute -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 1727 nsize: 4 1728 1729 test: 1730 suffix: p4est_convergence_test_2 1731 requires: p4est 1732 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 1733 1734 test: 1735 suffix: p4est_convergence_test_3 1736 requires: p4est 1737 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 1738 1739 test: 1740 suffix: p4est_convergence_test_4 1741 requires: p4est 1742 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 1743 timeoutfactor: 5 1744 1745 # Serial tests with GLVis visualization 1746 test: 1747 suffix: glvis_2d_tet_p1 1748 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 1749 test: 1750 suffix: glvis_2d_tet_p2 1751 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 1752 test: 1753 suffix: glvis_2d_hex_p1 1754 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 1755 test: 1756 suffix: glvis_2d_hex_p2 1757 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 1758 test: 1759 suffix: glvis_2d_hex_p2_p4est 1760 requires: p4est 1761 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 1762 test: 1763 suffix: glvis_2d_tet_p0 1764 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 1765 test: 1766 suffix: glvis_2d_hex_p0 1767 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 1768 1769 # PCHPDDM tests 1770 testset: 1771 nsize: 4 1772 requires: hpddm slepc !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES) 1773 args: -run_type test -run_test_check_ksp -quiet -petscspace_degree 1 -dm_distribute -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 1774 test: 1775 suffix: quad_singular_hpddm 1776 args: -dm_plex_box_faces 6,7 1777 test: 1778 requires: p4est 1779 suffix: p4est_singular_2d_hpddm 1780 args: -dm_plex_convert_type p4est -dm_forest_minimum_refinement 1 -dm_forest_initial_refinement 3 -dm_forest_maximum_refinement 3 1781 test: 1782 requires: p4est 1783 suffix: p4est_nc_singular_2d_hpddm 1784 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 1785 testset: 1786 nsize: 4 1787 requires: hpddm slepc triangle !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES) 1788 args: -run_type full -dm_distribute -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 1789 test: 1790 args: -pc_hpddm_coarse_mat_type baij -options_left no 1791 suffix: tri_hpddm_reuse_baij 1792 test: 1793 requires: !complex 1794 suffix: tri_hpddm_reuse 1795 testset: 1796 nsize: 4 1797 requires: hpddm slepc !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES) 1798 args: -run_type full -dm_distribute -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 1799 test: 1800 args: -pc_hpddm_coarse_mat_type baij -options_left no 1801 suffix: quad_hpddm_reuse_baij 1802 test: 1803 requires: !complex 1804 suffix: quad_hpddm_reuse 1805 testset: 1806 nsize: 4 1807 requires: hpddm slepc !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES) 1808 args: -run_type full -dm_distribute -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 1809 test: 1810 args: -pc_hpddm_coarse_mat_type baij -options_left no 1811 suffix: quad_hpddm_reuse_threshold_baij 1812 test: 1813 requires: !complex 1814 suffix: quad_hpddm_reuse_threshold 1815 testset: 1816 nsize: 4 1817 requires: hpddm slepc parmetis !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES) 1818 filter: sed -e "s/linear solver iterations=17/linear solver iterations=16/g" 1819 args: -run_type full -dm_distribute -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 1820 test: 1821 args: -pc_hpddm_coarse_mat_type baij -options_left no 1822 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" 1823 suffix: tri_parmetis_hpddm_baij 1824 test: 1825 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" 1826 requires: !complex 1827 suffix: tri_parmetis_hpddm 1828 1829 # 2D serial P1 tests for adaptive MG 1830 test: 1831 suffix: 2d_p1_adaptmg_0 1832 requires: triangle bamg 1833 args: -dm_refine_hierarchy 3 -dm_plex_box_faces 4,4 -bc_type dirichlet -petscspace_degree 1 \ 1834 -variable_coefficient checkerboard_0 -mat_petscspace_degree 0 -div 16 -k 3 \ 1835 -snes_max_it 1 -ksp_converged_reason \ 1836 -ksp_rtol 1e-8 -pc_type mg 1837 # -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 1838 test: 1839 suffix: 2d_p1_adaptmg_1 1840 requires: triangle bamg 1841 args: -dm_refine_hierarchy 3 -dm_plex_box_faces 4,4 -bc_type dirichlet -petscspace_degree 1 \ 1842 -variable_coefficient checkerboard_0 -mat_petscspace_degree 0 -div 16 -k 3 \ 1843 -snes_max_it 1 -ksp_converged_reason \ 1844 -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 \ 1845 -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 1846 1847 TEST*/ 1848