1 static char help[] = "Check that a DM can accurately represent and interpolate functions of a given polynomial order\n\n"; 2 3 #include <petscdmplex.h> 4 #include <petscdm.h> 5 #include <petscdmda.h> 6 #include <petscfe.h> 7 #include <petscds.h> 8 #include <petscksp.h> 9 #include <petscsnes.h> 10 11 typedef struct { 12 /* Domain and mesh definition */ 13 PetscBool useDA; /* Flag DMDA tensor product mesh */ 14 PetscBool shearCoords; /* Flag for shear transform */ 15 PetscBool nonaffineCoords; /* Flag for non-affine transform */ 16 /* Element definition */ 17 PetscInt qorder; /* Order of the quadrature */ 18 PetscInt numComponents; /* Number of field components */ 19 PetscFE fe; /* The finite element */ 20 /* Testing space */ 21 PetscInt porder; /* Order of polynomials to test */ 22 PetscBool convergence; /* Test for order of convergence */ 23 PetscBool convRefine; /* Test for convergence using refinement, otherwise use coarsening */ 24 PetscBool constraints; /* Test local constraints */ 25 PetscBool tree; /* Test tree routines */ 26 PetscBool testFEjacobian; /* Test finite element Jacobian assembly */ 27 PetscBool testFVgrad; /* Test finite difference gradient routine */ 28 PetscBool testInjector; /* Test finite element injection routines */ 29 PetscInt treeCell; /* Cell to refine in tree test */ 30 PetscReal constants[3]; /* Constant values for each dimension */ 31 } AppCtx; 32 33 /* u = 1 */ 34 PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx) 35 { 36 AppCtx *user = (AppCtx *) ctx; 37 PetscInt d; 38 for (d = 0; d < dim; ++d) u[d] = user->constants[d]; 39 return 0; 40 } 41 PetscErrorCode constantDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) 42 { 43 PetscInt d; 44 for (d = 0; d < dim; ++d) u[d] = 0.0; 45 return 0; 46 } 47 48 /* u = x */ 49 PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx) 50 { 51 PetscInt d; 52 for (d = 0; d < dim; ++d) u[d] = coords[d]; 53 return 0; 54 } 55 PetscErrorCode linearDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) 56 { 57 PetscInt d, e; 58 for (d = 0; d < dim; ++d) { 59 u[d] = 0.0; 60 for (e = 0; e < dim; ++e) u[d] += (d == e ? 1.0 : 0.0) * n[e]; 61 } 62 return 0; 63 } 64 65 /* u = x^2 or u = (x^2, xy) or u = (xy, yz, zx) */ 66 PetscErrorCode quadratic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx) 67 { 68 if (dim > 2) {u[0] = coords[0]*coords[1]; u[1] = coords[1]*coords[2]; u[2] = coords[2]*coords[0];} 69 else if (dim > 1) {u[0] = coords[0]*coords[0]; u[1] = coords[0]*coords[1];} 70 else if (dim > 0) {u[0] = coords[0]*coords[0];} 71 return 0; 72 } 73 PetscErrorCode quadraticDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) 74 { 75 if (dim > 2) {u[0] = coords[1]*n[0] + coords[0]*n[1]; u[1] = coords[2]*n[1] + coords[1]*n[2]; u[2] = coords[2]*n[0] + coords[0]*n[2];} 76 else if (dim > 1) {u[0] = 2.0*coords[0]*n[0]; u[1] = coords[1]*n[0] + coords[0]*n[1];} 77 else if (dim > 0) {u[0] = 2.0*coords[0]*n[0];} 78 return 0; 79 } 80 81 /* u = x^3 or u = (x^3, x^2y) or u = (x^2y, y^2z, z^2x) */ 82 PetscErrorCode cubic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx) 83 { 84 if (dim > 2) {u[0] = coords[0]*coords[0]*coords[1]; u[1] = coords[1]*coords[1]*coords[2]; u[2] = coords[2]*coords[2]*coords[0];} 85 else if (dim > 1) {u[0] = coords[0]*coords[0]*coords[0]; u[1] = coords[0]*coords[0]*coords[1];} 86 else if (dim > 0) {u[0] = coords[0]*coords[0]*coords[0];} 87 return 0; 88 } 89 PetscErrorCode cubicDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) 90 { 91 if (dim > 2) {u[0] = 2.0*coords[0]*coords[1]*n[0] + coords[0]*coords[0]*n[1]; u[1] = 2.0*coords[1]*coords[2]*n[1] + coords[1]*coords[1]*n[2]; u[2] = 2.0*coords[2]*coords[0]*n[2] + coords[2]*coords[2]*n[0];} 92 else if (dim > 1) {u[0] = 3.0*coords[0]*coords[0]*n[0]; u[1] = 2.0*coords[0]*coords[1]*n[0] + coords[0]*coords[0]*n[1];} 93 else if (dim > 0) {u[0] = 3.0*coords[0]*coords[0]*n[0];} 94 return 0; 95 } 96 97 /* u = tanh(x) */ 98 PetscErrorCode trig(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx) 99 { 100 PetscInt d; 101 for (d = 0; d < dim; ++d) u[d] = PetscTanhReal(coords[d] - 0.5); 102 return 0; 103 } 104 PetscErrorCode trigDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) 105 { 106 PetscInt d; 107 for (d = 0; d < dim; ++d) u[d] = 1.0/PetscSqr(PetscCoshReal(coords[d] - 0.5)) * n[d]; 108 return 0; 109 } 110 111 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 112 { 113 PetscInt n = 3; 114 PetscErrorCode ierr; 115 116 PetscFunctionBeginUser; 117 options->useDA = PETSC_FALSE; 118 options->shearCoords = PETSC_FALSE; 119 options->nonaffineCoords = PETSC_FALSE; 120 options->qorder = 0; 121 options->numComponents = PETSC_DEFAULT; 122 options->porder = 0; 123 options->convergence = PETSC_FALSE; 124 options->convRefine = PETSC_TRUE; 125 options->constraints = PETSC_FALSE; 126 options->tree = PETSC_FALSE; 127 options->treeCell = 0; 128 options->testFEjacobian = PETSC_FALSE; 129 options->testFVgrad = PETSC_FALSE; 130 options->testInjector = PETSC_FALSE; 131 options->constants[0] = 1.0; 132 options->constants[1] = 2.0; 133 options->constants[2] = 3.0; 134 135 ierr = PetscOptionsBegin(comm, "", "Projection Test Options", "DMPlex");CHKERRQ(ierr); 136 ierr = PetscOptionsBool("-use_da", "Flag for DMDA mesh", "ex3.c", options->useDA, &options->useDA, NULL);CHKERRQ(ierr); 137 ierr = PetscOptionsBool("-shear_coords", "Transform coordinates with a shear", "ex3.c", options->shearCoords, &options->shearCoords, NULL);CHKERRQ(ierr); 138 ierr = PetscOptionsBool("-non_affine_coords", "Transform coordinates with a non-affine transform", "ex3.c", options->nonaffineCoords, &options->nonaffineCoords, NULL);CHKERRQ(ierr); 139 ierr = PetscOptionsBoundedInt("-qorder", "The quadrature order", "ex3.c", options->qorder, &options->qorder, NULL,0);CHKERRQ(ierr); 140 ierr = PetscOptionsBoundedInt("-num_comp", "The number of field components", "ex3.c", options->numComponents, &options->numComponents, NULL,PETSC_DEFAULT);CHKERRQ(ierr); 141 ierr = PetscOptionsBoundedInt("-porder", "The order of polynomials to test", "ex3.c", options->porder, &options->porder, NULL,0);CHKERRQ(ierr); 142 ierr = PetscOptionsBool("-convergence", "Check the convergence rate", "ex3.c", options->convergence, &options->convergence, NULL);CHKERRQ(ierr); 143 ierr = PetscOptionsBool("-conv_refine", "Use refinement for the convergence rate", "ex3.c", options->convRefine, &options->convRefine, NULL);CHKERRQ(ierr); 144 ierr = PetscOptionsBool("-constraints", "Test local constraints (serial only)", "ex3.c", options->constraints, &options->constraints, NULL);CHKERRQ(ierr); 145 ierr = PetscOptionsBool("-tree", "Test tree routines", "ex3.c", options->tree, &options->tree, NULL);CHKERRQ(ierr); 146 ierr = PetscOptionsBoundedInt("-tree_cell", "cell to refine in tree test", "ex3.c", options->treeCell, &options->treeCell, NULL,0);CHKERRQ(ierr); 147 ierr = PetscOptionsBool("-test_fe_jacobian", "Test finite element Jacobian assembly", "ex3.c", options->testFEjacobian, &options->testFEjacobian, NULL);CHKERRQ(ierr); 148 ierr = PetscOptionsBool("-test_fv_grad", "Test finite volume gradient reconstruction", "ex3.c", options->testFVgrad, &options->testFVgrad, NULL);CHKERRQ(ierr); 149 ierr = PetscOptionsBool("-test_injector","Test finite element injection", "ex3.c", options->testInjector, &options->testInjector,NULL);CHKERRQ(ierr); 150 ierr = PetscOptionsRealArray("-constants","Set the constant values", "ex3.c", options->constants, &n,NULL);CHKERRQ(ierr); 151 ierr = PetscOptionsEnd();CHKERRQ(ierr); 152 153 PetscFunctionReturn(0); 154 } 155 156 static PetscErrorCode TransformCoordinates(DM dm, AppCtx *user) 157 { 158 PetscSection coordSection; 159 Vec coordinates; 160 PetscScalar *coords; 161 PetscInt vStart, vEnd, v; 162 PetscErrorCode ierr; 163 164 PetscFunctionBeginUser; 165 if (user->nonaffineCoords) { 166 /* x' = r^(1/p) (x/r), y' = r^(1/p) (y/r), z' = z */ 167 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 168 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 169 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 170 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 171 for (v = vStart; v < vEnd; ++v) { 172 PetscInt dof, off; 173 PetscReal p = 4.0, r; 174 175 ierr = PetscSectionGetDof(coordSection, v, &dof);CHKERRQ(ierr); 176 ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 177 switch (dof) { 178 case 2: 179 r = PetscSqr(PetscRealPart(coords[off+0])) + PetscSqr(PetscRealPart(coords[off+1])); 180 coords[off+0] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p)/(2*p))*coords[off+0]; 181 coords[off+1] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p)/(2*p))*coords[off+1]; 182 break; 183 case 3: 184 r = PetscSqr(PetscRealPart(coords[off+0])) + PetscSqr(PetscRealPart(coords[off+1])); 185 coords[off+0] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p)/(2*p))*coords[off+0]; 186 coords[off+1] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p)/(2*p))*coords[off+1]; 187 coords[off+2] = coords[off+2]; 188 break; 189 } 190 } 191 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 192 } 193 if (user->shearCoords) { 194 /* x' = x + m y + m z, y' = y + m z, z' = z */ 195 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 196 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 197 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 198 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 199 for (v = vStart; v < vEnd; ++v) { 200 PetscInt dof, off; 201 PetscReal m = 1.0; 202 203 ierr = PetscSectionGetDof(coordSection, v, &dof);CHKERRQ(ierr); 204 ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 205 switch (dof) { 206 case 2: 207 coords[off+0] = coords[off+0] + m*coords[off+1]; 208 coords[off+1] = coords[off+1]; 209 break; 210 case 3: 211 coords[off+0] = coords[off+0] + m*coords[off+1] + m*coords[off+2]; 212 coords[off+1] = coords[off+1] + m*coords[off+2]; 213 coords[off+2] = coords[off+2]; 214 break; 215 } 216 } 217 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 218 } 219 PetscFunctionReturn(0); 220 } 221 222 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 223 { 224 PetscInt dim = 2; 225 PetscBool simplex; 226 PetscErrorCode ierr; 227 228 PetscFunctionBeginUser; 229 if (user->useDA) { 230 switch (dim) { 231 case 2: 232 ierr = DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 2, 2, PETSC_DETERMINE, PETSC_DETERMINE, 1, 1, NULL, NULL, dm);CHKERRQ(ierr); 233 ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 234 ierr = DMSetUp(*dm);CHKERRQ(ierr); 235 ierr = DMDASetVertexCoordinates(*dm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0);CHKERRQ(ierr); 236 break; 237 default: 238 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot create structured mesh of dimension %d", dim); 239 } 240 ierr = PetscObjectSetName((PetscObject) *dm, "Hexahedral Mesh");CHKERRQ(ierr); 241 } else { 242 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 243 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 244 ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 245 246 ierr = DMGetDimension(*dm, &dim);CHKERRQ(ierr); 247 ierr = DMPlexIsSimplex(*dm, &simplex);CHKERRQ(ierr); 248 ierr = MPI_Bcast(&simplex, 1, MPIU_BOOL, 0, comm);CHKERRMPI(ierr); 249 if (user->tree) { 250 DM refTree, ncdm = NULL; 251 252 ierr = DMPlexCreateDefaultReferenceTree(comm,dim,simplex,&refTree);CHKERRQ(ierr); 253 ierr = DMPlexSetReferenceTree(*dm,refTree);CHKERRQ(ierr); 254 ierr = DMDestroy(&refTree);CHKERRQ(ierr); 255 ierr = DMPlexTreeRefineCell(*dm,user->treeCell,&ncdm);CHKERRQ(ierr); 256 if (ncdm) { 257 ierr = DMDestroy(dm);CHKERRQ(ierr); 258 *dm = ncdm; 259 ierr = DMPlexSetRefinementUniform(*dm, PETSC_FALSE);CHKERRQ(ierr); 260 } 261 ierr = DMViewFromOptions(*dm,NULL,"-tree_dm_view");CHKERRQ(ierr); 262 } else { 263 ierr = DMPlexSetRefinementUniform(*dm, PETSC_TRUE);CHKERRQ(ierr); 264 } 265 ierr = PetscObjectSetOptionsPrefix((PetscObject) *dm, "dist_");CHKERRQ(ierr); 266 ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 267 ierr = PetscObjectSetOptionsPrefix((PetscObject) *dm, NULL);CHKERRQ(ierr); 268 if (simplex) {ierr = PetscObjectSetName((PetscObject) *dm, "Simplicial Mesh");CHKERRQ(ierr);} 269 else {ierr = PetscObjectSetName((PetscObject) *dm, "Hexahedral Mesh");CHKERRQ(ierr);} 270 } 271 ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 272 ierr = TransformCoordinates(*dm, user);CHKERRQ(ierr); 273 ierr = DMViewFromOptions(*dm,NULL,"-dm_view");CHKERRQ(ierr); 274 PetscFunctionReturn(0); 275 } 276 277 static void simple_mass(PetscInt dim, PetscInt Nf, PetscInt NfAux, 278 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 279 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 280 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 281 { 282 PetscInt d, e; 283 for (d = 0, e = 0; d < dim; d++, e+=dim+1) { 284 g0[e] = 1.; 285 } 286 } 287 288 /* < \nabla v, 1/2(\nabla u + {\nabla u}^T) > */ 289 static void symmetric_gradient_inner_product(PetscInt dim, PetscInt Nf, PetscInt NfAux, 290 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 291 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 292 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar C[]) 293 { 294 PetscInt compI, compJ, d, e; 295 296 for (compI = 0; compI < dim; ++compI) { 297 for (compJ = 0; compJ < dim; ++compJ) { 298 for (d = 0; d < dim; ++d) { 299 for (e = 0; e < dim; e++) { 300 if (d == e && d == compI && d == compJ) { 301 C[((compI*dim+compJ)*dim+d)*dim+e] = 1.0; 302 } else if ((d == compJ && e == compI) || (d == e && compI == compJ)) { 303 C[((compI*dim+compJ)*dim+d)*dim+e] = 0.5; 304 } else { 305 C[((compI*dim+compJ)*dim+d)*dim+e] = 0.0; 306 } 307 } 308 } 309 } 310 } 311 } 312 313 static PetscErrorCode SetupSection(DM dm, AppCtx *user) 314 { 315 PetscErrorCode ierr; 316 317 PetscFunctionBeginUser; 318 if (user->constraints) { 319 /* test local constraints */ 320 DM coordDM; 321 PetscInt fStart, fEnd, f, vStart, vEnd, v; 322 PetscInt edgesx = 2, vertsx; 323 PetscInt edgesy = 2, vertsy; 324 PetscMPIInt size; 325 PetscInt numConst; 326 PetscSection aSec; 327 PetscInt *anchors; 328 PetscInt offset; 329 IS aIS; 330 MPI_Comm comm = PetscObjectComm((PetscObject)dm); 331 332 ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 333 if (size > 1) SETERRQ(comm,PETSC_ERR_SUP,"Local constraint test can only be performed in serial"); 334 335 /* we are going to test constraints by using them to enforce periodicity 336 * in one direction, and comparing to the existing method of enforcing 337 * periodicity */ 338 339 /* first create the coordinate section so that it does not clone the 340 * constraints */ 341 ierr = DMGetCoordinateDM(dm,&coordDM);CHKERRQ(ierr); 342 343 /* create the constrained-to-anchor section */ 344 ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr); 345 ierr = DMPlexGetDepthStratum(dm,1,&fStart,&fEnd);CHKERRQ(ierr); 346 ierr = PetscSectionCreate(PETSC_COMM_SELF,&aSec);CHKERRQ(ierr); 347 ierr = PetscSectionSetChart(aSec,PetscMin(fStart,vStart),PetscMax(fEnd,vEnd));CHKERRQ(ierr); 348 349 /* define the constraints */ 350 ierr = PetscOptionsGetInt(NULL,NULL,"-da_grid_x",&edgesx,NULL);CHKERRQ(ierr); 351 ierr = PetscOptionsGetInt(NULL,NULL,"-da_grid_y",&edgesy,NULL);CHKERRQ(ierr); 352 vertsx = edgesx + 1; 353 vertsy = edgesy + 1; 354 numConst = vertsy + edgesy; 355 ierr = PetscMalloc1(numConst,&anchors);CHKERRQ(ierr); 356 offset = 0; 357 for (v = vStart + edgesx; v < vEnd; v+= vertsx) { 358 ierr = PetscSectionSetDof(aSec,v,1);CHKERRQ(ierr); 359 anchors[offset++] = v - edgesx; 360 } 361 for (f = fStart + edgesx * vertsy + edgesx * edgesy; f < fEnd; f++) { 362 ierr = PetscSectionSetDof(aSec,f,1);CHKERRQ(ierr); 363 anchors[offset++] = f - edgesx * edgesy; 364 } 365 ierr = PetscSectionSetUp(aSec);CHKERRQ(ierr); 366 ierr = ISCreateGeneral(PETSC_COMM_SELF,numConst,anchors,PETSC_OWN_POINTER,&aIS);CHKERRQ(ierr); 367 368 ierr = DMPlexSetAnchors(dm,aSec,aIS);CHKERRQ(ierr); 369 ierr = PetscSectionDestroy(&aSec);CHKERRQ(ierr); 370 ierr = ISDestroy(&aIS);CHKERRQ(ierr); 371 } 372 ierr = DMSetNumFields(dm, 1);CHKERRQ(ierr); 373 ierr = DMSetField(dm, 0, NULL, (PetscObject) user->fe);CHKERRQ(ierr); 374 ierr = DMCreateDS(dm);CHKERRQ(ierr); 375 if (user->constraints) { 376 /* test getting local constraint matrix that matches section */ 377 PetscSection aSec; 378 IS aIS; 379 380 ierr = DMPlexGetAnchors(dm,&aSec,&aIS);CHKERRQ(ierr); 381 if (aSec) { 382 PetscDS ds; 383 PetscSection cSec, section; 384 PetscInt cStart, cEnd, c, numComp; 385 Mat cMat, mass; 386 Vec local; 387 const PetscInt *anchors; 388 389 ierr = DMGetLocalSection(dm,§ion);CHKERRQ(ierr); 390 /* this creates the matrix and preallocates the matrix structure: we 391 * just have to fill in the values */ 392 ierr = DMGetDefaultConstraints(dm,&cSec,&cMat);CHKERRQ(ierr); 393 ierr = PetscSectionGetChart(cSec,&cStart,&cEnd);CHKERRQ(ierr); 394 ierr = ISGetIndices(aIS,&anchors);CHKERRQ(ierr); 395 ierr = PetscFEGetNumComponents(user->fe, &numComp);CHKERRQ(ierr); 396 for (c = cStart; c < cEnd; c++) { 397 PetscInt cDof; 398 399 /* is this point constrained? (does it have an anchor?) */ 400 ierr = PetscSectionGetDof(aSec,c,&cDof);CHKERRQ(ierr); 401 if (cDof) { 402 PetscInt cOff, a, aDof, aOff, j; 403 if (cDof != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Found %d anchor points: should be just one",cDof); 404 405 /* find the anchor point */ 406 ierr = PetscSectionGetOffset(aSec,c,&cOff);CHKERRQ(ierr); 407 a = anchors[cOff]; 408 409 /* find the constrained dofs (row in constraint matrix) */ 410 ierr = PetscSectionGetDof(cSec,c,&cDof);CHKERRQ(ierr); 411 ierr = PetscSectionGetOffset(cSec,c,&cOff);CHKERRQ(ierr); 412 413 /* find the anchor dofs (column in constraint matrix) */ 414 ierr = PetscSectionGetDof(section,a,&aDof);CHKERRQ(ierr); 415 ierr = PetscSectionGetOffset(section,a,&aOff);CHKERRQ(ierr); 416 417 if (cDof != aDof) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Point and anchor have different number of dofs: %d, %d\n",cDof,aDof); 418 if (cDof % numComp) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Point dofs not divisible by field components: %d, %d\n",cDof,numComp); 419 420 /* put in a simple equality constraint */ 421 for (j = 0; j < cDof; j++) { 422 ierr = MatSetValue(cMat,cOff+j,aOff+j,1.,INSERT_VALUES);CHKERRQ(ierr); 423 } 424 } 425 } 426 ierr = MatAssemblyBegin(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 427 ierr = MatAssemblyEnd(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 428 ierr = ISRestoreIndices(aIS,&anchors);CHKERRQ(ierr); 429 430 /* Now that we have constructed the constraint matrix, any FE matrix 431 * that we construct will apply the constraints during construction */ 432 433 ierr = DMCreateMatrix(dm,&mass);CHKERRQ(ierr); 434 /* get a dummy local variable to serve as the solution */ 435 ierr = DMGetLocalVector(dm,&local);CHKERRQ(ierr); 436 ierr = DMGetDS(dm,&ds);CHKERRQ(ierr); 437 /* set the jacobian to be the mass matrix */ 438 ierr = PetscDSSetJacobian(ds, 0, 0, simple_mass, NULL, NULL, NULL);CHKERRQ(ierr); 439 /* build the mass matrix */ 440 ierr = DMPlexSNESComputeJacobianFEM(dm,local,mass,mass,NULL);CHKERRQ(ierr); 441 ierr = MatView(mass,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 442 ierr = MatDestroy(&mass);CHKERRQ(ierr); 443 ierr = DMRestoreLocalVector(dm,&local);CHKERRQ(ierr); 444 } 445 } 446 PetscFunctionReturn(0); 447 } 448 449 static PetscErrorCode TestFEJacobian(DM dm, AppCtx *user) 450 { 451 PetscErrorCode ierr; 452 453 PetscFunctionBeginUser; 454 if (!user->useDA) { 455 Vec local; 456 const Vec *vecs; 457 Mat E; 458 MatNullSpace sp; 459 PetscBool isNullSpace, hasConst; 460 PetscInt dim, n, i; 461 Vec res = NULL, localX, localRes; 462 PetscDS ds; 463 464 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 465 if (user->numComponents != dim) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "The number of components %d must be equal to the dimension %d for this test", user->numComponents, dim); 466 ierr = DMGetDS(dm,&ds);CHKERRQ(ierr); 467 ierr = PetscDSSetJacobian(ds,0,0,NULL,NULL,NULL,symmetric_gradient_inner_product);CHKERRQ(ierr); 468 ierr = DMCreateMatrix(dm,&E);CHKERRQ(ierr); 469 ierr = DMGetLocalVector(dm,&local);CHKERRQ(ierr); 470 ierr = DMPlexSNESComputeJacobianFEM(dm,local,E,E,NULL);CHKERRQ(ierr); 471 ierr = DMPlexCreateRigidBody(dm,0,&sp);CHKERRQ(ierr); 472 ierr = MatNullSpaceGetVecs(sp,&hasConst,&n,&vecs);CHKERRQ(ierr); 473 if (n) {ierr = VecDuplicate(vecs[0],&res);CHKERRQ(ierr);} 474 ierr = DMCreateLocalVector(dm,&localX);CHKERRQ(ierr); 475 ierr = DMCreateLocalVector(dm,&localRes);CHKERRQ(ierr); 476 for (i = 0; i < n; i++) { /* also test via matrix-free Jacobian application */ 477 PetscReal resNorm; 478 479 ierr = VecSet(localRes,0.);CHKERRQ(ierr); 480 ierr = VecSet(localX,0.);CHKERRQ(ierr); 481 ierr = VecSet(local,0.);CHKERRQ(ierr); 482 ierr = VecSet(res,0.);CHKERRQ(ierr); 483 ierr = DMGlobalToLocalBegin(dm,vecs[i],INSERT_VALUES,localX);CHKERRQ(ierr); 484 ierr = DMGlobalToLocalEnd(dm,vecs[i],INSERT_VALUES,localX);CHKERRQ(ierr); 485 ierr = DMSNESComputeJacobianAction(dm,local,localX,localRes,NULL);CHKERRQ(ierr); 486 ierr = DMLocalToGlobalBegin(dm,localRes,ADD_VALUES,res);CHKERRQ(ierr); 487 ierr = DMLocalToGlobalEnd(dm,localRes,ADD_VALUES,res);CHKERRQ(ierr); 488 ierr = VecNorm(res,NORM_2,&resNorm);CHKERRQ(ierr); 489 if (resNorm > PETSC_SMALL) { 490 ierr = PetscPrintf(PetscObjectComm((PetscObject)dm),"Symmetric gradient action null space vector %D residual: %E\n",i,resNorm);CHKERRQ(ierr); 491 } 492 } 493 ierr = VecDestroy(&localRes);CHKERRQ(ierr); 494 ierr = VecDestroy(&localX);CHKERRQ(ierr); 495 ierr = VecDestroy(&res);CHKERRQ(ierr); 496 ierr = MatNullSpaceTest(sp,E,&isNullSpace);CHKERRQ(ierr); 497 if (isNullSpace) { 498 ierr = PetscPrintf(PetscObjectComm((PetscObject)dm),"Symmetric gradient null space: PASS\n");CHKERRQ(ierr); 499 } else { 500 ierr = PetscPrintf(PetscObjectComm((PetscObject)dm),"Symmetric gradient null space: FAIL\n");CHKERRQ(ierr); 501 } 502 ierr = MatNullSpaceDestroy(&sp);CHKERRQ(ierr); 503 ierr = MatDestroy(&E);CHKERRQ(ierr); 504 ierr = DMRestoreLocalVector(dm,&local);CHKERRQ(ierr); 505 } 506 PetscFunctionReturn(0); 507 } 508 509 static PetscErrorCode TestInjector(DM dm, AppCtx *user) 510 { 511 DM refTree; 512 PetscMPIInt rank; 513 PetscErrorCode ierr; 514 515 PetscFunctionBegin; 516 ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr); 517 if (refTree) { 518 Mat inj; 519 520 ierr = DMPlexComputeInjectorReferenceTree(refTree,&inj);CHKERRQ(ierr); 521 ierr = PetscObjectSetName((PetscObject)inj,"Reference Tree Injector");CHKERRQ(ierr); 522 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 523 if (!rank) { 524 ierr = MatView(inj,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 525 } 526 ierr = MatDestroy(&inj);CHKERRQ(ierr); 527 } 528 PetscFunctionReturn(0); 529 } 530 531 static PetscErrorCode TestFVGrad(DM dm, AppCtx *user) 532 { 533 MPI_Comm comm; 534 DM dmRedist, dmfv, dmgrad, dmCell, refTree; 535 PetscFV fv; 536 PetscInt dim, nvecs, v, cStart, cEnd, cEndInterior; 537 PetscMPIInt size; 538 Vec cellgeom, grad, locGrad; 539 const PetscScalar *cgeom; 540 PetscReal allVecMaxDiff = 0., fvTol = 100. * PETSC_MACHINE_EPSILON; 541 PetscErrorCode ierr; 542 543 PetscFunctionBeginUser; 544 comm = PetscObjectComm((PetscObject)dm); 545 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 546 /* duplicate DM, give dup. a FV discretization */ 547 ierr = DMSetBasicAdjacency(dm,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr); 548 ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 549 dmRedist = NULL; 550 if (size > 1) { 551 ierr = DMPlexDistributeOverlap(dm,1,NULL,&dmRedist);CHKERRQ(ierr); 552 } 553 if (!dmRedist) { 554 dmRedist = dm; 555 ierr = PetscObjectReference((PetscObject)dmRedist);CHKERRQ(ierr); 556 } 557 ierr = PetscFVCreate(comm,&fv);CHKERRQ(ierr); 558 ierr = PetscFVSetType(fv,PETSCFVLEASTSQUARES);CHKERRQ(ierr); 559 ierr = PetscFVSetNumComponents(fv,user->numComponents);CHKERRQ(ierr); 560 ierr = PetscFVSetSpatialDimension(fv,dim);CHKERRQ(ierr); 561 ierr = PetscFVSetFromOptions(fv);CHKERRQ(ierr); 562 ierr = PetscFVSetUp(fv);CHKERRQ(ierr); 563 ierr = DMPlexConstructGhostCells(dmRedist,NULL,NULL,&dmfv);CHKERRQ(ierr); 564 ierr = DMDestroy(&dmRedist);CHKERRQ(ierr); 565 ierr = DMSetNumFields(dmfv,1);CHKERRQ(ierr); 566 ierr = DMSetField(dmfv, 0, NULL, (PetscObject) fv);CHKERRQ(ierr); 567 ierr = DMCreateDS(dmfv);CHKERRQ(ierr); 568 ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr); 569 if (refTree) {ierr = DMCopyDisc(dmfv,refTree);CHKERRQ(ierr);} 570 ierr = DMPlexGetGradientDM(dmfv, fv, &dmgrad);CHKERRQ(ierr); 571 ierr = DMPlexGetHeightStratum(dmfv,0,&cStart,&cEnd);CHKERRQ(ierr); 572 nvecs = dim * (dim+1) / 2; 573 ierr = DMPlexGetGeometryFVM(dmfv,NULL,&cellgeom,NULL);CHKERRQ(ierr); 574 ierr = VecGetDM(cellgeom,&dmCell);CHKERRQ(ierr); 575 ierr = VecGetArrayRead(cellgeom,&cgeom);CHKERRQ(ierr); 576 ierr = DMGetGlobalVector(dmgrad,&grad);CHKERRQ(ierr); 577 ierr = DMGetLocalVector(dmgrad,&locGrad);CHKERRQ(ierr); 578 ierr = DMPlexGetGhostCellStratum(dmgrad,&cEndInterior,NULL);CHKERRQ(ierr); 579 cEndInterior = (cEndInterior < 0) ? cEnd: cEndInterior; 580 for (v = 0; v < nvecs; v++) { 581 Vec locX; 582 PetscInt c; 583 PetscScalar trueGrad[3][3] = {{0.}}; 584 const PetscScalar *gradArray; 585 PetscReal maxDiff, maxDiffGlob; 586 587 ierr = DMGetLocalVector(dmfv,&locX);CHKERRQ(ierr); 588 /* get the local projection of the rigid body mode */ 589 for (c = cStart; c < cEnd; c++) { 590 PetscFVCellGeom *cg; 591 PetscScalar cx[3] = {0.,0.,0.}; 592 593 ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 594 if (v < dim) { 595 cx[v] = 1.; 596 } else { 597 PetscInt w = v - dim; 598 599 cx[(w + 1) % dim] = cg->centroid[(w + 2) % dim]; 600 cx[(w + 2) % dim] = -cg->centroid[(w + 1) % dim]; 601 } 602 ierr = DMPlexVecSetClosure(dmfv,NULL,locX,c,cx,INSERT_ALL_VALUES);CHKERRQ(ierr); 603 } 604 /* TODO: this isn't in any header */ 605 ierr = DMPlexReconstructGradientsFVM(dmfv,locX,grad);CHKERRQ(ierr); 606 ierr = DMGlobalToLocalBegin(dmgrad,grad,INSERT_VALUES,locGrad);CHKERRQ(ierr); 607 ierr = DMGlobalToLocalEnd(dmgrad,grad,INSERT_VALUES,locGrad);CHKERRQ(ierr); 608 ierr = VecGetArrayRead(locGrad,&gradArray);CHKERRQ(ierr); 609 /* compare computed gradient to exact gradient */ 610 if (v >= dim) { 611 PetscInt w = v - dim; 612 613 trueGrad[(w + 1) % dim][(w + 2) % dim] = 1.; 614 trueGrad[(w + 2) % dim][(w + 1) % dim] = -1.; 615 } 616 maxDiff = 0.; 617 for (c = cStart; c < cEndInterior; c++) { 618 PetscScalar *compGrad; 619 PetscInt i, j, k; 620 PetscReal FrobDiff = 0.; 621 622 ierr = DMPlexPointLocalRead(dmgrad, c, gradArray, &compGrad);CHKERRQ(ierr); 623 624 for (i = 0, k = 0; i < dim; i++) { 625 for (j = 0; j < dim; j++, k++) { 626 PetscScalar diff = compGrad[k] - trueGrad[i][j]; 627 FrobDiff += PetscRealPart(diff * PetscConj(diff)); 628 } 629 } 630 FrobDiff = PetscSqrtReal(FrobDiff); 631 maxDiff = PetscMax(maxDiff,FrobDiff); 632 } 633 ierr = MPI_Allreduce(&maxDiff,&maxDiffGlob,1,MPIU_REAL,MPIU_MAX,comm);CHKERRMPI(ierr); 634 allVecMaxDiff = PetscMax(allVecMaxDiff,maxDiffGlob); 635 ierr = VecRestoreArrayRead(locGrad,&gradArray);CHKERRQ(ierr); 636 ierr = DMRestoreLocalVector(dmfv,&locX);CHKERRQ(ierr); 637 } 638 if (allVecMaxDiff < fvTol) { 639 ierr = PetscPrintf(PetscObjectComm((PetscObject)dm),"Finite volume gradient reconstruction: PASS\n");CHKERRQ(ierr); 640 } else { 641 ierr = PetscPrintf(PetscObjectComm((PetscObject)dm),"Finite volume gradient reconstruction: FAIL at tolerance %g with max difference %g\n",fvTol,allVecMaxDiff);CHKERRQ(ierr); 642 } 643 ierr = DMRestoreLocalVector(dmgrad,&locGrad);CHKERRQ(ierr); 644 ierr = DMRestoreGlobalVector(dmgrad,&grad);CHKERRQ(ierr); 645 ierr = VecRestoreArrayRead(cellgeom,&cgeom);CHKERRQ(ierr); 646 ierr = DMDestroy(&dmfv);CHKERRQ(ierr); 647 ierr = PetscFVDestroy(&fv);CHKERRQ(ierr); 648 PetscFunctionReturn(0); 649 } 650 651 static PetscErrorCode ComputeError(DM dm, PetscErrorCode (**exactFuncs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), 652 PetscErrorCode (**exactFuncDers)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), 653 void **exactCtxs, PetscReal *error, PetscReal *errorDer, AppCtx *user) 654 { 655 Vec u; 656 PetscReal n[3] = {1.0, 1.0, 1.0}; 657 PetscErrorCode ierr; 658 659 PetscFunctionBeginUser; 660 ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr); 661 /* Project function into FE function space */ 662 ierr = DMProjectFunction(dm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, u);CHKERRQ(ierr); 663 ierr = VecViewFromOptions(u, NULL, "-projection_view");CHKERRQ(ierr); 664 /* Compare approximation to exact in L_2 */ 665 ierr = DMComputeL2Diff(dm, 0.0, exactFuncs, exactCtxs, u, error);CHKERRQ(ierr); 666 ierr = DMComputeL2GradientDiff(dm, 0.0, exactFuncDers, exactCtxs, u, n, errorDer);CHKERRQ(ierr); 667 ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr); 668 PetscFunctionReturn(0); 669 } 670 671 static PetscErrorCode CheckFunctions(DM dm, PetscInt order, AppCtx *user) 672 { 673 PetscErrorCode (*exactFuncs[1]) (PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 674 PetscErrorCode (*exactFuncDers[1]) (PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx); 675 void *exactCtxs[3]; 676 MPI_Comm comm; 677 PetscReal error, errorDer, tol = PETSC_SMALL; 678 PetscErrorCode ierr; 679 680 PetscFunctionBeginUser; 681 exactCtxs[0] = user; 682 exactCtxs[1] = user; 683 exactCtxs[2] = user; 684 ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 685 /* Setup functions to approximate */ 686 switch (order) { 687 case 0: 688 exactFuncs[0] = constant; 689 exactFuncDers[0] = constantDer; 690 break; 691 case 1: 692 exactFuncs[0] = linear; 693 exactFuncDers[0] = linearDer; 694 break; 695 case 2: 696 exactFuncs[0] = quadratic; 697 exactFuncDers[0] = quadraticDer; 698 break; 699 case 3: 700 exactFuncs[0] = cubic; 701 exactFuncDers[0] = cubicDer; 702 break; 703 default: 704 SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Could not determine functions to test for order %d", order); 705 } 706 ierr = ComputeError(dm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user);CHKERRQ(ierr); 707 /* Report result */ 708 if (error > tol) {ierr = PetscPrintf(comm, "Function tests FAIL for order %D at tolerance %g error %g\n", order, (double)tol,(double) error);CHKERRQ(ierr);} 709 else {ierr = PetscPrintf(comm, "Function tests pass for order %D at tolerance %g\n", order, (double)tol);CHKERRQ(ierr);} 710 if (errorDer > tol) {ierr = PetscPrintf(comm, "Function tests FAIL for order %D derivatives at tolerance %g error %g\n", order, (double)tol, (double)errorDer);CHKERRQ(ierr);} 711 else {ierr = PetscPrintf(comm, "Function tests pass for order %D derivatives at tolerance %g\n", order, (double)tol);CHKERRQ(ierr);} 712 PetscFunctionReturn(0); 713 } 714 715 static PetscErrorCode CheckInterpolation(DM dm, PetscBool checkRestrict, PetscInt order, AppCtx *user) 716 { 717 PetscErrorCode (*exactFuncs[1]) (PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx); 718 PetscErrorCode (*exactFuncDers[1]) (PetscInt, PetscReal, const PetscReal x[], const PetscReal n[], PetscInt, PetscScalar *u, void *ctx); 719 PetscReal n[3] = {1.0, 1.0, 1.0}; 720 void *exactCtxs[3]; 721 DM rdm, idm, fdm; 722 Mat Interp; 723 Vec iu, fu, scaling; 724 MPI_Comm comm; 725 PetscInt dim; 726 PetscReal error, errorDer, tol = PETSC_SMALL; 727 PetscErrorCode ierr; 728 729 PetscFunctionBeginUser; 730 exactCtxs[0] = user; 731 exactCtxs[1] = user; 732 exactCtxs[2] = user; 733 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 734 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 735 ierr = DMRefine(dm, comm, &rdm);CHKERRQ(ierr); 736 ierr = DMSetCoarseDM(rdm, dm);CHKERRQ(ierr); 737 ierr = DMPlexSetRegularRefinement(rdm, user->convRefine);CHKERRQ(ierr); 738 if (user->tree) { 739 DM refTree; 740 ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr); 741 ierr = DMPlexSetReferenceTree(rdm,refTree);CHKERRQ(ierr); 742 } 743 if (user->useDA) {ierr = DMDASetVertexCoordinates(rdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0);CHKERRQ(ierr);} 744 ierr = SetupSection(rdm, user);CHKERRQ(ierr); 745 /* Setup functions to approximate */ 746 switch (order) { 747 case 0: 748 exactFuncs[0] = constant; 749 exactFuncDers[0] = constantDer; 750 break; 751 case 1: 752 exactFuncs[0] = linear; 753 exactFuncDers[0] = linearDer; 754 break; 755 case 2: 756 exactFuncs[0] = quadratic; 757 exactFuncDers[0] = quadraticDer; 758 break; 759 case 3: 760 exactFuncs[0] = cubic; 761 exactFuncDers[0] = cubicDer; 762 break; 763 default: 764 SETERRQ2(comm, PETSC_ERR_ARG_OUTOFRANGE, "Could not determine functions to test for dimension %D order %D", dim, order); 765 } 766 idm = checkRestrict ? rdm : dm; 767 fdm = checkRestrict ? dm : rdm; 768 ierr = DMGetGlobalVector(idm, &iu);CHKERRQ(ierr); 769 ierr = DMGetGlobalVector(fdm, &fu);CHKERRQ(ierr); 770 ierr = DMSetApplicationContext(dm, user);CHKERRQ(ierr); 771 ierr = DMSetApplicationContext(rdm, user);CHKERRQ(ierr); 772 ierr = DMCreateInterpolation(dm, rdm, &Interp, &scaling);CHKERRQ(ierr); 773 /* Project function into initial FE function space */ 774 ierr = DMProjectFunction(idm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, iu);CHKERRQ(ierr); 775 /* Interpolate function into final FE function space */ 776 if (checkRestrict) {ierr = MatRestrict(Interp, iu, fu);CHKERRQ(ierr);ierr = VecPointwiseMult(fu, scaling, fu);CHKERRQ(ierr);} 777 else {ierr = MatInterpolate(Interp, iu, fu);CHKERRQ(ierr);} 778 /* Compare approximation to exact in L_2 */ 779 ierr = DMComputeL2Diff(fdm, 0.0, exactFuncs, exactCtxs, fu, &error);CHKERRQ(ierr); 780 ierr = DMComputeL2GradientDiff(fdm, 0.0, exactFuncDers, exactCtxs, fu, n, &errorDer);CHKERRQ(ierr); 781 /* Report result */ 782 if (error > tol) {ierr = PetscPrintf(comm, "Interpolation tests FAIL for order %D at tolerance %g error %g\n", order, (double)tol, (double)error);CHKERRQ(ierr);} 783 else {ierr = PetscPrintf(comm, "Interpolation tests pass for order %D at tolerance %g\n", order, (double)tol);CHKERRQ(ierr);} 784 if (errorDer > tol) {ierr = PetscPrintf(comm, "Interpolation tests FAIL for order %D derivatives at tolerance %g error %g\n", order, (double)tol, (double)errorDer);CHKERRQ(ierr);} 785 else {ierr = PetscPrintf(comm, "Interpolation tests pass for order %D derivatives at tolerance %g\n", order, (double)tol);CHKERRQ(ierr);} 786 ierr = DMRestoreGlobalVector(idm, &iu);CHKERRQ(ierr); 787 ierr = DMRestoreGlobalVector(fdm, &fu);CHKERRQ(ierr); 788 ierr = MatDestroy(&Interp);CHKERRQ(ierr); 789 ierr = VecDestroy(&scaling);CHKERRQ(ierr); 790 ierr = DMDestroy(&rdm);CHKERRQ(ierr); 791 PetscFunctionReturn(0); 792 } 793 794 static PetscErrorCode CheckConvergence(DM dm, PetscInt Nr, AppCtx *user) 795 { 796 DM odm = dm, rdm = NULL, cdm = NULL; 797 PetscErrorCode (*exactFuncs[1]) (PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) = {trig}; 798 PetscErrorCode (*exactFuncDers[1]) (PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) = {trigDer}; 799 void *exactCtxs[3]; 800 PetscInt r, c, cStart, cEnd; 801 PetscReal errorOld, errorDerOld, error, errorDer, rel, len, lenOld; 802 double p; 803 PetscErrorCode ierr; 804 805 PetscFunctionBeginUser; 806 if (!user->convergence) PetscFunctionReturn(0); 807 exactCtxs[0] = user; 808 exactCtxs[1] = user; 809 exactCtxs[2] = user; 810 ierr = PetscObjectReference((PetscObject) odm);CHKERRQ(ierr); 811 if (!user->convRefine) { 812 for (r = 0; r < Nr; ++r) { 813 ierr = DMRefine(odm, PetscObjectComm((PetscObject) dm), &rdm);CHKERRQ(ierr); 814 ierr = DMDestroy(&odm);CHKERRQ(ierr); 815 odm = rdm; 816 } 817 ierr = SetupSection(odm, user);CHKERRQ(ierr); 818 } 819 ierr = ComputeError(odm, exactFuncs, exactFuncDers, exactCtxs, &errorOld, &errorDerOld, user);CHKERRQ(ierr); 820 if (user->convRefine) { 821 for (r = 0; r < Nr; ++r) { 822 ierr = DMRefine(odm, PetscObjectComm((PetscObject) dm), &rdm);CHKERRQ(ierr); 823 if (user->useDA) {ierr = DMDASetVertexCoordinates(rdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0);CHKERRQ(ierr);} 824 ierr = SetupSection(rdm, user);CHKERRQ(ierr); 825 ierr = ComputeError(rdm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user);CHKERRQ(ierr); 826 p = PetscLog2Real(errorOld/error); 827 ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Function convergence rate at refinement %D: %.2f\n", r, (double)p);CHKERRQ(ierr); 828 p = PetscLog2Real(errorDerOld/errorDer); 829 ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Derivative convergence rate at refinement %D: %.2f\n", r, (double)p);CHKERRQ(ierr); 830 ierr = DMDestroy(&odm);CHKERRQ(ierr); 831 odm = rdm; 832 errorOld = error; 833 errorDerOld = errorDer; 834 } 835 } else { 836 /* ierr = ComputeLongestEdge(dm, &lenOld);CHKERRQ(ierr); */ 837 ierr = DMPlexGetHeightStratum(odm, 0, &cStart, &cEnd);CHKERRQ(ierr); 838 lenOld = cEnd - cStart; 839 for (c = 0; c < Nr; ++c) { 840 ierr = DMCoarsen(odm, PetscObjectComm((PetscObject) dm), &cdm);CHKERRQ(ierr); 841 if (user->useDA) {ierr = DMDASetVertexCoordinates(cdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0);CHKERRQ(ierr);} 842 ierr = SetupSection(cdm, user);CHKERRQ(ierr); 843 ierr = ComputeError(cdm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user);CHKERRQ(ierr); 844 /* ierr = ComputeLongestEdge(cdm, &len);CHKERRQ(ierr); */ 845 ierr = DMPlexGetHeightStratum(cdm, 0, &cStart, &cEnd);CHKERRQ(ierr); 846 len = cEnd - cStart; 847 rel = error/errorOld; 848 p = PetscLogReal(rel) / PetscLogReal(lenOld / len); 849 ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Function convergence rate at coarsening %D: %.2f\n", c, (double)p);CHKERRQ(ierr); 850 rel = errorDer/errorDerOld; 851 p = PetscLogReal(rel) / PetscLogReal(lenOld / len); 852 ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Derivative convergence rate at coarsening %D: %.2f\n", c, (double)p);CHKERRQ(ierr); 853 ierr = DMDestroy(&odm);CHKERRQ(ierr); 854 odm = cdm; 855 errorOld = error; 856 errorDerOld = errorDer; 857 lenOld = len; 858 } 859 } 860 ierr = DMDestroy(&odm);CHKERRQ(ierr); 861 PetscFunctionReturn(0); 862 } 863 864 int main(int argc, char **argv) 865 { 866 DM dm; 867 AppCtx user; /* user-defined work context */ 868 PetscInt dim = 2; 869 PetscBool simplex = PETSC_FALSE; 870 PetscErrorCode ierr; 871 872 ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 873 ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 874 ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 875 if (!user.useDA) { 876 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 877 ierr = DMPlexIsSimplex(dm, &simplex);CHKERRQ(ierr); 878 } 879 user.numComponents = user.numComponents < 0 ? dim : user.numComponents; 880 ierr = PetscFECreateDefault(PETSC_COMM_WORLD, dim, user.numComponents, simplex, NULL, user.qorder, &user.fe);CHKERRQ(ierr); 881 ierr = SetupSection(dm, &user);CHKERRQ(ierr); 882 if (user.testFEjacobian) {ierr = TestFEJacobian(dm, &user);CHKERRQ(ierr);} 883 if (user.testFVgrad) {ierr = TestFVGrad(dm, &user);CHKERRQ(ierr);} 884 if (user.testInjector) {ierr = TestInjector(dm, &user);CHKERRQ(ierr);} 885 ierr = CheckFunctions(dm, user.porder, &user);CHKERRQ(ierr); 886 { 887 PetscDualSpace dsp; 888 PetscInt k; 889 890 ierr = PetscFEGetDualSpace(user.fe, &dsp);CHKERRQ(ierr); 891 ierr = PetscDualSpaceGetDeRahm(dsp, &k);CHKERRQ(ierr); 892 if (dim == 2 && simplex == PETSC_TRUE && user.tree == PETSC_FALSE && k == 0) { 893 ierr = CheckInterpolation(dm, PETSC_FALSE, user.porder, &user);CHKERRQ(ierr); 894 ierr = CheckInterpolation(dm, PETSC_TRUE, user.porder, &user);CHKERRQ(ierr); 895 } 896 } 897 ierr = CheckConvergence(dm, 3, &user);CHKERRQ(ierr); 898 ierr = PetscFEDestroy(&user.fe);CHKERRQ(ierr); 899 ierr = DMDestroy(&dm);CHKERRQ(ierr); 900 ierr = PetscFinalize(); 901 return ierr; 902 } 903 904 /*TEST 905 906 test: 907 suffix: 1 908 requires: triangle 909 910 # 2D P_1 on a triangle 911 test: 912 suffix: p1_2d_0 913 requires: triangle 914 args: -petscspace_degree 1 -qorder 1 -convergence 915 test: 916 suffix: p1_2d_1 917 requires: triangle 918 args: -petscspace_degree 1 -qorder 1 -porder 1 919 test: 920 suffix: p1_2d_2 921 requires: triangle 922 args: -petscspace_degree 1 -qorder 1 -porder 2 923 test: 924 suffix: p1_2d_3 925 requires: triangle pragmatic 926 args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0 927 filter: grep -v DEBUG 928 test: 929 suffix: p1_2d_4 930 requires: triangle pragmatic 931 args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0 932 test: 933 suffix: p1_2d_5 934 requires: triangle pragmatic 935 args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0 936 937 # 3D P_1 on a tetrahedron 938 test: 939 suffix: p1_3d_0 940 requires: ctetgen 941 args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -convergence 942 test: 943 suffix: p1_3d_1 944 requires: ctetgen 945 args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -porder 1 946 test: 947 suffix: p1_3d_2 948 requires: ctetgen 949 args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -porder 2 950 test: 951 suffix: p1_3d_3 952 requires: ctetgen pragmatic 953 args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0 954 filter: grep -v DEBUG 955 test: 956 suffix: p1_3d_4 957 requires: ctetgen pragmatic 958 args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0 959 test: 960 suffix: p1_3d_5 961 requires: ctetgen pragmatic 962 args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0 963 964 # 2D P_2 on a triangle 965 test: 966 suffix: p2_2d_0 967 requires: triangle 968 args: -petscspace_degree 2 -qorder 2 -convergence 969 test: 970 suffix: p2_2d_1 971 requires: triangle 972 args: -petscspace_degree 2 -qorder 2 -porder 1 973 test: 974 suffix: p2_2d_2 975 requires: triangle 976 args: -petscspace_degree 2 -qorder 2 -porder 2 977 test: 978 suffix: p2_2d_3 979 requires: triangle pragmatic 980 args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -convergence -conv_refine 0 981 filter: grep -v DEBUG 982 test: 983 suffix: p2_2d_4 984 requires: triangle pragmatic 985 args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 1 -conv_refine 0 986 test: 987 suffix: p2_2d_5 988 requires: triangle pragmatic 989 args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 2 -conv_refine 0 990 991 # 3D P_2 on a tetrahedron 992 test: 993 suffix: p2_3d_0 994 requires: ctetgen 995 args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -convergence 996 test: 997 suffix: p2_3d_1 998 requires: ctetgen 999 args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -porder 1 1000 test: 1001 suffix: p2_3d_2 1002 requires: ctetgen 1003 args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -porder 2 1004 test: 1005 suffix: p2_3d_3 1006 requires: ctetgen pragmatic 1007 args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -convergence -conv_refine 0 1008 filter: grep -v DEBUG 1009 test: 1010 suffix: p2_3d_4 1011 requires: ctetgen pragmatic 1012 args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 1 -conv_refine 0 1013 test: 1014 suffix: p2_3d_5 1015 requires: ctetgen pragmatic 1016 args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 2 -conv_refine 0 1017 1018 # 2D Q_1 on a quadrilaterial DA 1019 test: 1020 suffix: q1_2d_da_0 1021 requires: mpi_type_get_envelope broken 1022 args: -use_da 1 -petscspace_degree 1 -qorder 1 -convergence 1023 test: 1024 suffix: q1_2d_da_1 1025 requires: mpi_type_get_envelope broken 1026 args: -use_da 1 -petscspace_degree 1 -qorder 1 -porder 1 1027 test: 1028 suffix: q1_2d_da_2 1029 requires: mpi_type_get_envelope broken 1030 args: -use_da 1 -petscspace_degree 1 -qorder 1 -porder 2 1031 1032 # 2D Q_1 on a quadrilaterial Plex 1033 test: 1034 suffix: q1_2d_plex_0 1035 args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -convergence 1036 test: 1037 suffix: q1_2d_plex_1 1038 args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 1 1039 test: 1040 suffix: q1_2d_plex_2 1041 args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 2 1042 test: 1043 suffix: q1_2d_plex_3 1044 args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 1 -shear_coords 1045 test: 1046 suffix: q1_2d_plex_4 1047 args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 2 -shear_coords 1048 test: 1049 suffix: q1_2d_plex_5 1050 args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 0 -non_affine_coords -convergence 1051 test: 1052 suffix: q1_2d_plex_6 1053 args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 1 -non_affine_coords -convergence 1054 test: 1055 suffix: q1_2d_plex_7 1056 args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 2 -non_affine_coords -convergence 1057 1058 # 2D Q_2 on a quadrilaterial 1059 test: 1060 suffix: q2_2d_plex_0 1061 requires: mpi_type_get_envelope 1062 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -convergence 1063 test: 1064 suffix: q2_2d_plex_1 1065 requires: mpi_type_get_envelope 1066 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 1067 test: 1068 suffix: q2_2d_plex_2 1069 requires: mpi_type_get_envelope 1070 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 1071 test: 1072 suffix: q2_2d_plex_3 1073 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 -shear_coords 1074 test: 1075 suffix: q2_2d_plex_4 1076 requires: mpi_type_get_envelope 1077 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 -shear_coords 1078 test: 1079 suffix: q2_2d_plex_5 1080 requires: mpi_type_get_envelope 1081 args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 0 -non_affine_coords -convergence 1082 test: 1083 suffix: q2_2d_plex_6 1084 requires: mpi_type_get_envelope 1085 args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 1 -non_affine_coords -convergence 1086 test: 1087 suffix: q2_2d_plex_7 1088 requires: mpi_type_get_envelope 1089 args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 2 -non_affine_coords -convergence 1090 1091 # 2D P_3 on a triangle 1092 test: 1093 suffix: p3_2d_0 1094 requires: triangle !single 1095 args: -petscspace_degree 3 -qorder 3 -convergence 1096 test: 1097 suffix: p3_2d_1 1098 requires: triangle !single 1099 args: -petscspace_degree 3 -qorder 3 -porder 1 1100 test: 1101 suffix: p3_2d_2 1102 requires: triangle !single 1103 args: -petscspace_degree 3 -qorder 3 -porder 2 1104 test: 1105 suffix: p3_2d_3 1106 requires: triangle !single 1107 args: -petscspace_degree 3 -qorder 3 -porder 3 1108 test: 1109 suffix: p3_2d_4 1110 requires: triangle pragmatic 1111 args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -convergence -conv_refine 0 1112 filter: grep -v DEBUG 1113 test: 1114 suffix: p3_2d_5 1115 requires: triangle pragmatic 1116 args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -porder 1 -conv_refine 0 1117 test: 1118 suffix: p3_2d_6 1119 requires: triangle pragmatic 1120 args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -porder 3 -conv_refine 0 1121 1122 # 2D Q_3 on a quadrilaterial 1123 test: 1124 suffix: q3_2d_0 1125 requires: mpi_type_get_envelope !single 1126 args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -convergence 1127 test: 1128 suffix: q3_2d_1 1129 requires: mpi_type_get_envelope !single 1130 args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 1 1131 test: 1132 suffix: q3_2d_2 1133 requires: mpi_type_get_envelope !single 1134 args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 2 1135 test: 1136 suffix: q3_2d_3 1137 requires: mpi_type_get_envelope !single 1138 args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 3 1139 1140 # 2D P_1disc on a triangle/quadrilateral 1141 test: 1142 suffix: p1d_2d_0 1143 requires: triangle 1144 args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -convergence 1145 test: 1146 suffix: p1d_2d_1 1147 requires: triangle 1148 args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 1 1149 test: 1150 suffix: p1d_2d_2 1151 requires: triangle 1152 args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 2 1153 test: 1154 suffix: p1d_2d_3 1155 requires: triangle 1156 args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -convergence 1157 filter: sed -e "s/convergence rate at refinement 0: 2/convergence rate at refinement 0: 1.9/g" 1158 test: 1159 suffix: p1d_2d_4 1160 requires: triangle 1161 args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 1 1162 test: 1163 suffix: p1d_2d_5 1164 requires: triangle 1165 args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 2 1166 1167 # 2D BDM_1 on a triangle 1168 test: 1169 suffix: bdm1_2d_0 1170 requires: triangle 1171 args: -petscspace_degree 1 -petscdualspace_type bdm \ 1172 -num_comp 2 -qorder 1 -convergence 1173 test: 1174 suffix: bdm1_2d_1 1175 requires: triangle 1176 args: -petscspace_degree 1 -petscdualspace_type bdm \ 1177 -num_comp 2 -qorder 1 -porder 1 1178 test: 1179 suffix: bdm1_2d_2 1180 requires: triangle 1181 args: -petscspace_degree 1 -petscdualspace_type bdm \ 1182 -num_comp 2 -qorder 1 -porder 2 1183 1184 # 2D BDM_1 on a quadrilateral 1185 test: 1186 suffix: bdm1q_2d_0 1187 requires: triangle 1188 args: -petscspace_degree 1 -petscdualspace_type bdm \ 1189 -petscdualspace_lagrange_tensor 1 \ 1190 -dm_plex_simplex 0 -num_comp 2 -qorder 1 -convergence 1191 test: 1192 suffix: bdm1q_2d_1 1193 requires: triangle 1194 args: -petscspace_degree 1 -petscdualspace_type bdm \ 1195 -petscdualspace_lagrange_tensor 1 \ 1196 -dm_plex_simplex 0 -num_comp 2 -qorder 1 -porder 1 1197 test: 1198 suffix: bdm1q_2d_2 1199 requires: triangle 1200 args: -petscspace_degree 1 -petscdualspace_type bdm \ 1201 -petscdualspace_lagrange_tensor 1 \ 1202 -dm_plex_simplex 0 -num_comp 2 -qorder 1 -porder 2 1203 1204 # Test high order quadrature 1205 test: 1206 suffix: p1_quad_2 1207 requires: triangle 1208 args: -petscspace_degree 1 -qorder 2 -porder 1 1209 test: 1210 suffix: p1_quad_5 1211 requires: triangle 1212 args: -petscspace_degree 1 -qorder 5 -porder 1 1213 test: 1214 suffix: p2_quad_3 1215 requires: triangle 1216 args: -petscspace_degree 2 -qorder 3 -porder 2 1217 test: 1218 suffix: p2_quad_5 1219 requires: triangle 1220 args: -petscspace_degree 2 -qorder 5 -porder 2 1221 test: 1222 suffix: q1_quad_2 1223 requires: mpi_type_get_envelope 1224 args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 2 -porder 1 1225 test: 1226 suffix: q1_quad_5 1227 requires: mpi_type_get_envelope 1228 args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 5 -porder 1 1229 test: 1230 suffix: q2_quad_3 1231 requires: mpi_type_get_envelope 1232 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 3 -porder 1 1233 test: 1234 suffix: q2_quad_5 1235 requires: mpi_type_get_envelope 1236 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 5 -porder 1 1237 1238 # Nonconforming tests 1239 test: 1240 suffix: constraints 1241 args: -dm_coord_space 0 -dm_plex_simplex 0 -petscspace_type tensor -petscspace_degree 1 -qorder 0 -constraints 1242 test: 1243 suffix: nonconforming_tensor_2 1244 nsize: 4 1245 args: -dist_dm_distribute -test_fe_jacobian -test_injector -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_max_projection_height 1 -petscspace_type tensor -petscspace_degree 2 -qorder 2 -dm_view ascii::ASCII_INFO_DETAIL 1246 test: 1247 suffix: nonconforming_tensor_3 1248 nsize: 4 1249 args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_plex_max_projection_height 2 -petscspace_type tensor -petscspace_degree 1 -qorder 1 -dm_view ascii::ASCII_INFO_DETAIL 1250 test: 1251 suffix: nonconforming_tensor_2_fv 1252 nsize: 4 1253 args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -dm_plex_simplex 0 -num_comp 2 1254 test: 1255 suffix: nonconforming_tensor_3_fv 1256 nsize: 4 1257 args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -num_comp 3 1258 test: 1259 suffix: nonconforming_tensor_2_hi 1260 requires: !single 1261 nsize: 4 1262 args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_max_projection_height 1 -petscspace_type tensor -petscspace_degree 4 -qorder 4 1263 test: 1264 suffix: nonconforming_tensor_3_hi 1265 requires: !single skip 1266 nsize: 4 1267 args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_plex_max_projection_height 2 -petscspace_type tensor -petscspace_degree 4 -qorder 4 1268 test: 1269 suffix: nonconforming_simplex_2 1270 requires: triangle 1271 nsize: 4 1272 args: -dist_dm_distribute -test_fe_jacobian -test_injector -petscpartitioner_type simple -tree -dm_plex_max_projection_height 1 -petscspace_degree 2 -qorder 2 -dm_view ascii::ASCII_INFO_DETAIL 1273 test: 1274 suffix: nonconforming_simplex_2_hi 1275 requires: triangle !single 1276 nsize: 4 1277 args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_max_projection_height 1 -petscspace_degree 4 -qorder 4 1278 test: 1279 suffix: nonconforming_simplex_2_fv 1280 requires: triangle 1281 nsize: 4 1282 args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -num_comp 2 1283 test: 1284 suffix: nonconforming_simplex_3 1285 requires: ctetgen 1286 nsize: 4 1287 args: -dist_dm_distribute -test_fe_jacobian -test_injector -petscpartitioner_type simple -tree -dm_plex_dim 3 -dm_plex_max_projection_height 2 -petscspace_degree 2 -qorder 2 -dm_view ascii::ASCII_INFO_DETAIL 1288 test: 1289 suffix: nonconforming_simplex_3_hi 1290 requires: ctetgen skip 1291 nsize: 4 1292 args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_dim 3 -dm_plex_max_projection_height 2 -petscspace_degree 4 -qorder 4 1293 test: 1294 suffix: nonconforming_simplex_3_fv 1295 requires: ctetgen 1296 nsize: 4 1297 args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -dm_plex_dim 3 -num_comp 3 1298 1299 TEST*/ 1300 1301 /* 1302 # 2D Q_2 on a quadrilaterial Plex 1303 test: 1304 suffix: q2_2d_plex_0 1305 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -convergence 1306 test: 1307 suffix: q2_2d_plex_1 1308 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 1309 test: 1310 suffix: q2_2d_plex_2 1311 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 1312 test: 1313 suffix: q2_2d_plex_3 1314 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 -shear_coords 1315 test: 1316 suffix: q2_2d_plex_4 1317 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 -shear_coords 1318 test: 1319 suffix: q2_2d_plex_5 1320 args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 0 -non_affine_coords 1321 test: 1322 suffix: q2_2d_plex_6 1323 args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 1 -non_affine_coords 1324 test: 1325 suffix: q2_2d_plex_7 1326 args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 2 -non_affine_coords 1327 1328 test: 1329 suffix: p1d_2d_6 1330 requires: pragmatic 1331 args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0 1332 test: 1333 suffix: p1d_2d_7 1334 requires: pragmatic 1335 args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0 1336 test: 1337 suffix: p1d_2d_8 1338 requires: pragmatic 1339 args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0 1340 */ 1341