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