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