1 static char help[] = "Tests L2 projection with DMSwarm using delta function particles.\n"; 2 3 #include <petscdmplex.h> 4 #include <petscfe.h> 5 #include <petscdmswarm.h> 6 #include <petscds.h> 7 #include <petscksp.h> 8 #include <petsc/private/petscfeimpl.h> 9 typedef struct { 10 char meshFilename[PETSC_MAX_PATH_LEN]; /* Name of the mesh filename if any */ 11 PetscReal L[3]; /* Dimensions of the mesh bounding box */ 12 PetscInt particlesPerCell; /* The number of partices per cell */ 13 PetscReal particleRelDx; /* Relative particle position perturbation compared to average cell diameter h */ 14 PetscReal meshRelDx; /* Relative vertex position perturbation compared to average cell diameter h */ 15 PetscInt k; /* Mode number for test function */ 16 PetscReal momentTol; /* Tolerance for checking moment conservation */ 17 PetscBool useBlockDiagPrec; /* Use the block diagonal of the normal equations as a preconditioner */ 18 PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *); 19 } AppCtx; 20 21 /* const char *const ex2FunctionTypes[] = {"linear","x2_x4","sin","ex2FunctionTypes","EX2_FUNCTION_",0}; */ 22 23 static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx) 24 { 25 AppCtx *ctx = (AppCtx *) a_ctx; 26 PetscInt d; 27 28 u[0] = 0.0; 29 for (d = 0; d < dim; ++d) u[0] += x[d]/(ctx->L[d]); 30 return 0; 31 } 32 33 static PetscErrorCode x2_x4(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx) 34 { 35 AppCtx *ctx = (AppCtx *) a_ctx; 36 PetscInt d; 37 38 u[0] = 1; 39 for (d = 0; d < dim; ++d) u[0] *= PetscSqr(x[d])*PetscSqr(ctx->L[d]) - PetscPowRealInt(x[d], 4); 40 return 0; 41 } 42 43 static PetscErrorCode sinx(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx) 44 { 45 AppCtx *ctx = (AppCtx *) a_ctx; 46 47 u[0] = PetscSinScalar(2*PETSC_PI*ctx->k*x[0]/(ctx->L[0])); 48 return 0; 49 } 50 51 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 52 { 53 char fstring[PETSC_MAX_PATH_LEN] = "linear"; 54 PetscBool flag; 55 PetscErrorCode ierr; 56 57 PetscFunctionBeginUser; 58 options->particlesPerCell = 1; 59 options->k = 1; 60 options->particleRelDx = 1.e-20; 61 options->meshRelDx = 1.e-20; 62 options->momentTol = 100.*PETSC_MACHINE_EPSILON; 63 options->useBlockDiagPrec = PETSC_FALSE; 64 ierr = PetscStrcpy(options->meshFilename, "");CHKERRQ(ierr); 65 66 ierr = PetscOptionsBegin(comm, "", "L2 Projection Options", "DMPLEX");CHKERRQ(ierr); 67 ierr = PetscOptionsString("-mesh", "Name of the mesh filename if any", "ex2.c", options->meshFilename, options->meshFilename, sizeof(options->meshFilename), NULL);CHKERRQ(ierr); 68 ierr = PetscOptionsInt("-k", "Mode number of test", "ex2.c", options->k, &options->k, NULL);CHKERRQ(ierr); 69 ierr = PetscOptionsInt("-particlesPerCell", "Number of particles per cell", "ex2.c", options->particlesPerCell, &options->particlesPerCell, NULL);CHKERRQ(ierr); 70 ierr = PetscOptionsReal("-particle_perturbation", "Relative perturbation of particles (0,1)", "ex2.c", options->particleRelDx, &options->particleRelDx, NULL);CHKERRQ(ierr); 71 ierr = PetscOptionsReal("-mesh_perturbation", "Relative perturbation of mesh points (0,1)", "ex2.c", options->meshRelDx, &options->meshRelDx, NULL);CHKERRQ(ierr); 72 ierr = PetscOptionsString("-function", "Name of test function", "ex2.c", fstring, fstring, sizeof(fstring), NULL);CHKERRQ(ierr); 73 ierr = PetscStrcmp(fstring, "linear", &flag);CHKERRQ(ierr); 74 if (flag) { 75 options->func = linear; 76 } else { 77 ierr = PetscStrcmp(fstring, "sin", &flag);CHKERRQ(ierr); 78 if (flag) { 79 options->func = sinx; 80 } else { 81 ierr = PetscStrcmp(fstring, "x2_x4", &flag);CHKERRQ(ierr); 82 options->func = x2_x4; 83 if (!flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown function %s",fstring); 84 } 85 } 86 ierr = PetscOptionsBool("-block_diag_prec", "Use the block diagonal of the normal equations to precondition the particle projection", "ex2.c", options->useBlockDiagPrec, &options->useBlockDiagPrec, NULL);CHKERRQ(ierr); 87 ierr = PetscOptionsEnd(); 88 89 PetscFunctionReturn(0); 90 } 91 92 static PetscErrorCode PerturbVertices(DM dm, AppCtx *user) 93 { 94 PetscRandom rnd; 95 PetscReal interval = user->meshRelDx; 96 Vec coordinates; 97 PetscScalar *coords; 98 PetscReal *hh, low[3], high[3]; 99 PetscInt d, cdim, cEnd, N, p, bs; 100 PetscErrorCode ierr; 101 102 PetscFunctionBeginUser; 103 ierr = DMGetBoundingBox(dm, low, high);CHKERRQ(ierr); 104 ierr = DMPlexGetHeightStratum(dm, 0, NULL, &cEnd);CHKERRQ(ierr); 105 ierr = PetscRandomCreate(PetscObjectComm((PetscObject) dm), &rnd);CHKERRQ(ierr); 106 ierr = PetscRandomSetInterval(rnd, -interval, interval);CHKERRQ(ierr); 107 ierr = PetscRandomSetFromOptions(rnd);CHKERRQ(ierr); 108 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 109 ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr); 110 ierr = PetscCalloc1(cdim,&hh);CHKERRQ(ierr); 111 for (d = 0; d < cdim; ++d) hh[d] = (user->L[d])/PetscPowReal(cEnd, 1./cdim); 112 ierr = VecGetLocalSize(coordinates, &N);CHKERRQ(ierr); 113 ierr = VecGetBlockSize(coordinates, &bs);CHKERRQ(ierr); 114 if (bs != cdim) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Coordinate vector has wrong block size %D != %D", bs, cdim); 115 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 116 for (p = 0; p < N; p += cdim) { 117 PetscScalar *coord = &coords[p], value; 118 119 for (d = 0; d < cdim; ++d) { 120 ierr = PetscRandomGetValue(rnd, &value);CHKERRQ(ierr); 121 coord[d] = PetscMax(low[d], PetscMin(high[d], PetscRealPart(coord[d] + value*hh[d]))); 122 } 123 } 124 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 125 ierr = PetscRandomDestroy(&rnd);CHKERRQ(ierr); 126 ierr = PetscFree(hh);CHKERRQ(ierr); 127 PetscFunctionReturn(0); 128 } 129 130 static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user) 131 { 132 PetscReal low[3], high[3]; 133 PetscInt cdim, d; 134 PetscBool flg; 135 PetscErrorCode ierr; 136 137 PetscFunctionBeginUser; 138 ierr = PetscStrcmp(user->meshFilename, "", &flg);CHKERRQ(ierr); 139 if (flg) { 140 ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr); 141 } else { 142 ierr = DMPlexCreateFromFile(comm, user->meshFilename, PETSC_TRUE, dm);CHKERRQ(ierr); 143 } 144 ierr = DMGetCoordinateDim(*dm, &cdim);CHKERRQ(ierr); 145 ierr = DMGetBoundingBox(*dm, low, high);CHKERRQ(ierr); 146 for (d = 0; d < cdim; ++d) user->L[d] = high[d] - low[d]; 147 { 148 DM distributedMesh = NULL; 149 150 ierr = DMPlexDistribute(*dm, 0, NULL, &distributedMesh);CHKERRQ(ierr); 151 if (distributedMesh) { 152 ierr = DMDestroy(dm);CHKERRQ(ierr); 153 *dm = distributedMesh; 154 } 155 } 156 ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr); /* needed for periodic */ 157 ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 158 ierr = PerturbVertices(*dm, user);CHKERRQ(ierr); 159 ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr); 160 ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 161 PetscFunctionReturn(0); 162 } 163 164 static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux, 165 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 166 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 167 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 168 { 169 g0[0] = 1.0; 170 } 171 172 static PetscErrorCode CreateFEM(DM dm, AppCtx *user) 173 { 174 PetscFE fe; 175 PetscDS ds; 176 DMPolytopeType ct; 177 PetscBool simplex; 178 PetscInt dim, cStart; 179 PetscErrorCode ierr; 180 181 PetscFunctionBeginUser; 182 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 183 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr); 184 ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr); 185 simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 186 ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, simplex, NULL, -1, &fe);CHKERRQ(ierr); 187 ierr = PetscObjectSetName((PetscObject) fe, "fe");CHKERRQ(ierr); 188 ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 189 ierr = DMCreateDS(dm);CHKERRQ(ierr); 190 ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 191 /* Setup to form mass matrix */ 192 ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 193 ierr = PetscDSSetJacobian(ds, 0, 0, identity, NULL, NULL, NULL);CHKERRQ(ierr); 194 PetscFunctionReturn(0); 195 } 196 197 static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user) 198 { 199 PetscRandom rnd, rndp; 200 PetscReal interval = user->particleRelDx; 201 DMPolytopeType ct; 202 PetscBool simplex; 203 PetscScalar value, *vals; 204 PetscReal *centroid, *coords, *xi0, *v0, *J, *invJ, detJ; 205 PetscInt *cellid; 206 PetscInt Ncell, Np = user->particlesPerCell, p, cStart, c, dim, d; 207 PetscErrorCode ierr; 208 209 PetscFunctionBeginUser; 210 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 211 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr); 212 ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr); 213 simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 214 ierr = DMCreate(PetscObjectComm((PetscObject) dm), sw);CHKERRQ(ierr); 215 ierr = DMSetType(*sw, DMSWARM);CHKERRQ(ierr); 216 ierr = DMSetDimension(*sw, dim);CHKERRQ(ierr); 217 218 ierr = PetscRandomCreate(PetscObjectComm((PetscObject) dm), &rnd);CHKERRQ(ierr); 219 ierr = PetscRandomSetInterval(rnd, -1.0, 1.0);CHKERRQ(ierr); 220 ierr = PetscRandomSetFromOptions(rnd);CHKERRQ(ierr); 221 ierr = PetscRandomCreate(PetscObjectComm((PetscObject) dm), &rndp);CHKERRQ(ierr); 222 ierr = PetscRandomSetInterval(rndp, -interval, interval);CHKERRQ(ierr); 223 ierr = PetscRandomSetFromOptions(rndp);CHKERRQ(ierr); 224 225 ierr = DMSwarmSetType(*sw, DMSWARM_PIC);CHKERRQ(ierr); 226 ierr = DMSwarmSetCellDM(*sw, dm);CHKERRQ(ierr); 227 ierr = DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR);CHKERRQ(ierr); 228 ierr = DMSwarmFinalizeFieldRegister(*sw);CHKERRQ(ierr); 229 ierr = DMPlexGetHeightStratum(dm, 0, NULL, &Ncell);CHKERRQ(ierr); 230 ierr = DMSwarmSetLocalSizes(*sw, Ncell * Np, 0);CHKERRQ(ierr); 231 ierr = DMSetFromOptions(*sw);CHKERRQ(ierr); 232 ierr = DMSwarmGetField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 233 ierr = DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr); 234 ierr = DMSwarmGetField(*sw, "w_q", NULL, NULL, (void **) &vals);CHKERRQ(ierr); 235 236 ierr = PetscMalloc5(dim, ¢roid, dim, &xi0, dim, &v0, dim*dim, &J, dim*dim, &invJ);CHKERRQ(ierr); 237 for (c = 0; c < Ncell; ++c) { 238 if (Np == 1) { 239 ierr = DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL);CHKERRQ(ierr); 240 cellid[c] = c; 241 for (d = 0; d < dim; ++d) coords[c*dim+d] = centroid[d]; 242 } else { 243 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); /* affine */ 244 for (d = 0; d < dim; ++d) xi0[d] = -1.0; 245 for (p = 0; p < Np; ++p) { 246 const PetscInt n = c*Np + p; 247 PetscReal sum = 0.0, refcoords[3]; 248 249 cellid[n] = c; 250 for (d = 0; d < dim; ++d) {ierr = PetscRandomGetValue(rnd, &value);CHKERRQ(ierr); refcoords[d] = PetscRealPart(value); sum += refcoords[d];} 251 if (simplex && sum > 0.0) for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim)*sum; 252 CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n*dim]); 253 } 254 } 255 } 256 ierr = PetscFree5(centroid, xi0, v0, J, invJ);CHKERRQ(ierr); 257 for (c = 0; c < Ncell; ++c) { 258 for (p = 0; p < Np; ++p) { 259 const PetscInt n = c*Np + p; 260 261 for (d = 0; d < dim; ++d) {ierr = PetscRandomGetValue(rndp, &value);CHKERRQ(ierr); coords[n*dim+d] += PetscRealPart(value);} 262 user->func(dim, 0.0, &coords[n*dim], 1, &vals[c], user); 263 } 264 } 265 ierr = DMSwarmRestoreField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 266 ierr = DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr); 267 ierr = DMSwarmRestoreField(*sw, "w_q", NULL, NULL, (void **) &vals);CHKERRQ(ierr); 268 ierr = PetscRandomDestroy(&rnd);CHKERRQ(ierr); 269 ierr = PetscRandomDestroy(&rndp);CHKERRQ(ierr); 270 ierr = PetscObjectSetName((PetscObject) *sw, "Particles");CHKERRQ(ierr); 271 ierr = DMViewFromOptions(*sw, NULL, "-sw_view");CHKERRQ(ierr); 272 PetscFunctionReturn(0); 273 } 274 275 static PetscErrorCode computeParticleMoments(DM sw, PetscReal moments[3], AppCtx *user) 276 { 277 DM dm; 278 const PetscReal *coords; 279 const PetscScalar *w; 280 PetscReal mom[3] = {0.0, 0.0, 0.0}; 281 PetscInt cell, cStart, cEnd, dim; 282 PetscErrorCode ierr; 283 284 PetscFunctionBeginUser; 285 ierr = DMGetDimension(sw, &dim);CHKERRQ(ierr); 286 ierr = DMSwarmGetCellDM(sw, &dm);CHKERRQ(ierr); 287 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 288 ierr = DMSwarmSortGetAccess(sw);CHKERRQ(ierr); 289 ierr = DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 290 ierr = DMSwarmGetField(sw, "w_q", NULL, NULL, (void **) &w);CHKERRQ(ierr); 291 for (cell = cStart; cell < cEnd; ++cell) { 292 PetscInt *pidx; 293 PetscInt Np, p, d; 294 295 ierr = DMSwarmSortGetPointsPerCell(sw, cell, &Np, &pidx);CHKERRQ(ierr); 296 for (p = 0; p < Np; ++p) { 297 const PetscInt idx = pidx[p]; 298 const PetscReal *c = &coords[idx*dim]; 299 300 mom[0] += PetscRealPart(w[idx]); 301 mom[1] += PetscRealPart(w[idx]) * c[0]; 302 for (d = 0; d < dim; ++d) mom[2] += PetscRealPart(w[idx]) * c[d]*c[d]; 303 } 304 ierr = PetscFree(pidx);CHKERRQ(ierr); 305 } 306 ierr = DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 307 ierr = DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **) &w);CHKERRQ(ierr); 308 ierr = DMSwarmSortRestoreAccess(sw);CHKERRQ(ierr); 309 ierr = MPI_Allreduce(mom, moments, 3, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject) sw));CHKERRQ(ierr); 310 PetscFunctionReturn(0); 311 } 312 313 static void f0_1(PetscInt dim, PetscInt Nf, PetscInt NfAux, 314 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 315 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 316 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 317 { 318 f0[0] = u[0]; 319 } 320 321 static void f0_x(PetscInt dim, PetscInt Nf, PetscInt NfAux, 322 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 323 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 324 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 325 { 326 f0[0] = x[0]*u[0]; 327 } 328 329 static void f0_r2(PetscInt dim, PetscInt Nf, PetscInt NfAux, 330 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 331 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 332 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 333 { 334 PetscInt d; 335 336 f0[0] = 0.0; 337 for (d = 0; d < dim; ++d) f0[0] += PetscSqr(x[d])*u[0]; 338 } 339 340 static PetscErrorCode computeFEMMoments(DM dm, Vec u, PetscReal moments[3], AppCtx *user) 341 { 342 PetscDS prob; 343 PetscErrorCode ierr; 344 PetscScalar mom; 345 346 PetscFunctionBeginUser; 347 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 348 ierr = PetscDSSetObjective(prob, 0, &f0_1);CHKERRQ(ierr); 349 ierr = DMPlexComputeIntegralFEM(dm, u, &mom, user);CHKERRQ(ierr); 350 moments[0] = PetscRealPart(mom); 351 ierr = PetscDSSetObjective(prob, 0, &f0_x);CHKERRQ(ierr); 352 ierr = DMPlexComputeIntegralFEM(dm, u, &mom, user);CHKERRQ(ierr); 353 moments[1] = PetscRealPart(mom); 354 ierr = PetscDSSetObjective(prob, 0, &f0_r2);CHKERRQ(ierr); 355 ierr = DMPlexComputeIntegralFEM(dm, u, &mom, user);CHKERRQ(ierr); 356 moments[2] = PetscRealPart(mom); 357 PetscFunctionReturn(0); 358 } 359 360 static PetscErrorCode TestL2ProjectionParticlesToField(DM dm, DM sw, AppCtx *user) 361 { 362 MPI_Comm comm; 363 KSP ksp; 364 Mat M; /* FEM mass matrix */ 365 Mat M_p; /* Particle mass matrix */ 366 Vec f, rhs, fhat; /* Particle field f, \int phi_i f, FEM field */ 367 PetscReal pmoments[3]; /* \int f, \int x f, \int r^2 f */ 368 PetscReal fmoments[3]; /* \int \hat f, \int x \hat f, \int r^2 \hat f */ 369 PetscInt m; 370 PetscErrorCode ierr; 371 372 PetscFunctionBeginUser; 373 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 374 ierr = KSPCreate(comm, &ksp);CHKERRQ(ierr); 375 ierr = KSPSetOptionsPrefix(ksp, "ptof_");CHKERRQ(ierr); 376 ierr = DMGetGlobalVector(dm, &fhat);CHKERRQ(ierr); 377 ierr = DMGetGlobalVector(dm, &rhs);CHKERRQ(ierr); 378 379 ierr = DMCreateMassMatrix(sw, dm, &M_p);CHKERRQ(ierr); 380 ierr = MatViewFromOptions(M_p, NULL, "-M_p_view");CHKERRQ(ierr); 381 382 /* make particle weight vector */ 383 ierr = DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f);CHKERRQ(ierr); 384 385 /* create matrix RHS vector */ 386 ierr = MatMultTranspose(M_p, f, rhs);CHKERRQ(ierr); 387 ierr = DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f);CHKERRQ(ierr); 388 ierr = PetscObjectSetName((PetscObject) rhs,"rhs");CHKERRQ(ierr); 389 ierr = VecViewFromOptions(rhs, NULL, "-rhs_view");CHKERRQ(ierr); 390 391 ierr = DMCreateMatrix(dm, &M);CHKERRQ(ierr); 392 ierr = DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user);CHKERRQ(ierr); 393 ierr = MatViewFromOptions(M, NULL, "-M_view");CHKERRQ(ierr); 394 ierr = KSPSetOperators(ksp, M, M);CHKERRQ(ierr); 395 ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 396 ierr = KSPSolve(ksp, rhs, fhat);CHKERRQ(ierr); 397 ierr = PetscObjectSetName((PetscObject) fhat,"fhat");CHKERRQ(ierr); 398 ierr = VecViewFromOptions(fhat, NULL, "-fhat_view");CHKERRQ(ierr); 399 400 /* Check moments of field */ 401 ierr = computeParticleMoments(sw, pmoments, user);CHKERRQ(ierr); 402 ierr = computeFEMMoments(dm, fhat, fmoments, user);CHKERRQ(ierr); 403 ierr = PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2]);CHKERRQ(ierr); 404 for (m = 0; m < 3; ++m) { 405 if (PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]) > user->momentTol) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "Moment %D error too large %g > %g", m, PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]), user->momentTol); 406 } 407 408 ierr = KSPDestroy(&ksp);CHKERRQ(ierr); 409 ierr = MatDestroy(&M);CHKERRQ(ierr); 410 ierr = MatDestroy(&M_p);CHKERRQ(ierr); 411 ierr = DMRestoreGlobalVector(dm, &fhat);CHKERRQ(ierr); 412 ierr = DMRestoreGlobalVector(dm, &rhs);CHKERRQ(ierr); 413 414 PetscFunctionReturn(0); 415 } 416 417 418 static PetscErrorCode TestL2ProjectionFieldToParticles(DM dm, DM sw, AppCtx *user) 419 { 420 421 MPI_Comm comm; 422 KSP ksp; 423 Mat M; /* FEM mass matrix */ 424 Mat M_p, PM_p; /* Particle mass matrix M_p, and the preconditioning matrix, e.g. M_p M^T_p */ 425 Vec f, rhs, fhat; /* Particle field f, \int phi_i fhat, FEM field */ 426 PetscReal pmoments[3]; /* \int f, \int x f, \int r^2 f */ 427 PetscReal fmoments[3]; /* \int \hat f, \int x \hat f, \int r^2 \hat f */ 428 PetscInt m; 429 PetscErrorCode ierr; 430 431 PetscFunctionBeginUser; 432 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 433 434 ierr = KSPCreate(comm, &ksp);CHKERRQ(ierr); 435 ierr = KSPSetOptionsPrefix(ksp, "ftop_");CHKERRQ(ierr); 436 ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 437 438 ierr = DMGetGlobalVector(dm, &fhat);CHKERRQ(ierr); 439 ierr = DMGetGlobalVector(dm, &rhs);CHKERRQ(ierr); 440 441 ierr = DMCreateMassMatrix(sw, dm, &M_p);CHKERRQ(ierr); 442 ierr = MatViewFromOptions(M_p, NULL, "-M_p_view");CHKERRQ(ierr); 443 444 /* make particle weight vector */ 445 ierr = DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f);CHKERRQ(ierr); 446 447 /* create matrix RHS vector, in this case the FEM field fhat with the coefficients vector #alpha */ 448 ierr = PetscObjectSetName((PetscObject) rhs,"rhs");CHKERRQ(ierr); 449 ierr = VecViewFromOptions(rhs, NULL, "-rhs_view");CHKERRQ(ierr); 450 ierr = DMCreateMatrix(dm, &M);CHKERRQ(ierr); 451 ierr = DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user);CHKERRQ(ierr); 452 ierr = MatViewFromOptions(M, NULL, "-M_view");CHKERRQ(ierr); 453 ierr = MatMultTranspose(M, fhat, rhs);CHKERRQ(ierr); 454 if (user->useBlockDiagPrec) {ierr = DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p);CHKERRQ(ierr);} 455 else {ierr = PetscObjectReference((PetscObject) M_p);CHKERRQ(ierr); PM_p = M_p;} 456 457 ierr = KSPSetOperators(ksp, M_p, PM_p);CHKERRQ(ierr); 458 ierr = KSPSolveTranspose(ksp, rhs, f);CHKERRQ(ierr); 459 ierr = PetscObjectSetName((PetscObject) fhat,"fhat");CHKERRQ(ierr); 460 ierr = VecViewFromOptions(fhat, NULL, "-fhat_view");CHKERRQ(ierr); 461 462 ierr = DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f);CHKERRQ(ierr); 463 464 /* Check moments */ 465 ierr = computeParticleMoments(sw, pmoments, user);CHKERRQ(ierr); 466 ierr = computeFEMMoments(dm, fhat, fmoments, user);CHKERRQ(ierr); 467 ierr = PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2]);CHKERRQ(ierr); 468 for (m = 0; m < 3; ++m) { 469 if (PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]) > user->momentTol) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "Moment %D error too large %g > %g", m, PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]), user->momentTol); 470 } 471 ierr = KSPDestroy(&ksp);CHKERRQ(ierr); 472 ierr = MatDestroy(&M);CHKERRQ(ierr); 473 ierr = MatDestroy(&M_p);CHKERRQ(ierr); 474 ierr = MatDestroy(&PM_p);CHKERRQ(ierr); 475 ierr = DMRestoreGlobalVector(dm, &fhat);CHKERRQ(ierr); 476 ierr = DMRestoreGlobalVector(dm, &rhs);CHKERRQ(ierr); 477 PetscFunctionReturn(0); 478 } 479 480 /* Interpolate the gradient of an FEM (FVM) field. Code repurposed from DMPlexComputeGradientClementInterpolant */ 481 static PetscErrorCode InterpolateGradient(DM dm, Vec locX, Vec locC){ 482 483 DM_Plex *mesh = (DM_Plex *) dm->data; 484 PetscInt debug = mesh->printFEM; 485 DM dmC; 486 PetscSection section; 487 PetscQuadrature quad = NULL; 488 PetscScalar *interpolant, *gradsum; 489 PetscFEGeom fegeom; 490 PetscReal *coords; 491 const PetscReal *quadPoints, *quadWeights; 492 PetscInt dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, vStart, vEnd, v, field, fieldOffset; 493 PetscErrorCode ierr; 494 495 PetscFunctionBegin; 496 ierr = VecGetDM(locC, &dmC);CHKERRQ(ierr); 497 ierr = VecSet(locC, 0.0);CHKERRQ(ierr); 498 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 499 ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr); 500 fegeom.dimEmbed = coordDim; 501 ierr = DMGetLocalSection(dm, §ion);CHKERRQ(ierr); 502 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 503 for (field = 0; field < numFields; ++field) { 504 PetscObject obj; 505 PetscClassId id; 506 PetscInt Nc; 507 508 ierr = DMGetField(dm, field, NULL, &obj);CHKERRQ(ierr); 509 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 510 if (id == PETSCFE_CLASSID) { 511 PetscFE fe = (PetscFE) obj; 512 513 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 514 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 515 } else if (id == PETSCFV_CLASSID) { 516 PetscFV fv = (PetscFV) obj; 517 518 ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr); 519 ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 520 } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 521 numComponents += Nc; 522 } 523 ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 524 if ((qNc != 1) && (qNc != numComponents)) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Quadrature components %D != %D field components", qNc, numComponents); 525 ierr = PetscMalloc6(coordDim*numComponents*2,&gradsum,coordDim*numComponents,&interpolant,coordDim*Nq,&coords,Nq,&fegeom.detJ,coordDim*coordDim*Nq,&fegeom.J,coordDim*coordDim*Nq,&fegeom.invJ);CHKERRQ(ierr); 526 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 527 ierr = DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 528 for (v = vStart; v < vEnd; ++v) { 529 PetscInt *star = NULL; 530 PetscInt starSize, st, d, fc; 531 532 ierr = PetscArrayzero(gradsum, coordDim*numComponents);CHKERRQ(ierr); 533 ierr = DMPlexGetTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 534 for (st = 0; st < starSize*2; st += 2) { 535 const PetscInt cell = star[st]; 536 PetscScalar *grad = &gradsum[coordDim*numComponents]; 537 PetscScalar *x = NULL; 538 539 if ((cell < cStart) || (cell >= cEnd)) continue; 540 ierr = DMPlexComputeCellGeometryFEM(dm, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ);CHKERRQ(ierr); 541 ierr = DMPlexVecGetClosure(dm, NULL, locX, cell, NULL, &x);CHKERRQ(ierr); 542 for (field = 0, fieldOffset = 0; field < numFields; ++field) { 543 PetscObject obj; 544 PetscClassId id; 545 PetscInt Nb, Nc, q, qc = 0; 546 547 ierr = PetscArrayzero(grad, coordDim*numComponents);CHKERRQ(ierr); 548 ierr = DMGetField(dm, field, NULL, &obj);CHKERRQ(ierr); 549 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 550 if (id == PETSCFE_CLASSID) {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);} 551 else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;} 552 else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 553 for (q = 0; q < Nq; ++q) { 554 if (fegeom.detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, quadrature points %D", (double)fegeom.detJ[q], cell, q); 555 if (ierr) { 556 PetscErrorCode ierr2; 557 ierr2 = DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x);CHKERRQ(ierr2); 558 ierr2 = DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr2); 559 ierr2 = PetscFree6(gradsum,interpolant,coords,fegeom.detJ,fegeom.J,fegeom.invJ);CHKERRQ(ierr2); 560 CHKERRQ(ierr); 561 } 562 if (id == PETSCFE_CLASSID) {ierr = PetscFEInterpolateGradient_Static((PetscFE) obj, &x[fieldOffset], &fegeom, q, interpolant);CHKERRQ(ierr);} 563 else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 564 for (fc = 0; fc < Nc; ++fc) { 565 const PetscReal wt = quadWeights[q*qNc+qc+fc]; 566 567 for (d = 0; d < coordDim; ++d) grad[fc*coordDim+d] += interpolant[fc*dim+d]*wt*fegeom.detJ[q]; 568 } 569 } 570 fieldOffset += Nb; 571 } 572 ierr = DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x);CHKERRQ(ierr); 573 for (fc = 0; fc < numComponents; ++fc) { 574 for (d = 0; d < coordDim; ++d) { 575 gradsum[fc*coordDim+d] += grad[fc*coordDim+d]; 576 } 577 } 578 if (debug) { 579 ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D gradient: [", cell);CHKERRQ(ierr); 580 for (fc = 0; fc < numComponents; ++fc) { 581 for (d = 0; d < coordDim; ++d) { 582 if (fc || d > 0) {ierr = PetscPrintf(PETSC_COMM_SELF, ", ");CHKERRQ(ierr);} 583 ierr = PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(grad[fc*coordDim+d]));CHKERRQ(ierr); 584 } 585 } 586 ierr = PetscPrintf(PETSC_COMM_SELF, "]\n");CHKERRQ(ierr); 587 } 588 } 589 ierr = DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 590 ierr = DMPlexVecSetClosure(dmC, NULL, locC, v, gradsum, INSERT_VALUES);CHKERRQ(ierr); 591 } 592 ierr = PetscFree6(gradsum,interpolant,coords,fegeom.detJ,fegeom.J,fegeom.invJ);CHKERRQ(ierr); 593 PetscFunctionReturn(0); 594 } 595 596 static PetscErrorCode TestFieldGradientProjection(DM dm, DM sw, AppCtx *user) 597 { 598 599 MPI_Comm comm; 600 KSP ksp; 601 Mat M; /* FEM mass matrix */ 602 Mat M_p; /* Particle mass matrix */ 603 Vec f, rhs, fhat, grad; /* Particle field f, \int phi_i f, FEM field */ 604 PetscReal pmoments[3]; /* \int f, \int x f, \int r^2 f */ 605 PetscReal fmoments[3]; /* \int \hat f, \int x \hat f, \int r^2 \hat f */ 606 PetscInt m; 607 PetscErrorCode ierr; 608 609 PetscFunctionBeginUser; 610 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 611 ierr = KSPCreate(comm, &ksp);CHKERRQ(ierr); 612 ierr = KSPSetOptionsPrefix(ksp, "ptof_");CHKERRQ(ierr); 613 ierr = DMGetGlobalVector(dm, &fhat);CHKERRQ(ierr); 614 ierr = DMGetGlobalVector(dm, &rhs);CHKERRQ(ierr); 615 ierr = DMGetGlobalVector(dm, &grad);CHKERRQ(ierr); 616 617 ierr = DMCreateMassMatrix(sw, dm, &M_p);CHKERRQ(ierr); 618 ierr = MatViewFromOptions(M_p, NULL, "-M_p_view");CHKERRQ(ierr); 619 620 /* make particle weight vector */ 621 ierr = DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f);CHKERRQ(ierr); 622 623 /* create matrix RHS vector */ 624 ierr = MatMultTranspose(M_p, f, rhs);CHKERRQ(ierr); 625 ierr = DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f);CHKERRQ(ierr); 626 ierr = PetscObjectSetName((PetscObject) rhs,"rhs");CHKERRQ(ierr); 627 ierr = VecViewFromOptions(rhs, NULL, "-rhs_view");CHKERRQ(ierr); 628 629 ierr = DMCreateMatrix(dm, &M);CHKERRQ(ierr); 630 ierr = DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user);CHKERRQ(ierr); 631 632 ierr = InterpolateGradient(dm, fhat, grad);CHKERRQ(ierr); 633 634 ierr = MatViewFromOptions(M, NULL, "-M_view");CHKERRQ(ierr); 635 ierr = KSPSetOperators(ksp, M, M);CHKERRQ(ierr); 636 ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 637 ierr = KSPSolve(ksp, rhs, grad);CHKERRQ(ierr); 638 ierr = PetscObjectSetName((PetscObject) fhat,"fhat");CHKERRQ(ierr); 639 ierr = VecViewFromOptions(fhat, NULL, "-fhat_view");CHKERRQ(ierr); 640 641 /* Check moments of field */ 642 ierr = computeParticleMoments(sw, pmoments, user);CHKERRQ(ierr); 643 ierr = computeFEMMoments(dm, grad, fmoments, user);CHKERRQ(ierr); 644 ierr = PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2]);CHKERRQ(ierr); 645 for (m = 0; m < 3; ++m) { 646 if (PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]) > user->momentTol) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "Moment %D error too large %g > %g", m, PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]), user->momentTol); 647 } 648 649 ierr = KSPDestroy(&ksp);CHKERRQ(ierr); 650 ierr = MatDestroy(&M);CHKERRQ(ierr); 651 ierr = MatDestroy(&M_p);CHKERRQ(ierr); 652 ierr = DMRestoreGlobalVector(dm, &fhat);CHKERRQ(ierr); 653 ierr = DMRestoreGlobalVector(dm, &rhs);CHKERRQ(ierr); 654 ierr = DMRestoreGlobalVector(dm, &grad);CHKERRQ(ierr); 655 656 PetscFunctionReturn(0); 657 } 658 659 int main (int argc, char * argv[]) { 660 MPI_Comm comm; 661 DM dm, sw; 662 AppCtx user; 663 PetscErrorCode ierr; 664 665 ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 666 comm = PETSC_COMM_WORLD; 667 ierr = ProcessOptions(comm, &user);CHKERRQ(ierr); 668 ierr = CreateMesh(comm, &dm, &user);CHKERRQ(ierr); 669 ierr = CreateFEM(dm, &user);CHKERRQ(ierr); 670 ierr = CreateParticles(dm, &sw, &user);CHKERRQ(ierr); 671 ierr = TestL2ProjectionParticlesToField(dm, sw, &user);CHKERRQ(ierr); 672 ierr = TestL2ProjectionFieldToParticles(dm, sw, &user);CHKERRQ(ierr); 673 ierr = TestFieldGradientProjection(dm, sw, &user);CHKERRQ(ierr); 674 ierr = DMDestroy(&dm);CHKERRQ(ierr); 675 ierr = DMDestroy(&sw);CHKERRQ(ierr); 676 ierr = PetscFinalize(); 677 return ierr; 678 } 679 680 681 /*TEST 682 683 # Swarm does not handle complex or quad 684 build: 685 requires: !complex double 686 687 test: 688 suffix: proj_tri_0 689 requires: triangle 690 args: -dm_plex_box_faces 1,1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 691 filter: grep -v marker | grep -v atomic | grep -v usage 692 693 test: 694 suffix: proj_tri_2_faces 695 requires: triangle 696 args: -dm_plex_box_faces 2,2 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 697 filter: grep -v marker | grep -v atomic | grep -v usage 698 699 test: 700 suffix: proj_quad_0 701 requires: triangle 702 args: -dm_plex_box_simplex 0 -dm_plex_box_faces 1,1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 703 filter: grep -v marker | grep -v atomic | grep -v usage 704 705 test: 706 suffix: proj_quad_2_faces 707 requires: triangle 708 args: -dm_plex_box_simplex 0 -dm_plex_box_faces 2,2 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 709 filter: grep -v marker | grep -v atomic | grep -v usage 710 711 test: 712 suffix: proj_tri_5P 713 requires: triangle 714 args: -dm_plex_box_faces 1,1 -particlesPerCell 5 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 715 filter: grep -v marker | grep -v atomic | grep -v usage 716 717 test: 718 suffix: proj_quad_5P 719 requires: triangle 720 args: -dm_plex_box_simplex 0 -dm_plex_box_faces 1,1 -particlesPerCell 5 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 721 filter: grep -v marker | grep -v atomic | grep -v usage 722 723 test: 724 suffix: proj_tri_mdx 725 requires: triangle 726 args: -dm_plex_box_faces 1,1 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 727 filter: grep -v marker | grep -v atomic | grep -v usage 728 729 test: 730 suffix: proj_tri_mdx_5P 731 requires: triangle 732 args: -dm_plex_box_faces 1,1 -particlesPerCell 5 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 733 filter: grep -v marker | grep -v atomic | grep -v usage 734 735 test: 736 suffix: proj_tri_3d 737 requires: ctetgen 738 args: -dm_plex_box_dim 3 -dm_plex_box_faces 1,1,1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 739 filter: grep -v marker | grep -v atomic | grep -v usage 740 741 test: 742 suffix: proj_tri_3d_2_faces 743 requires: ctetgen 744 args: -dm_plex_box_dim 3 -dm_plex_box_faces 2,2,2 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 745 filter: grep -v marker | grep -v atomic | grep -v usage 746 747 test: 748 suffix: proj_tri_3d_5P 749 requires: ctetgen 750 args: -dm_plex_box_dim 3 -dm_plex_box_faces 1,1,1 -particlesPerCell 5 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 751 filter: grep -v marker | grep -v atomic | grep -v usage 752 753 test: 754 suffix: proj_tri_3d_mdx 755 requires: ctetgen 756 args: -dm_plex_box_dim 3 -dm_plex_box_faces 1,1,1 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 757 filter: grep -v marker | grep -v atomic | grep -v usage 758 759 test: 760 suffix: proj_tri_3d_mdx_5P 761 requires: ctetgen 762 args: -dm_plex_box_dim 3 -dm_plex_box_faces 1,1,1 -particlesPerCell 5 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 763 filter: grep -v marker | grep -v atomic | grep -v usage 764 765 test: 766 suffix: proj_tri_3d_mdx_2_faces 767 requires: ctetgen 768 args: -dm_plex_box_dim 3 -dm_plex_box_faces 2,2,2 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 769 filter: grep -v marker | grep -v atomic | grep -v usage 770 771 test: 772 suffix: proj_tri_3d_mdx_5P_2_faces 773 requires: ctetgen 774 args: -dm_plex_box_dim 3 -dm_plex_box_faces 2,2,2 -particlesPerCell 5 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 775 filter: grep -v marker | grep -v atomic | grep -v usage 776 777 test: 778 suffix: proj_quad_lsqr_scale 779 requires: !complex 780 args: -dm_plex_box_simplex 0 -dm_plex_box_faces 4,4 \ 781 -petscspace_degree 2 -petscfe_default_quadrature_order 3 \ 782 -particlesPerCell 200 \ 783 -ptof_pc_type lu \ 784 -ftop_ksp_rtol 1e-17 -ftop_ksp_type lsqr -ftop_pc_type none 785 786 test: 787 suffix: proj_quad_lsqr_prec_scale 788 requires: !complex 789 args: -dm_plex_box_simplex 0 -dm_plex_box_faces 4,4 \ 790 -petscspace_degree 2 -petscfe_default_quadrature_order 3 \ 791 -particlesPerCell 200 \ 792 -ptof_pc_type lu \ 793 -ftop_ksp_rtol 1e-17 -ftop_ksp_type lsqr -ftop_pc_type lu -ftop_pc_factor_shift_type nonzero -block_diag_prec 794 795 TEST*/ 796