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