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 463 PetscFunctionBeginUser; 464 options->runType = RUN_FULL; 465 options->bcType = DIRICHLET; 466 options->variableCoefficient = COEFF_NONE; 467 options->fieldBC = PETSC_FALSE; 468 options->jacobianMF = PETSC_FALSE; 469 options->showInitial = PETSC_FALSE; 470 options->showSolution = PETSC_FALSE; 471 options->restart = PETSC_FALSE; 472 options->quiet = PETSC_FALSE; 473 options->nonzInit = PETSC_FALSE; 474 options->bdIntegral = PETSC_FALSE; 475 options->checkksp = PETSC_FALSE; 476 options->div = 4; 477 options->k = 1; 478 options->kgrid = NULL; 479 options->rand = PETSC_FALSE; 480 481 PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX"); 482 run = options->runType; 483 PetscCall(PetscOptionsEList("-run_type", "The run type", "ex12.c", runTypes, 4, runTypes[options->runType], &run, NULL)); 484 options->runType = (RunType) run; 485 bc = options->bcType; 486 PetscCall(PetscOptionsEList("-bc_type","Type of boundary condition","ex12.c",bcTypes,3,bcTypes[options->bcType],&bc,NULL)); 487 options->bcType = (BCType) bc; 488 coeff = options->variableCoefficient; 489 PetscCall(PetscOptionsEList("-variable_coefficient","Type of variable coefficent","ex12.c",coeffTypes,8,coeffTypes[options->variableCoefficient],&coeff,NULL)); 490 options->variableCoefficient = (CoeffType) coeff; 491 492 PetscCall(PetscOptionsBool("-field_bc", "Use a field representation for the BC", "ex12.c", options->fieldBC, &options->fieldBC, NULL)); 493 PetscCall(PetscOptionsBool("-jacobian_mf", "Calculate the action of the Jacobian on the fly", "ex12.c", options->jacobianMF, &options->jacobianMF, NULL)); 494 PetscCall(PetscOptionsBool("-show_initial", "Output the initial guess for verification", "ex12.c", options->showInitial, &options->showInitial, NULL)); 495 PetscCall(PetscOptionsBool("-show_solution", "Output the solution for verification", "ex12.c", options->showSolution, &options->showSolution, NULL)); 496 PetscCall(PetscOptionsBool("-restart", "Read in the mesh and solution from a file", "ex12.c", options->restart, &options->restart, NULL)); 497 PetscCall(PetscOptionsBool("-quiet", "Don't print any vecs", "ex12.c", options->quiet, &options->quiet, NULL)); 498 PetscCall(PetscOptionsBool("-nonzero_initial_guess", "nonzero initial guess", "ex12.c", options->nonzInit, &options->nonzInit, NULL)); 499 PetscCall(PetscOptionsBool("-bd_integral", "Compute the integral of the solution on the boundary", "ex12.c", options->bdIntegral, &options->bdIntegral, NULL)); 500 if (options->runType == RUN_TEST) { 501 PetscCall(PetscOptionsBool("-run_test_check_ksp", "Check solution of KSP", "ex12.c", options->checkksp, &options->checkksp, NULL)); 502 } 503 PetscCall(PetscOptionsInt("-div", "The number of division for the checkerboard coefficient", "ex12.c", options->div, &options->div, NULL)); 504 PetscCall(PetscOptionsInt("-k", "The exponent for the checkerboard coefficient", "ex12.c", options->k, &options->k, NULL)); 505 PetscCall(PetscOptionsBool("-k_random", "Assign random k values to checkerboard", "ex12.c", options->rand, &options->rand, NULL)); 506 PetscOptionsEnd(); 507 PetscFunctionReturn(0); 508 } 509 510 static PetscErrorCode CreateBCLabel(DM dm, const char name[]) 511 { 512 DM plex; 513 DMLabel label; 514 515 PetscFunctionBeginUser; 516 PetscCall(DMCreateLabel(dm, name)); 517 PetscCall(DMGetLabel(dm, name, &label)); 518 PetscCall(DMConvert(dm, DMPLEX, &plex)); 519 PetscCall(DMPlexMarkBoundaryFaces(plex, 1, label)); 520 PetscCall(DMDestroy(&plex)); 521 PetscFunctionReturn(0); 522 } 523 524 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 525 { 526 PetscFunctionBeginUser; 527 PetscCall(DMCreate(comm, dm)); 528 PetscCall(DMSetType(*dm, DMPLEX)); 529 PetscCall(DMSetFromOptions(*dm)); 530 { 531 char convType[256]; 532 PetscBool flg; 533 534 PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX"); 535 PetscCall(PetscOptionsFList("-dm_plex_convert_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg)); 536 PetscOptionsEnd(); 537 if (flg) { 538 DM dmConv; 539 540 PetscCall(DMConvert(*dm,convType,&dmConv)); 541 if (dmConv) { 542 PetscCall(DMDestroy(dm)); 543 *dm = dmConv; 544 } 545 PetscCall(DMSetFromOptions(*dm)); 546 PetscCall(DMSetUp(*dm)); 547 } 548 } 549 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 550 if (user->rand) { 551 PetscRandom r; 552 PetscReal val; 553 PetscInt dim, N, i; 554 555 PetscCall(DMGetDimension(*dm, &dim)); 556 N = PetscPowInt(user->div, dim); 557 PetscCall(PetscMalloc1(N, &user->kgrid)); 558 PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &r)); 559 PetscCall(PetscRandomSetFromOptions(r)); 560 PetscCall(PetscRandomSetInterval(r, 0.0, user->k)); 561 PetscCall(PetscRandomSetSeed(r, 1973)); 562 PetscCall(PetscRandomSeed(r)); 563 for (i = 0; i < N; ++i) { 564 PetscCall(PetscRandomGetValueReal(r, &val)); 565 user->kgrid[i] = 1 + (PetscInt) val; 566 } 567 PetscCall(PetscRandomDestroy(&r)); 568 } 569 PetscFunctionReturn(0); 570 } 571 572 static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 573 { 574 PetscDS ds; 575 DMLabel label; 576 PetscWeakForm wf; 577 const PetscReal *L; 578 const PetscInt id = 1; 579 PetscInt bd, dim; 580 581 PetscFunctionBeginUser; 582 PetscCall(DMGetDS(dm, &ds)); 583 PetscCall(DMGetDimension(dm, &dim)); 584 PetscCall(DMGetPeriodicity(dm, NULL, NULL, &L)); 585 switch (user->variableCoefficient) { 586 case COEFF_NONE: 587 if (L && L[0]) { 588 if (L && L[1]) { 589 PetscCall(PetscDSSetResidual(ds, 0, f0_xytrig_u, f1_u)); 590 PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 591 } else { 592 PetscCall(PetscDSSetResidual(ds, 0, f0_xtrig_u, f1_u)); 593 PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 594 } 595 } else { 596 PetscCall(PetscDSSetResidual(ds, 0, f0_u, f1_u)); 597 PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 598 } 599 break; 600 case COEFF_ANALYTIC: 601 PetscCall(PetscDSSetResidual(ds, 0, f0_analytic_u, f1_analytic_u)); 602 PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_analytic_uu)); 603 break; 604 case COEFF_FIELD: 605 PetscCall(PetscDSSetResidual(ds, 0, f0_analytic_u, f1_field_u)); 606 PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_field_uu)); 607 break; 608 case COEFF_NONLINEAR: 609 PetscCall(PetscDSSetResidual(ds, 0, f0_analytic_nonlinear_u, f1_analytic_nonlinear_u)); 610 PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_analytic_nonlinear_uu)); 611 break; 612 case COEFF_BALL: 613 PetscCall(PetscDSSetResidual(ds, 0, f0_ball_u, f1_u)); 614 PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 615 break; 616 case COEFF_CROSS: 617 switch (dim) { 618 case 2: 619 PetscCall(PetscDSSetResidual(ds, 0, f0_cross_u_2d, f1_u)); 620 break; 621 case 3: 622 PetscCall(PetscDSSetResidual(ds, 0, f0_cross_u_3d, f1_u)); 623 break; 624 default: 625 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %" PetscInt_FMT, dim); 626 } 627 PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 628 break; 629 case COEFF_CHECKERBOARD_0: 630 PetscCall(PetscDSSetResidual(ds, 0, f0_checkerboard_0_u, f1_field_u)); 631 PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_field_uu)); 632 break; 633 default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid variable coefficient type %d", user->variableCoefficient); 634 } 635 switch (dim) { 636 case 2: 637 switch (user->variableCoefficient) { 638 case COEFF_BALL: 639 user->exactFuncs[0] = ball_u_2d;break; 640 case COEFF_CROSS: 641 user->exactFuncs[0] = cross_u_2d;break; 642 case COEFF_CHECKERBOARD_0: 643 user->exactFuncs[0] = zero;break; 644 default: 645 if (L && L[0]) { 646 if (L && L[1]) { 647 user->exactFuncs[0] = xytrig_u_2d; 648 } else { 649 user->exactFuncs[0] = xtrig_u_2d; 650 } 651 } else { 652 user->exactFuncs[0] = quadratic_u_2d; 653 user->exactFields[0] = quadratic_u_field_2d; 654 } 655 } 656 if (user->bcType == NEUMANN) { 657 PetscCall(DMGetLabel(dm, "boundary", &label)); 658 PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "wall", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 659 PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 660 PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_bd_u, 0, NULL)); 661 } 662 break; 663 case 3: 664 switch (user->variableCoefficient) { 665 case COEFF_BALL: 666 user->exactFuncs[0] = ball_u_3d;break; 667 case COEFF_CROSS: 668 user->exactFuncs[0] = cross_u_3d;break; 669 default: 670 user->exactFuncs[0] = quadratic_u_3d; 671 user->exactFields[0] = quadratic_u_field_3d; 672 } 673 if (user->bcType == NEUMANN) { 674 PetscCall(DMGetLabel(dm, "boundary", &label)); 675 PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "wall", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 676 PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 677 PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_bd_u, 0, NULL)); 678 } 679 break; 680 default: 681 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %" PetscInt_FMT, dim); 682 } 683 /* Setup constants */ 684 switch (user->variableCoefficient) { 685 case COEFF_CHECKERBOARD_0: 686 { 687 PetscScalar constants[2]; 688 689 constants[0] = user->div; 690 constants[1] = user->k; 691 PetscCall(PetscDSSetConstants(ds, 2, constants)); 692 } 693 break; 694 default: break; 695 } 696 PetscCall(PetscDSSetExactSolution(ds, 0, user->exactFuncs[0], user)); 697 /* Setup Boundary Conditions */ 698 if (user->bcType == DIRICHLET) { 699 PetscCall(DMGetLabel(dm, "marker", &label)); 700 if (!label) { 701 /* Right now, p4est cannot create labels immediately */ 702 PetscCall(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)); 703 } else { 704 PetscCall(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)); 705 } 706 } 707 PetscFunctionReturn(0); 708 } 709 710 static PetscErrorCode SetupMaterial(DM dm, DM dmAux, AppCtx *user) 711 { 712 PetscErrorCode (*matFuncs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx) = {nu_2d}; 713 void *ctx[1]; 714 Vec nu; 715 716 PetscFunctionBegin; 717 ctx[0] = user; 718 if (user->variableCoefficient == COEFF_CHECKERBOARD_0) {matFuncs[0] = checkerboardCoeff;} 719 PetscCall(DMCreateLocalVector(dmAux, &nu)); 720 PetscCall(PetscObjectSetName((PetscObject) nu, "Coefficient")); 721 PetscCall(DMProjectFunctionLocal(dmAux, 0.0, matFuncs, ctx, INSERT_ALL_VALUES, nu)); 722 PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, nu)); 723 PetscCall(VecDestroy(&nu)); 724 PetscFunctionReturn(0); 725 } 726 727 static PetscErrorCode SetupBC(DM dm, DM dmAux, AppCtx *user) 728 { 729 PetscErrorCode (*bcFuncs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx); 730 Vec uexact; 731 PetscInt dim; 732 733 PetscFunctionBegin; 734 PetscCall(DMGetDimension(dm, &dim)); 735 if (dim == 2) bcFuncs[0] = quadratic_u_2d; 736 else bcFuncs[0] = quadratic_u_3d; 737 PetscCall(DMCreateLocalVector(dmAux, &uexact)); 738 PetscCall(DMProjectFunctionLocal(dmAux, 0.0, bcFuncs, NULL, INSERT_ALL_VALUES, uexact)); 739 PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, uexact)); 740 PetscCall(VecDestroy(&uexact)); 741 PetscFunctionReturn(0); 742 } 743 744 static PetscErrorCode SetupAuxDM(DM dm, PetscFE feAux, AppCtx *user) 745 { 746 DM dmAux, coordDM; 747 748 PetscFunctionBegin; 749 /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */ 750 PetscCall(DMGetCoordinateDM(dm, &coordDM)); 751 if (!feAux) PetscFunctionReturn(0); 752 PetscCall(DMClone(dm, &dmAux)); 753 PetscCall(DMSetCoordinateDM(dmAux, coordDM)); 754 PetscCall(DMSetField(dmAux, 0, NULL, (PetscObject) feAux)); 755 PetscCall(DMCreateDS(dmAux)); 756 if (user->fieldBC) PetscCall(SetupBC(dm, dmAux, user)); 757 else PetscCall(SetupMaterial(dm, dmAux, user)); 758 PetscCall(DMDestroy(&dmAux)); 759 PetscFunctionReturn(0); 760 } 761 762 static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 763 { 764 DM plex, cdm = dm; 765 PetscFE fe, feAux = NULL; 766 PetscBool simplex; 767 PetscInt dim; 768 MPI_Comm comm; 769 770 PetscFunctionBeginUser; 771 PetscCall(DMGetDimension(dm, &dim)); 772 PetscCall(DMConvert(dm, DMPLEX, &plex)); 773 PetscCall(DMPlexIsSimplex(plex, &simplex)); 774 PetscCall(DMDestroy(&plex)); 775 PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 776 PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, -1, &fe)); 777 PetscCall(PetscObjectSetName((PetscObject) fe, "potential")); 778 if (user->variableCoefficient == COEFF_FIELD || user->variableCoefficient == COEFF_CHECKERBOARD_0) { 779 PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "mat_", -1, &feAux)); 780 PetscCall(PetscObjectSetName((PetscObject) feAux, "coefficient")); 781 PetscCall(PetscFECopyQuadrature(fe, feAux)); 782 } else if (user->fieldBC) { 783 PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "bc_", -1, &feAux)); 784 PetscCall(PetscFECopyQuadrature(fe, feAux)); 785 } 786 /* Set discretization and boundary conditions for each mesh */ 787 PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe)); 788 PetscCall(DMCreateDS(dm)); 789 PetscCall(SetupProblem(dm, user)); 790 while (cdm) { 791 PetscCall(SetupAuxDM(cdm, feAux, user)); 792 if (user->bcType == DIRICHLET) { 793 PetscBool hasLabel; 794 795 PetscCall(DMHasLabel(cdm, "marker", &hasLabel)); 796 if (!hasLabel) PetscCall(CreateBCLabel(cdm, "marker")); 797 } 798 PetscCall(DMCopyDisc(dm, cdm)); 799 PetscCall(DMGetCoarseDM(cdm, &cdm)); 800 } 801 PetscCall(PetscFEDestroy(&fe)); 802 PetscCall(PetscFEDestroy(&feAux)); 803 PetscFunctionReturn(0); 804 } 805 806 int main(int argc, char **argv) 807 { 808 DM dm; /* Problem specification */ 809 SNES snes; /* nonlinear solver */ 810 Vec u; /* solution vector */ 811 Mat A,J; /* Jacobian matrix */ 812 MatNullSpace nullSpace; /* May be necessary for Neumann conditions */ 813 AppCtx user; /* user-defined work context */ 814 JacActionCtx userJ; /* context for Jacobian MF action */ 815 PetscReal error = 0.0; /* L_2 error in the solution */ 816 817 PetscFunctionBeginUser; 818 PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 819 PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 820 PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 821 PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 822 PetscCall(SNESSetDM(snes, dm)); 823 PetscCall(DMSetApplicationContext(dm, &user)); 824 825 PetscCall(PetscMalloc2(1, &user.exactFuncs, 1, &user.exactFields)); 826 PetscCall(SetupDiscretization(dm, &user)); 827 828 PetscCall(DMCreateGlobalVector(dm, &u)); 829 PetscCall(PetscObjectSetName((PetscObject) u, "potential")); 830 831 PetscCall(DMCreateMatrix(dm, &J)); 832 if (user.jacobianMF) { 833 PetscInt M, m, N, n; 834 835 PetscCall(MatGetSize(J, &M, &N)); 836 PetscCall(MatGetLocalSize(J, &m, &n)); 837 PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 838 PetscCall(MatSetSizes(A, m, n, M, N)); 839 PetscCall(MatSetType(A, MATSHELL)); 840 PetscCall(MatSetUp(A)); 841 #if 0 842 PetscCall(MatShellSetOperation(A, MATOP_MULT, (void (*)(void))FormJacobianAction)); 843 #endif 844 845 userJ.dm = dm; 846 userJ.J = J; 847 userJ.user = &user; 848 849 PetscCall(DMCreateLocalVector(dm, &userJ.u)); 850 if (user.fieldBC) PetscCall(DMProjectFieldLocal(dm, 0.0, userJ.u, user.exactFields, INSERT_BC_VALUES, userJ.u)); 851 else PetscCall(DMProjectFunctionLocal(dm, 0.0, user.exactFuncs, NULL, INSERT_BC_VALUES, userJ.u)); 852 PetscCall(MatShellSetContext(A, &userJ)); 853 } else { 854 A = J; 855 } 856 857 nullSpace = NULL; 858 if (user.bcType != DIRICHLET) { 859 PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_TRUE, 0, NULL, &nullSpace)); 860 PetscCall(MatSetNullSpace(A, nullSpace)); 861 } 862 863 PetscCall(DMPlexSetSNESLocalFEM(dm,&user,&user,&user)); 864 PetscCall(SNESSetJacobian(snes, A, J, NULL, NULL)); 865 866 PetscCall(SNESSetFromOptions(snes)); 867 868 if (user.fieldBC) PetscCall(DMProjectField(dm, 0.0, u, user.exactFields, INSERT_ALL_VALUES, u)); 869 else PetscCall(DMProjectFunction(dm, 0.0, user.exactFuncs, NULL, INSERT_ALL_VALUES, u)); 870 if (user.restart) { 871 #if defined(PETSC_HAVE_HDF5) 872 PetscViewer viewer; 873 char filename[PETSC_MAX_PATH_LEN]; 874 875 PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_filename", filename, sizeof(filename), NULL)); 876 PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer)); 877 PetscCall(PetscViewerSetType(viewer, PETSCVIEWERHDF5)); 878 PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_READ)); 879 PetscCall(PetscViewerFileSetName(viewer, filename)); 880 PetscCall(PetscViewerHDF5PushGroup(viewer, "/fields")); 881 PetscCall(VecLoad(u, viewer)); 882 PetscCall(PetscViewerHDF5PopGroup(viewer)); 883 PetscCall(PetscViewerDestroy(&viewer)); 884 #endif 885 } 886 if (user.showInitial) { 887 Vec lv; 888 PetscCall(DMGetLocalVector(dm, &lv)); 889 PetscCall(DMGlobalToLocalBegin(dm, u, INSERT_VALUES, lv)); 890 PetscCall(DMGlobalToLocalEnd(dm, u, INSERT_VALUES, lv)); 891 PetscCall(DMPrintLocalVec(dm, "Local function", 1.0e-10, lv)); 892 PetscCall(DMRestoreLocalVector(dm, &lv)); 893 } 894 if (user.runType == RUN_FULL || user.runType == RUN_EXACT) { 895 PetscErrorCode (*initialGuess[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx) = {zero}; 896 897 if (user.nonzInit) initialGuess[0] = ecks; 898 if (user.runType == RUN_FULL) { 899 PetscCall(DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u)); 900 } 901 PetscCall(VecViewFromOptions(u, NULL, "-guess_vec_view")); 902 PetscCall(SNESSolve(snes, NULL, u)); 903 PetscCall(SNESGetSolution(snes, &u)); 904 PetscCall(SNESGetDM(snes, &dm)); 905 906 if (user.showSolution) { 907 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Solution\n")); 908 PetscCall(VecChop(u, 3.0e-9)); 909 PetscCall(VecView(u, PETSC_VIEWER_STDOUT_WORLD)); 910 } 911 } else if (user.runType == RUN_PERF) { 912 Vec r; 913 PetscReal res = 0.0; 914 915 PetscCall(SNESGetFunction(snes, &r, NULL, NULL)); 916 PetscCall(SNESComputeFunction(snes, u, r)); 917 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n")); 918 PetscCall(VecChop(r, 1.0e-10)); 919 PetscCall(VecNorm(r, NORM_2, &res)); 920 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", (double)res)); 921 } else { 922 Vec r; 923 PetscReal res = 0.0, tol = 1.0e-11; 924 925 /* Check discretization error */ 926 PetscCall(SNESGetFunction(snes, &r, NULL, NULL)); 927 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n")); 928 if (!user.quiet) PetscCall(VecView(u, PETSC_VIEWER_STDOUT_WORLD)); 929 PetscCall(DMComputeL2Diff(dm, 0.0, user.exactFuncs, NULL, u, &error)); 930 if (error < tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: < %2.1e\n", (double)tol)); 931 else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g\n", (double)error)); 932 /* Check residual */ 933 PetscCall(SNESComputeFunction(snes, u, r)); 934 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n")); 935 PetscCall(VecChop(r, 1.0e-10)); 936 if (!user.quiet) PetscCall(VecView(r, PETSC_VIEWER_STDOUT_WORLD)); 937 PetscCall(VecNorm(r, NORM_2, &res)); 938 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", (double)res)); 939 /* Check Jacobian */ 940 { 941 Vec b; 942 943 PetscCall(SNESComputeJacobian(snes, u, A, A)); 944 PetscCall(VecDuplicate(u, &b)); 945 PetscCall(VecSet(r, 0.0)); 946 PetscCall(SNESComputeFunction(snes, r, b)); 947 PetscCall(MatMult(A, u, r)); 948 PetscCall(VecAXPY(r, 1.0, b)); 949 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Au - b = Au + F(0)\n")); 950 PetscCall(VecChop(r, 1.0e-10)); 951 if (!user.quiet) PetscCall(VecView(r, PETSC_VIEWER_STDOUT_WORLD)); 952 PetscCall(VecNorm(r, NORM_2, &res)); 953 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", (double)res)); 954 /* check solver */ 955 if (user.checkksp) { 956 KSP ksp; 957 958 if (nullSpace) PetscCall(MatNullSpaceRemove(nullSpace, u)); 959 PetscCall(SNESComputeJacobian(snes, u, A, J)); 960 PetscCall(MatMult(A, u, b)); 961 PetscCall(SNESGetKSP(snes, &ksp)); 962 PetscCall(KSPSetOperators(ksp, A, J)); 963 PetscCall(KSPSolve(ksp, b, r)); 964 PetscCall(VecAXPY(r, -1.0, u)); 965 PetscCall(VecNorm(r, NORM_2, &res)); 966 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP Error: %g\n", (double)res)); 967 } 968 PetscCall(VecDestroy(&b)); 969 } 970 } 971 PetscCall(VecViewFromOptions(u, NULL, "-vec_view")); 972 { 973 Vec nu; 974 975 PetscCall(DMGetAuxiliaryVec(dm, NULL, 0, 0, &nu)); 976 if (nu) PetscCall(VecViewFromOptions(nu, NULL, "-coeff_view")); 977 } 978 979 if (user.bdIntegral) { 980 DMLabel label; 981 PetscInt id = 1; 982 PetscScalar bdInt = 0.0; 983 PetscReal exact = 3.3333333333; 984 985 PetscCall(DMGetLabel(dm, "marker", &label)); 986 PetscCall(DMPlexComputeBdIntegral(dm, u, label, 1, &id, bd_integral_2d, &bdInt, NULL)); 987 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Solution boundary integral: %.4g\n", (double) PetscAbsScalar(bdInt))); 988 PetscCheck(PetscAbsReal(PetscAbsScalar(bdInt) - exact) <= PETSC_SQRT_MACHINE_EPSILON,PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Invalid boundary integral %g != %g", (double) PetscAbsScalar(bdInt), (double)exact); 989 } 990 991 PetscCall(MatNullSpaceDestroy(&nullSpace)); 992 if (user.jacobianMF) PetscCall(VecDestroy(&userJ.u)); 993 if (A != J) PetscCall(MatDestroy(&A)); 994 PetscCall(MatDestroy(&J)); 995 PetscCall(VecDestroy(&u)); 996 PetscCall(SNESDestroy(&snes)); 997 PetscCall(DMDestroy(&dm)); 998 PetscCall(PetscFree2(user.exactFuncs, user.exactFields)); 999 PetscCall(PetscFree(user.kgrid)); 1000 PetscCall(PetscFinalize()); 1001 return 0; 1002 } 1003 1004 /*TEST 1005 # 2D serial P1 test 0-4 1006 test: 1007 suffix: 2d_p1_0 1008 requires: triangle 1009 args: -run_type test -bc_type dirichlet -dm_plex_interpolate 0 -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1010 1011 test: 1012 suffix: 2d_p1_1 1013 requires: triangle 1014 args: -run_type test -bc_type dirichlet -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1015 1016 test: 1017 suffix: 2d_p1_2 1018 requires: triangle 1019 args: -run_type test -dm_refine_volume_limit_pre 0.0625 -bc_type dirichlet -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1020 1021 test: 1022 suffix: 2d_p1_neumann_0 1023 requires: triangle 1024 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 1025 1026 test: 1027 suffix: 2d_p1_neumann_1 1028 requires: triangle 1029 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 1030 1031 # 2D serial P2 test 5-8 1032 test: 1033 suffix: 2d_p2_0 1034 requires: triangle 1035 args: -run_type test -bc_type dirichlet -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1036 1037 test: 1038 suffix: 2d_p2_1 1039 requires: triangle 1040 args: -run_type test -dm_refine_volume_limit_pre 0.0625 -bc_type dirichlet -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1041 1042 test: 1043 suffix: 2d_p2_neumann_0 1044 requires: triangle 1045 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 1046 1047 test: 1048 suffix: 2d_p2_neumann_1 1049 requires: triangle 1050 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 1051 1052 test: 1053 suffix: bd_int_0 1054 requires: triangle 1055 args: -run_type test -bc_type dirichlet -petscspace_degree 2 -bd_integral -dm_view -quiet 1056 1057 test: 1058 suffix: bd_int_1 1059 requires: triangle 1060 args: -run_type test -dm_refine 2 -bc_type dirichlet -petscspace_degree 2 -bd_integral -dm_view -quiet 1061 1062 # 3D serial P1 test 9-12 1063 test: 1064 suffix: 3d_p1_0 1065 requires: ctetgen 1066 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 1067 1068 test: 1069 suffix: 3d_p1_1 1070 requires: ctetgen 1071 args: -run_type test -dm_plex_dim 3 -bc_type dirichlet -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 -dm_view 1072 1073 test: 1074 suffix: 3d_p1_2 1075 requires: ctetgen 1076 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 1077 1078 test: 1079 suffix: 3d_p1_neumann_0 1080 requires: ctetgen 1081 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 1082 1083 # Analytic variable coefficient 13-20 1084 test: 1085 suffix: 13 1086 requires: triangle 1087 args: -run_type test -variable_coefficient analytic -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1088 test: 1089 suffix: 14 1090 requires: triangle 1091 args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient analytic -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1092 test: 1093 suffix: 15 1094 requires: triangle 1095 args: -run_type test -variable_coefficient analytic -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1096 test: 1097 suffix: 16 1098 requires: triangle 1099 args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient analytic -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1100 test: 1101 suffix: 17 1102 requires: ctetgen 1103 args: -run_type test -dm_plex_dim 3 -variable_coefficient analytic -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1104 1105 test: 1106 suffix: 18 1107 requires: ctetgen 1108 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 1109 1110 test: 1111 suffix: 19 1112 requires: ctetgen 1113 args: -run_type test -dm_plex_dim 3 -variable_coefficient analytic -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1114 1115 test: 1116 suffix: 20 1117 requires: ctetgen 1118 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 1119 1120 # P1 variable coefficient 21-28 1121 test: 1122 suffix: 21 1123 requires: triangle 1124 args: -run_type test -variable_coefficient field -petscspace_degree 1 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1125 1126 test: 1127 suffix: 22 1128 requires: triangle 1129 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 1130 1131 test: 1132 suffix: 23 1133 requires: triangle 1134 args: -run_type test -variable_coefficient field -petscspace_degree 2 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1135 1136 test: 1137 suffix: 24 1138 requires: triangle 1139 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 1140 1141 test: 1142 suffix: 25 1143 requires: ctetgen 1144 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 1145 1146 test: 1147 suffix: 26 1148 requires: ctetgen 1149 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 1150 1151 test: 1152 suffix: 27 1153 requires: ctetgen 1154 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 1155 1156 test: 1157 suffix: 28 1158 requires: ctetgen 1159 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 1160 1161 # P0 variable coefficient 29-36 1162 test: 1163 suffix: 29 1164 requires: triangle 1165 args: -run_type test -variable_coefficient field -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1166 1167 test: 1168 suffix: 30 1169 requires: triangle 1170 args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient field -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1171 1172 test: 1173 suffix: 31 1174 requires: triangle 1175 args: -run_type test -variable_coefficient field -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1176 1177 test: 1178 requires: triangle 1179 suffix: 32 1180 args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient field -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1181 1182 test: 1183 requires: ctetgen 1184 suffix: 33 1185 args: -run_type test -dm_plex_dim 3 -variable_coefficient field -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 1186 1187 test: 1188 suffix: 34 1189 requires: ctetgen 1190 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 1191 1192 test: 1193 suffix: 35 1194 requires: ctetgen 1195 args: -run_type test -dm_plex_dim 3 -variable_coefficient field -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1196 1197 test: 1198 suffix: 36 1199 requires: ctetgen 1200 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 1201 1202 # Full solve 39-44 1203 test: 1204 suffix: 39 1205 requires: triangle !single 1206 args: -run_type full -dm_refine_volume_limit_pre 0.015625 -petscspace_degree 2 -pc_type gamg -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 -snes_rtol 1.0e-6 -ksp_rtol 1.0e-7 -ksp_monitor -ksp_converged_reason -snes_monitor_short -snes_converged_reason ::ascii_info_detail 1207 test: 1208 suffix: 40 1209 requires: triangle !single 1210 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 1211 test: 1212 suffix: 41 1213 requires: triangle !single 1214 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 1215 test: 1216 suffix: 42 1217 requires: triangle !single 1218 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 1219 test: 1220 suffix: 43 1221 requires: triangle !single 1222 nsize: 2 1223 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 1224 1225 test: 1226 suffix: 44 1227 requires: triangle !single 1228 nsize: 2 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 -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 1230 1231 # These tests use a loose tolerance just to exercise the PtAP operations for MATIS and multiple PCBDDC setup calls inside PCMG 1232 testset: 1233 requires: triangle !single 1234 nsize: 3 1235 args: -run_type full -petscspace_degree 1 -dm_mat_type is -pc_type mg -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 1236 test: 1237 suffix: gmg_bddc 1238 filter: sed -e "s/CONVERGED_FNORM_RELATIVE iterations 3/CONVERGED_FNORM_RELATIVE iterations 4/g" 1239 args: -mg_levels_pc_type jacobi 1240 test: 1241 filter: sed -e "s/iterations [0-4]/iterations 4/g" 1242 suffix: gmg_bddc_lev 1243 args: -mg_levels_pc_type bddc 1244 1245 # Restarting 1246 testset: 1247 suffix: restart 1248 requires: hdf5 triangle !complex 1249 args: -run_type test -bc_type dirichlet -petscspace_degree 1 1250 test: 1251 args: -dm_view hdf5:sol.h5 -vec_view hdf5:sol.h5::append 1252 test: 1253 args: -dm_plex_filename sol.h5 -dm_plex_name box -restart 1254 1255 # Periodicity 1256 test: 1257 suffix: periodic_0 1258 requires: triangle 1259 args: -run_type full -bc_type dirichlet -petscspace_degree 1 -snes_converged_reason ::ascii_info_detail 1260 1261 test: 1262 requires: !complex 1263 suffix: periodic_1 1264 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 1265 1266 # 2D serial P1 test with field bc 1267 test: 1268 suffix: field_bc_2d_p1_0 1269 requires: triangle 1270 args: -run_type test -bc_type dirichlet -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1271 1272 test: 1273 suffix: field_bc_2d_p1_1 1274 requires: triangle 1275 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 1276 1277 test: 1278 suffix: field_bc_2d_p1_neumann_0 1279 requires: triangle 1280 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 1281 1282 test: 1283 suffix: field_bc_2d_p1_neumann_1 1284 requires: triangle 1285 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 1286 1287 # 3D serial P1 test with field bc 1288 test: 1289 suffix: field_bc_3d_p1_0 1290 requires: ctetgen 1291 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 1292 1293 test: 1294 suffix: field_bc_3d_p1_1 1295 requires: ctetgen 1296 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 1297 1298 test: 1299 suffix: field_bc_3d_p1_neumann_0 1300 requires: ctetgen 1301 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 1302 1303 test: 1304 suffix: field_bc_3d_p1_neumann_1 1305 requires: ctetgen 1306 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 1307 1308 # 2D serial P2 test with field bc 1309 test: 1310 suffix: field_bc_2d_p2_0 1311 requires: triangle 1312 args: -run_type test -bc_type dirichlet -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1 1313 1314 test: 1315 suffix: field_bc_2d_p2_1 1316 requires: triangle 1317 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 1318 1319 test: 1320 suffix: field_bc_2d_p2_neumann_0 1321 requires: triangle 1322 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 1323 1324 test: 1325 suffix: field_bc_2d_p2_neumann_1 1326 requires: triangle 1327 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 1328 1329 # 3D serial P2 test with field bc 1330 test: 1331 suffix: field_bc_3d_p2_0 1332 requires: ctetgen 1333 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 1334 1335 test: 1336 suffix: field_bc_3d_p2_1 1337 requires: ctetgen 1338 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 1339 1340 test: 1341 suffix: field_bc_3d_p2_neumann_0 1342 requires: ctetgen 1343 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 1344 1345 test: 1346 suffix: field_bc_3d_p2_neumann_1 1347 requires: ctetgen 1348 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 1349 1350 # Full solve simplex: Convergence 1351 test: 1352 suffix: 3d_p1_conv 1353 requires: ctetgen 1354 args: -run_type full -dm_plex_dim 3 -dm_refine 1 -bc_type dirichlet -petscspace_degree 1 \ 1355 -snes_convergence_estimate -convest_num_refine 1 -pc_type lu 1356 1357 # Full solve simplex: PCBDDC 1358 test: 1359 suffix: tri_bddc 1360 requires: triangle !single 1361 nsize: 5 1362 args: -run_type full -petscpartitioner_type simple -dm_refine 2 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -dm_mat_type is -pc_type bddc -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0 1363 1364 # Full solve simplex: PCBDDC 1365 test: 1366 suffix: tri_parmetis_bddc 1367 requires: triangle !single parmetis 1368 nsize: 4 1369 args: -run_type full -petscpartitioner_type parmetis -dm_refine 2 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -dm_mat_type is -pc_type bddc -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0 1370 1371 testset: 1372 args: -run_type full -dm_plex_simplex 0 -dm_plex_box_faces 3,3 -petscpartitioner_type simple -dm_refine 2 -bc_type dirichlet -petscspace_degree 2 -dm_mat_type is -pc_type bddc -ksp_type gmres -snes_monitor_short -ksp_monitor_short -snes_view -petscspace_poly_tensor -pc_bddc_corner_selection -ksp_rtol 1.e-9 -pc_bddc_use_edges 0 1373 nsize: 5 1374 output_file: output/ex12_quad_bddc.out 1375 filter: sed -e "s/aijcusparse/aij/g" -e "s/aijviennacl/aij/g" -e "s/factorization: cusparse/factorization: petsc/g" 1376 test: 1377 requires: !single 1378 suffix: quad_bddc 1379 test: 1380 requires: !single cuda 1381 suffix: quad_bddc_cuda 1382 args: -matis_localmat_type aijcusparse -pc_bddc_dirichlet_pc_factor_mat_solver_type cusparse -pc_bddc_neumann_pc_factor_mat_solver_type cusparse 1383 test: 1384 requires: !single viennacl 1385 suffix: quad_bddc_viennacl 1386 args: -matis_localmat_type aijviennacl 1387 1388 # Full solve simplex: ASM 1389 test: 1390 suffix: tri_q2q1_asm_lu 1391 requires: triangle !single 1392 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 1393 1394 test: 1395 suffix: tri_q2q1_msm_lu 1396 requires: triangle !single 1397 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 1398 1399 test: 1400 suffix: tri_q2q1_asm_sor 1401 requires: triangle !single 1402 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 1403 1404 test: 1405 suffix: tri_q2q1_msm_sor 1406 requires: triangle !single 1407 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 1408 1409 # Full solve simplex: FAS 1410 test: 1411 suffix: fas_newton_0 1412 requires: triangle !single 1413 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 1414 1415 test: 1416 suffix: fas_newton_1 1417 requires: triangle !single 1418 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 1419 filter: sed -e "s/total number of linear solver iterations=14/total number of linear solver iterations=15/g" 1420 1421 test: 1422 suffix: fas_ngs_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 ngs -fas_levels_1_snes_monitor_short 1425 1426 # These two tests are broken because DMPlexComputeInjectorFEM() only works for regularly refined meshes 1427 test: 1428 suffix: fas_newton_coarse_0 1429 requires: pragmatic triangle 1430 TODO: broken 1431 args: -run_type full -variable_coefficient nonlinear -petscspace_degree 1 \ 1432 -dm_refine 2 -dm_coarsen_hierarchy 1 -dm_plex_hash_location -dm_adaptor pragmatic \ 1433 -snes_type fas -snes_fas_levels 2 -snes_converged_reason ::ascii_info_detail -snes_monitor_short -snes_view \ 1434 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -fas_coarse_snes_linesearch_type basic \ 1435 -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 1436 1437 test: 1438 suffix: mg_newton_coarse_0 1439 requires: triangle pragmatic 1440 TODO: broken 1441 args: -run_type full -petscspace_degree 1 \ 1442 -dm_refine 3 -dm_coarsen_hierarchy 3 -dm_plex_hash_location -dm_adaptor pragmatic \ 1443 -snes_atol 1.0e-8 -snes_rtol 0.0 -snes_monitor_short -snes_converged_reason ::ascii_info_detail -snes_view \ 1444 -ksp_type richardson -ksp_atol 1.0e-8 -ksp_rtol 0.0 -ksp_norm_type unpreconditioned -ksp_monitor_true_residual \ 1445 -pc_type mg -pc_mg_levels 4 \ 1446 -mg_levels_ksp_type gmres -mg_levels_pc_type ilu -mg_levels_ksp_max_it 10 1447 1448 # Full solve tensor 1449 test: 1450 suffix: tensor_plex_2d 1451 args: -run_type test -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -dm_refine_hierarchy 2 1452 1453 test: 1454 suffix: tensor_p4est_2d 1455 requires: p4est 1456 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 1457 1458 test: 1459 suffix: tensor_plex_3d 1460 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 1461 1462 test: 1463 suffix: tensor_p4est_3d 1464 requires: p4est 1465 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 1466 1467 test: 1468 suffix: p4est_test_q2_conformal_serial 1469 requires: p4est 1470 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 1471 1472 test: 1473 suffix: p4est_test_q2_conformal_parallel 1474 requires: p4est 1475 nsize: 7 1476 args: -run_type test -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -petscpartitioner_type simple 1477 1478 test: 1479 suffix: p4est_test_q2_conformal_parallel_parmetis 1480 requires: parmetis p4est 1481 nsize: 4 1482 args: -run_type test -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -petscpartitioner_type parmetis 1483 1484 test: 1485 suffix: p4est_test_q2_nonconformal_serial 1486 requires: p4est 1487 filter: grep -v "CG or CGNE: variant" 1488 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 1489 1490 test: 1491 suffix: p4est_test_q2_nonconformal_parallel 1492 requires: p4est 1493 filter: grep -v "CG or CGNE: variant" 1494 nsize: 7 1495 args: -run_type test -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple 1496 1497 test: 1498 suffix: p4est_test_q2_nonconformal_parallel_parmetis 1499 requires: parmetis p4est 1500 nsize: 4 1501 args: -run_type test -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type parmetis 1502 1503 test: 1504 suffix: p4est_exact_q2_conformal_serial 1505 requires: p4est !single !complex !__float128 1506 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 1507 1508 test: 1509 suffix: p4est_exact_q2_conformal_parallel 1510 requires: p4est !single !complex !__float128 1511 nsize: 4 1512 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 1513 1514 test: 1515 suffix: p4est_exact_q2_conformal_parallel_parmetis 1516 requires: parmetis p4est !single 1517 nsize: 4 1518 args: -run_type exact -petscspace_degree 2 -fas_levels_snes_atol 1.e-10 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -petscpartitioner_type parmetis 1519 1520 test: 1521 suffix: p4est_exact_q2_nonconformal_serial 1522 requires: p4est 1523 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 1524 1525 test: 1526 suffix: p4est_exact_q2_nonconformal_parallel 1527 requires: p4est 1528 nsize: 7 1529 args: -run_type exact -petscspace_degree 2 -fas_levels_snes_atol 1.e-10 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple 1530 1531 test: 1532 suffix: p4est_exact_q2_nonconformal_parallel_parmetis 1533 requires: parmetis p4est 1534 nsize: 4 1535 args: -run_type exact -petscspace_degree 2 -fas_levels_snes_atol 1.e-10 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type parmetis 1536 1537 test: 1538 suffix: p4est_full_q2_nonconformal_serial 1539 requires: p4est !single 1540 filter: grep -v "variant HERMITIAN" 1541 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 1542 1543 test: 1544 suffix: p4est_full_q2_nonconformal_parallel 1545 requires: p4est !single 1546 filter: grep -v "variant HERMITIAN" 1547 nsize: 7 1548 args: -run_type full -petscspace_degree 2 -snes_max_it 20 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type jacobi -fas_coarse_ksp_type cg -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type jacobi -fas_levels_ksp_type cg -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple 1549 1550 test: 1551 suffix: p4est_full_q2_nonconformal_parallel_bddcfas 1552 requires: p4est !single 1553 filter: grep -v "variant HERMITIAN" 1554 nsize: 7 1555 args: -run_type full -petscspace_degree 2 -snes_max_it 20 -snes_type fas -snes_fas_levels 3 -dm_mat_type is -fas_coarse_pc_type bddc -fas_coarse_ksp_type cg -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type bddc -fas_levels_ksp_type cg -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple 1556 1557 test: 1558 suffix: p4est_full_q2_nonconformal_parallel_bddc 1559 requires: p4est !single 1560 filter: grep -v "variant HERMITIAN" 1561 nsize: 7 1562 args: -run_type full -petscspace_degree 2 -snes_max_it 20 -snes_type newtonls -dm_mat_type is -pc_type bddc -ksp_type cg -snes_monitor_short -snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple 1563 1564 test: 1565 TODO: broken 1566 suffix: p4est_fas_q2_conformal_serial 1567 requires: p4est !complex !__float128 1568 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 1569 1570 test: 1571 TODO: broken 1572 suffix: p4est_fas_q2_nonconformal_serial 1573 requires: p4est 1574 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 1575 1576 test: 1577 suffix: fas_newton_0_p4est 1578 requires: p4est !single !__float128 1579 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 1580 1581 # Full solve simplicial AMR 1582 test: 1583 suffix: tri_p1_adapt_init_pragmatic 1584 requires: pragmatic 1585 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 1586 1587 test: 1588 suffix: tri_p2_adapt_init_pragmatic 1589 requires: pragmatic 1590 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 1591 1592 test: 1593 suffix: tri_p1_adapt_init_mmg 1594 requires: mmg 1595 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 1596 1597 test: 1598 suffix: tri_p2_adapt_init_mmg 1599 requires: mmg 1600 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 1601 1602 test: 1603 suffix: tri_p1_adapt_seq_pragmatic 1604 requires: pragmatic 1605 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 1606 1607 test: 1608 suffix: tri_p2_adapt_seq_pragmatic 1609 requires: pragmatic 1610 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 1611 1612 test: 1613 suffix: tri_p1_adapt_seq_mmg 1614 requires: mmg 1615 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 1616 1617 test: 1618 suffix: tri_p2_adapt_seq_mmg 1619 requires: mmg 1620 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 1621 1622 test: 1623 suffix: tri_p1_adapt_analytic_pragmatic 1624 requires: pragmatic 1625 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 1626 1627 test: 1628 suffix: tri_p2_adapt_analytic_pragmatic 1629 requires: pragmatic 1630 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 1631 1632 test: 1633 suffix: tri_p1_adapt_analytic_mmg 1634 requires: mmg 1635 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 1636 1637 test: 1638 suffix: tri_p2_adapt_analytic_mmg 1639 requires: mmg 1640 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 1641 1642 test: 1643 suffix: tri_p1_adapt_uniform_pragmatic 1644 requires: pragmatic tetgen 1645 nsize: 2 1646 args: -run_type full -dm_plex_box_faces 8,8,8 -bc_type dirichlet -petscspace_degree 1 -variable_coefficient none -snes_converged_reason ::ascii_info_detail -ksp_type cg -pc_type sor -snes_adapt_sequence 3 -adaptor_target_num 400 -dm_plex_metric_h_max 0.5 -dm_plex_dim 3 -dm_adaptor pragmatic 1647 timeoutfactor: 2 1648 1649 test: 1650 suffix: tri_p2_adapt_uniform_pragmatic 1651 requires: pragmatic tetgen 1652 nsize: 2 1653 args: -run_type full -dm_plex_box_faces 8,8,8 -bc_type dirichlet -petscspace_degree 2 -variable_coefficient none -snes_converged_reason ::ascii_info_detail -ksp_type cg -pc_type sor -snes_adapt_sequence 1 -adaptor_target_num 400 -dm_plex_metric_h_max 0.5 -dm_plex_dim 3 -dm_adaptor pragmatic 1654 timeoutfactor: 1 1655 1656 test: 1657 suffix: tri_p1_adapt_uniform_mmg 1658 requires: mmg tetgen 1659 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 1660 timeoutfactor: 2 1661 1662 test: 1663 suffix: tri_p2_adapt_uniform_mmg 1664 requires: mmg tetgen 1665 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 1666 timeoutfactor: 1 1667 1668 test: 1669 suffix: tri_p1_adapt_uniform_parmmg 1670 requires: parmmg tetgen 1671 nsize: 2 1672 args: -run_type full -dm_plex_box_faces 8,8,8 -bc_type dirichlet -petscspace_degree 1 -variable_coefficient none -snes_converged_reason ::ascii_info_detail -ksp_type cg -pc_type sor -snes_adapt_sequence 3 -adaptor_target_num 400 -dm_plex_metric_h_max 0.5 -dm_plex_dim 3 -dm_adaptor parmmg 1673 timeoutfactor: 2 1674 1675 test: 1676 suffix: tri_p2_adapt_uniform_parmmg 1677 requires: parmmg tetgen 1678 nsize: 2 1679 args: -run_type full -dm_plex_box_faces 8,8,8 -bc_type dirichlet -petscspace_degree 2 -variable_coefficient none -snes_converged_reason ::ascii_info_detail -ksp_type cg -pc_type sor -snes_adapt_sequence 1 -adaptor_target_num 400 -dm_plex_metric_h_max 0.5 -dm_plex_dim 3 -dm_adaptor parmmg 1680 timeoutfactor: 1 1681 1682 # Full solve tensor AMR 1683 test: 1684 suffix: quad_q1_adapt_0 1685 requires: p4est 1686 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 1687 filter: grep -v DM_ 1688 1689 test: 1690 suffix: amr_0 1691 nsize: 5 1692 args: -run_type test -petscpartitioner_type simple -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -dm_refine 1 1693 1694 test: 1695 suffix: amr_1 1696 requires: p4est !complex 1697 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 1698 1699 test: 1700 suffix: p4est_solve_bddc 1701 requires: p4est !complex 1702 args: -run_type full -variable_coefficient nonlinear -nonzero_initial_guess 1 -petscspace_degree 2 -snes_max_it 20 -snes_type newtonls -dm_mat_type is -pc_type bddc -ksp_type cg -snes_monitor_short -ksp_monitor -snes_linesearch_type bt -snes_converged_reason -snes_view -dm_plex_simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple -pc_bddc_detect_disconnected 1703 nsize: 4 1704 1705 test: 1706 suffix: p4est_solve_fas 1707 requires: p4est 1708 args: -run_type full -variable_coefficient nonlinear -nonzero_initial_guess 1 -petscspace_degree 2 -snes_max_it 10 -snes_type fas -snes_linesearch_type bt -snes_fas_levels 3 -fas_coarse_snes_type newtonls -fas_coarse_snes_linesearch_type basic -fas_coarse_ksp_type cg -fas_coarse_pc_type jacobi -fas_coarse_snes_monitor_short -fas_levels_snes_max_it 4 -fas_levels_snes_type newtonls -fas_levels_snes_linesearch_type bt -fas_levels_ksp_type cg -fas_levels_pc_type jacobi -fas_levels_snes_monitor_short -fas_levels_cycle_snes_linesearch_type bt -snes_monitor_short -snes_converged_reason -snes_view -dm_plex_simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash 1709 nsize: 4 1710 TODO: identical machine two runs produce slightly different solver trackers 1711 1712 test: 1713 suffix: p4est_convergence_test_1 1714 requires: p4est 1715 args: -quiet -run_type test -petscspace_degree 1 -dm_plex_simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 2 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash 1716 nsize: 4 1717 1718 test: 1719 suffix: p4est_convergence_test_2 1720 requires: p4est 1721 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 1722 1723 test: 1724 suffix: p4est_convergence_test_3 1725 requires: p4est 1726 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 1727 1728 test: 1729 suffix: p4est_convergence_test_4 1730 requires: p4est 1731 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 1732 timeoutfactor: 5 1733 1734 # Serial tests with GLVis visualization 1735 test: 1736 suffix: glvis_2d_tet_p1 1737 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 1738 test: 1739 suffix: glvis_2d_tet_p2 1740 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 1741 test: 1742 suffix: glvis_2d_hex_p1 1743 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 1744 test: 1745 suffix: glvis_2d_hex_p2 1746 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 1747 test: 1748 suffix: glvis_2d_hex_p2_p4est 1749 requires: p4est 1750 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 1751 test: 1752 suffix: glvis_2d_tet_p0 1753 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 1754 test: 1755 suffix: glvis_2d_hex_p0 1756 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 1757 1758 # PCHPDDM tests 1759 testset: 1760 nsize: 4 1761 requires: hpddm slepc !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES) 1762 args: -run_type test -run_test_check_ksp -quiet -petscspace_degree 1 -petscpartitioner_type simple -bc_type none -dm_plex_simplex 0 -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_nev 2 -pc_hpddm_coarse_p 1 -pc_hpddm_coarse_pc_type svd -ksp_rtol 1.e-10 -pc_hpddm_levels_1_st_pc_factor_shift_type INBLOCKS -ksp_converged_reason 1763 test: 1764 suffix: quad_singular_hpddm 1765 args: -dm_plex_box_faces 6,7 1766 test: 1767 requires: p4est 1768 suffix: p4est_singular_2d_hpddm 1769 args: -dm_plex_convert_type p4est -dm_forest_minimum_refinement 1 -dm_forest_initial_refinement 3 -dm_forest_maximum_refinement 3 1770 test: 1771 requires: p4est 1772 suffix: p4est_nc_singular_2d_hpddm 1773 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 1774 testset: 1775 nsize: 4 1776 requires: hpddm slepc triangle !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES) 1777 args: -run_type full -petscpartitioner_type simple -dm_refine 2 -bc_type dirichlet -petscspace_degree 2 -ksp_type gmres -ksp_gmres_restart 100 -pc_type hpddm -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0 -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_nev 4 -pc_hpddm_coarse_p 2 -pc_hpddm_coarse_pc_type redundant -ksp_rtol 1.e-1 1778 test: 1779 args: -pc_hpddm_coarse_mat_type baij -options_left no 1780 suffix: tri_hpddm_reuse_baij 1781 test: 1782 requires: !complex 1783 suffix: tri_hpddm_reuse 1784 testset: 1785 nsize: 4 1786 requires: hpddm slepc !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES) 1787 args: -run_type full -petscpartitioner_type simple -dm_plex_box_faces 7,5 -dm_refine 2 -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 2 -ksp_type gmres -ksp_gmres_restart 100 -pc_type hpddm -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0 -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_nev 4 -pc_hpddm_coarse_p 2 -pc_hpddm_coarse_pc_type redundant -ksp_rtol 1.e-1 1788 test: 1789 args: -pc_hpddm_coarse_mat_type baij -options_left no 1790 suffix: quad_hpddm_reuse_baij 1791 test: 1792 requires: !complex 1793 suffix: quad_hpddm_reuse 1794 testset: 1795 nsize: 4 1796 requires: hpddm slepc !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES) 1797 args: -run_type full -petscpartitioner_type simple -dm_plex_box_faces 7,5 -dm_refine 2 -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -pc_type hpddm -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0 -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_threshold 0.1 -pc_hpddm_coarse_p 2 -pc_hpddm_coarse_pc_type redundant -ksp_rtol 1.e-1 1798 test: 1799 args: -pc_hpddm_coarse_mat_type baij -options_left no 1800 suffix: quad_hpddm_reuse_threshold_baij 1801 test: 1802 requires: !complex 1803 suffix: quad_hpddm_reuse_threshold 1804 testset: 1805 nsize: 4 1806 requires: hpddm slepc parmetis !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES) 1807 filter: sed -e "s/linear solver iterations=17/linear solver iterations=16/g" 1808 args: -run_type full -petscpartitioner_type parmetis -dm_refine 3 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -pc_type hpddm -snes_monitor_short -snes_converged_reason ::ascii_info_detail -snes_view -show_solution 0 -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type icc -pc_hpddm_levels_1_eps_nev 20 -pc_hpddm_coarse_p 2 -pc_hpddm_coarse_pc_type redundant -ksp_rtol 1.e-10 -dm_plex_filename ${PETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_boundary_label marker -pc_hpddm_levels_1_sub_pc_factor_levels 3 -variable_coefficient ball -dm_plex_gmsh_periodic 0 1809 test: 1810 args: -pc_hpddm_coarse_mat_type baij -options_left no 1811 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" 1812 suffix: tri_parmetis_hpddm_baij 1813 test: 1814 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" 1815 requires: !complex 1816 suffix: tri_parmetis_hpddm 1817 1818 # 2D serial P1 tests for adaptive MG 1819 test: 1820 suffix: 2d_p1_adaptmg_0 1821 requires: triangle 1822 args: -petscpartitioner_type simple -dm_refine_hierarchy 3 -dm_plex_box_faces 4,4 -bc_type dirichlet -petscspace_degree 1 \ 1823 -variable_coefficient checkerboard_0 -mat_petscspace_degree 0 -div 16 -k 3 \ 1824 -snes_max_it 1 -ksp_converged_reason \ 1825 -ksp_rtol 1e-8 -pc_type mg 1826 test: 1827 suffix: 2d_p1_adaptmg_1 1828 requires: triangle bamg todo 1829 args: -petscpartitioner_type simple -dm_refine_hierarchy 3 -dm_plex_box_faces 4,4 -bc_type dirichlet -petscspace_degree 1 \ 1830 -variable_coefficient checkerboard_0 -mat_petscspace_degree 0 -div 16 -k 3 \ 1831 -snes_max_it 1 -ksp_converged_reason \ 1832 -ksp_rtol 1e-8 -pc_type mg -pc_mg_galerkin -pc_mg_adapt_interp_coarse_space eigenvector -pc_mg_adapt_interp_n 1 \ 1833 -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 1834 test: 1835 suffix: 2d_p1_adaptmg_gdsw 1836 requires: triangle 1837 nsize: 4 1838 args: -petscpartitioner_type simple -dm_refine 3 -dm_plex_box_faces 4,4 -bc_type dirichlet -petscspace_degree 1 \ 1839 -variable_coefficient checkerboard_0 -mat_petscspace_degree 0 -div 16 -k 3 \ 1840 -snes_max_it 1 -ksp_converged_reason \ 1841 -ksp_rtol 1e-8 -pc_type mg -pc_mg_galerkin -pc_mg_adapt_interp_coarse_space gdsw -pc_mg_levels 2 -mg_levels_pc_type asm -dm_mat_type {{aij is}} 1842 1843 test: 1844 suffix: 2d_p1_adaptmg_agdsw 1845 requires: triangle mumps 1846 nsize: 4 1847 args: -petscpartitioner_type simple -dm_refine 3 -dm_plex_box_faces 4,4 -bc_type dirichlet -petscspace_degree 1 \ 1848 -variable_coefficient checkerboard_0 -mat_petscspace_degree 0 -div 16 -k 3 \ 1849 -snes_max_it 1 -ksp_converged_reason \ 1850 -ksp_rtol 1e-8 -pc_type mg -pc_mg_galerkin -pc_mg_adapt_interp_coarse_space gdsw -pc_mg_levels 2 -mg_levels_pc_type asm -dm_mat_type is -mg_levels_gdsw_tolerance 0.1 -mg_levels_gdsw_pseudo_pc_type qr 1851 1852 TEST*/ 1853