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