1 static char help[] = "Particle Basis Landau Example using nonlinear solve + Implicit Midpoint-like time stepping."; 2 3 /* TODO 4 5 1) SNES is sensitive to epsilon (but not to h). Should we do continuation in it? 6 7 2) Put this timestepper in library, maybe by changing DG 8 9 3) Add monitor to visualize distributions 10 11 */ 12 13 /* References 14 [1] https://arxiv.org/abs/1910.03080v2 15 */ 16 17 #include <petscdmplex.h> 18 #include <petscdmswarm.h> 19 #include <petscts.h> 20 #include <petscviewer.h> 21 #include <petscmath.h> 22 23 typedef struct { 24 /* Velocity space grid and functions */ 25 PetscInt N; /* The number of partices per spatial cell */ 26 PetscReal L; /* Velocity space is [-L, L]^d */ 27 PetscReal h; /* Spacing for grid 2L / N^{1/d} */ 28 PetscReal epsilon; /* gaussian regularization parameter */ 29 PetscReal momentTol; /* Tolerance for checking moment conservation */ 30 } AppCtx; 31 32 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) { 33 PetscFunctionBeginUser; 34 options->N = 1; 35 options->momentTol = 100.0 * PETSC_MACHINE_EPSILON; 36 options->L = 1.0; 37 options->h = -1.0; 38 options->epsilon = -1.0; 39 40 PetscOptionsBegin(comm, "", "Collision Options", "DMPLEX"); 41 PetscCall(PetscOptionsInt("-N", "Number of particles per spatial cell", "ex27.c", options->N, &options->N, NULL)); 42 PetscCall(PetscOptionsReal("-L", "Velocity-space extent", "ex27.c", options->L, &options->L, NULL)); 43 PetscCall(PetscOptionsReal("-h", "Velocity-space resolution", "ex27.c", options->h, &options->h, NULL)); 44 PetscCall(PetscOptionsReal("-epsilon", "Mollifier regularization parameter", "ex27.c", options->epsilon, &options->epsilon, NULL)); 45 PetscOptionsEnd(); 46 PetscFunctionReturn(0); 47 } 48 49 static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user) { 50 PetscFunctionBeginUser; 51 PetscCall(DMCreate(comm, dm)); 52 PetscCall(DMSetType(*dm, DMPLEX)); 53 PetscCall(DMSetFromOptions(*dm)); 54 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 55 PetscFunctionReturn(0); 56 } 57 58 static PetscErrorCode SetInitialCoordinates(DM sw) { 59 AppCtx *user; 60 PetscRandom rnd, rndv; 61 DM dm; 62 DMPolytopeType ct; 63 PetscBool simplex; 64 PetscReal *centroid, *coords, *velocity, *xi0, *v0, *J, *invJ, detJ, *vals; 65 PetscInt dim, d, cStart, cEnd, c, Np, p; 66 67 PetscFunctionBeginUser; 68 /* Randomization for coordinates */ 69 PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)sw), &rnd)); 70 PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0)); 71 PetscCall(PetscRandomSetFromOptions(rnd)); 72 PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)sw), &rndv)); 73 PetscCall(PetscRandomSetInterval(rndv, -1., 1.)); 74 PetscCall(PetscRandomSetFromOptions(rndv)); 75 PetscCall(DMGetApplicationContext(sw, &user)); 76 Np = user->N; 77 PetscCall(DMGetDimension(sw, &dim)); 78 PetscCall(DMSwarmGetCellDM(sw, &dm)); 79 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 80 PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 81 simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; 82 PetscCall(PetscMalloc5(dim, ¢roid, dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ)); 83 for (d = 0; d < dim; ++d) xi0[d] = -1.0; 84 PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 85 PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&velocity)); 86 PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&vals)); 87 for (c = cStart; c < cEnd; ++c) { 88 if (Np == 1) { 89 PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL)); 90 for (d = 0; d < dim; ++d) { coords[c * dim + d] = centroid[d]; } 91 vals[c] = 1.0; 92 } else { 93 PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); /* affine */ 94 for (p = 0; p < Np; ++p) { 95 const PetscInt n = c * Np + p; 96 PetscReal sum = 0.0, refcoords[3]; 97 98 for (d = 0; d < dim; ++d) { 99 PetscCall(PetscRandomGetValueReal(rnd, &refcoords[d])); 100 sum += refcoords[d]; 101 } 102 if (simplex && sum > 0.0) 103 for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim) * sum; 104 vals[n] = 1.0; 105 PetscCall(DMPlexReferenceToCoordinates(dm, c, 1, refcoords, &coords[n * dim])); 106 } 107 } 108 } 109 /* Random velocity IC */ 110 for (c = cStart; c < cEnd; ++c) { 111 for (p = 0; p < Np; ++p) { 112 for (d = 0; d < dim; ++d) { 113 PetscReal v_val; 114 115 PetscCall(PetscRandomGetValueReal(rndv, &v_val)); 116 velocity[p * dim + d] = v_val; 117 } 118 } 119 } 120 PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 121 PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&velocity)); 122 PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&vals)); 123 PetscCall(PetscFree5(centroid, xi0, v0, J, invJ)); 124 PetscCall(PetscRandomDestroy(&rnd)); 125 PetscCall(PetscRandomDestroy(&rndv)); 126 PetscFunctionReturn(0); 127 } 128 129 /* Get velocities from swarm and place in solution vector */ 130 static PetscErrorCode SetInitialConditions(DM dmSw, Vec u) { 131 DM dm; 132 AppCtx *user; 133 PetscReal *velocity; 134 PetscScalar *initialConditions; 135 PetscInt dim, d, cStart, cEnd, c, Np, p, n; 136 137 PetscFunctionBeginUser; 138 PetscCall(VecGetLocalSize(u, &n)); 139 PetscCall(DMGetApplicationContext(dmSw, &user)); 140 Np = user->N; 141 PetscCall(DMSwarmGetCellDM(dmSw, &dm)); 142 PetscCall(DMGetDimension(dm, &dim)); 143 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 144 PetscCall(DMSwarmGetField(dmSw, "velocity", NULL, NULL, (void **)&velocity)); 145 PetscCall(VecGetArray(u, &initialConditions)); 146 for (c = cStart; c < cEnd; ++c) { 147 for (p = 0; p < Np; ++p) { 148 const PetscInt n = c * Np + p; 149 for (d = 0; d < dim; d++) { initialConditions[n * dim + d] = velocity[n * dim + d]; } 150 } 151 } 152 PetscCall(VecRestoreArray(u, &initialConditions)); 153 PetscCall(DMSwarmRestoreField(dmSw, "velocity", NULL, NULL, (void **)&velocity)); 154 PetscFunctionReturn(0); 155 } 156 157 static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user) { 158 PetscInt *cellid; 159 PetscInt dim, cStart, cEnd, c, Np = user->N, p; 160 PetscBool view = PETSC_FALSE; 161 162 PetscFunctionBeginUser; 163 PetscCall(DMGetDimension(dm, &dim)); 164 PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw)); 165 PetscCall(DMSetType(*sw, DMSWARM)); 166 PetscCall(DMSetDimension(*sw, dim)); 167 /* h = 2L/n and N = n^d */ 168 if (user->h < 0.) user->h = 2. * user->L / PetscPowReal(user->N, 1. / dim); 169 /* From Section 4 in [1], \epsilon = 0.64 h^.98 */ 170 if (user->epsilon < 0.) user->epsilon = 0.64 * pow(user->h, 1.98); 171 PetscCall(PetscOptionsGetBool(NULL, NULL, "-param_view", &view, NULL)); 172 if (view) { PetscCall(PetscPrintf(PETSC_COMM_SELF, "N: %" PetscInt_FMT " L: %g h: %g eps: %g\n", user->N, (double)user->L, (double)user->h, (double)user->epsilon)); } 173 PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC)); 174 PetscCall(DMSwarmSetCellDM(*sw, dm)); 175 PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "velocity", dim, PETSC_REAL)); 176 PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR)); 177 PetscCall(DMSwarmFinalizeFieldRegister(*sw)); 178 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 179 PetscCall(DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0)); 180 PetscCall(DMSetFromOptions(*sw)); 181 PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 182 for (c = cStart; c < cEnd; ++c) { 183 for (p = 0; p < Np; ++p) { 184 const PetscInt n = c * Np + p; 185 cellid[n] = c; 186 } 187 } 188 PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 189 PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles")); 190 PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view")); 191 PetscFunctionReturn(0); 192 } 193 194 /* Internal dmplex function, same as found in dmpleximpl.h */ 195 static void DMPlex_WaxpyD_Internal(PetscInt dim, PetscReal a, const PetscReal *x, const PetscReal *y, PetscReal *w) { 196 PetscInt d; 197 198 for (d = 0; d < dim; ++d) w[d] = a * x[d] + y[d]; 199 } 200 201 /* Internal dmplex function, same as found in dmpleximpl.h */ 202 static PetscReal DMPlex_DotD_Internal(PetscInt dim, const PetscScalar *x, const PetscReal *y) { 203 PetscReal sum = 0.0; 204 PetscInt d; 205 206 for (d = 0; d < dim; ++d) sum += PetscRealPart(x[d]) * y[d]; 207 return sum; 208 } 209 210 /* Internal dmplex function, same as found in dmpleximpl.h */ 211 static void DMPlex_MultAdd2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[]) { 212 PetscScalar z[2]; 213 z[0] = x[0]; 214 z[1] = x[ldx]; 215 y[0] += A[0] * z[0] + A[1] * z[1]; 216 y[ldx] += A[2] * z[0] + A[3] * z[1]; 217 (void)PetscLogFlops(6.0); 218 } 219 220 /* Internal dmplex function, same as found in dmpleximpl.h to avoid private includes. */ 221 static void DMPlex_MultAdd3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[]) { 222 PetscScalar z[3]; 223 z[0] = x[0]; 224 z[1] = x[ldx]; 225 z[2] = x[ldx * 2]; 226 y[0] += A[0] * z[0] + A[1] * z[1] + A[2] * z[2]; 227 y[ldx] += A[3] * z[0] + A[4] * z[1] + A[5] * z[2]; 228 y[ldx * 2] += A[6] * z[0] + A[7] * z[1] + A[8] * z[2]; 229 (void)PetscLogFlops(15.0); 230 } 231 232 /* 233 Gaussian - The Gaussian function G(x) 234 235 Input Parameters: 236 + dim - The number of dimensions, or size of x 237 . mu - The mean, or center 238 . sigma - The standard deviation, or width 239 - x - The evaluation point of the function 240 241 Output Parameter: 242 . ret - The value G(x) 243 */ 244 static PetscReal Gaussian(PetscInt dim, const PetscReal mu[], PetscReal sigma, const PetscReal x[]) { 245 PetscReal arg = 0.0; 246 PetscInt d; 247 248 for (d = 0; d < dim; ++d) arg += PetscSqr(x[d] - mu[d]); 249 return PetscPowReal(2.0 * PETSC_PI * sigma, -dim / 2.0) * PetscExpReal(-arg / (2.0 * sigma)); 250 } 251 252 /* 253 ComputeGradS - Compute grad_v dS_eps/df 254 255 Input Parameters: 256 + dim - The dimension 257 . Np - The number of particles 258 . vp - The velocity v_p of the particle at which we evaluate 259 . velocity - The velocity field for all particles 260 . epsilon - The regularization strength 261 262 Output Parameter: 263 . integral - The output grad_v dS_eps/df (v_p) 264 265 Note: 266 This comes from (3.6) in [1], and we are computing 267 $ \nabla_v S_p = \grad \psi_\epsilon(v_p - v) log \sum_q \psi_\epsilon(v - v_q) 268 which is discretized by using a one-point quadrature in each box l at its center v^c_l 269 $ \sum_l h^d \nabla\psi_\epsilon(v_p - v^c_l) \log\left( \sum_q w_q \psi_\epsilon(v^c_l - v_q) \right) 270 where h^d is the volume of each box. 271 */ 272 static PetscErrorCode ComputeGradS(PetscInt dim, PetscInt Np, const PetscReal vp[], const PetscReal velocity[], PetscReal integral[], AppCtx *ctx) { 273 PetscReal vc_l[3], L = ctx->L, h = ctx->h, epsilon = ctx->epsilon, init = 0.5 * h - L; 274 PetscInt nx = roundf(2. * L / h); 275 PetscInt ny = dim > 1 ? nx : 1; 276 PetscInt nz = dim > 2 ? nx : 1; 277 PetscInt i, j, k, d, q, dbg = 0; 278 279 PetscFunctionBeginHot; 280 for (d = 0; d < dim; ++d) integral[d] = 0.0; 281 for (k = 0, vc_l[2] = init; k < nz; ++k, vc_l[2] += h) { 282 for (j = 0, vc_l[1] = init; j < ny; ++j, vc_l[1] += h) { 283 for (i = 0, vc_l[0] = init; i < nx; ++i, vc_l[0] += h) { 284 PetscReal sum = 0.0; 285 286 if (dbg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "(%" PetscInt_FMT " %" PetscInt_FMT ") vc_l: %g %g\n", i, j, (double)vc_l[0], (double)vc_l[1])); 287 /* \log \sum_k \psi(v - v_k) */ 288 for (q = 0; q < Np; ++q) sum += Gaussian(dim, &velocity[q * dim], epsilon, vc_l); 289 sum = PetscLogReal(sum); 290 for (d = 0; d < dim; ++d) integral[d] += (-1. / (epsilon)) * PetscAbsReal(vp[d] - vc_l[d]) * (Gaussian(dim, vp, epsilon, vc_l)) * sum; 291 } 292 } 293 } 294 PetscFunctionReturn(0); 295 } 296 297 /* Q = 1/|xi| (I - xi xi^T / |xi|^2), xi = vp - vq */ 298 static PetscErrorCode QCompute(PetscInt dim, const PetscReal vp[], const PetscReal vq[], PetscReal Q[]) { 299 PetscReal xi[3], xi2, xi3, mag; 300 PetscInt d, e; 301 302 PetscFunctionBeginHot; 303 DMPlex_WaxpyD_Internal(dim, -1.0, vq, vp, xi); 304 xi2 = DMPlex_DotD_Internal(dim, xi, xi); 305 mag = PetscSqrtReal(xi2); 306 xi3 = xi2 * mag; 307 for (d = 0; d < dim; ++d) { 308 for (e = 0; e < dim; ++e) { Q[d * dim + e] = -xi[d] * xi[e] / xi3; } 309 Q[d * dim + d] += 1. / mag; 310 } 311 PetscFunctionReturn(0); 312 } 313 314 static PetscErrorCode RHSFunctionParticles(TS ts, PetscReal t, Vec U, Vec R, void *ctx) { 315 AppCtx *user = (AppCtx *)ctx; 316 PetscInt dbg = 0; 317 DM sw; /* Particles */ 318 Vec sol; /* Solution vector at current time */ 319 const PetscScalar *u; /* input solution vector */ 320 PetscScalar *r; 321 PetscReal *velocity; 322 PetscInt dim, Np, p, q; 323 324 PetscFunctionBeginUser; 325 PetscCall(VecZeroEntries(R)); 326 PetscCall(TSGetDM(ts, &sw)); 327 PetscCall(DMGetDimension(sw, &dim)); 328 PetscCall(VecGetLocalSize(U, &Np)); 329 PetscCall(TSGetSolution(ts, &sol)); 330 PetscCall(VecGetArray(sol, &velocity)); 331 PetscCall(VecGetArray(R, &r)); 332 PetscCall(VecGetArrayRead(U, &u)); 333 Np /= dim; 334 if (dbg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Part ppr x y\n")); 335 for (p = 0; p < Np; ++p) { 336 PetscReal gradS_p[3] = {0., 0., 0.}; 337 338 PetscCall(ComputeGradS(dim, Np, &velocity[p * dim], velocity, gradS_p, user)); 339 for (q = 0; q < Np; ++q) { 340 PetscReal gradS_q[3] = {0., 0., 0.}, GammaS[3] = {0., 0., 0.}, Q[9]; 341 342 if (q == p) continue; 343 PetscCall(ComputeGradS(dim, Np, &velocity[q * dim], velocity, gradS_q, user)); 344 DMPlex_WaxpyD_Internal(dim, -1.0, gradS_q, gradS_p, GammaS); 345 PetscCall(QCompute(dim, &u[p * dim], &u[q * dim], Q)); 346 switch (dim) { 347 case 2: DMPlex_MultAdd2DReal_Internal(Q, 1, GammaS, &r[p * dim]); break; 348 case 3: DMPlex_MultAdd3DReal_Internal(Q, 1, GammaS, &r[p * dim]); break; 349 default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Do not support dimension %" PetscInt_FMT, dim); 350 } 351 } 352 if (dbg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final %4" PetscInt_FMT " %10.8lf %10.8lf\n", p, r[p * dim + 0], r[p * dim + 1])); 353 } 354 PetscCall(VecRestoreArrayRead(U, &u)); 355 PetscCall(VecRestoreArray(R, &r)); 356 PetscCall(VecRestoreArray(sol, &velocity)); 357 PetscCall(VecViewFromOptions(R, NULL, "-residual_view")); 358 PetscFunctionReturn(0); 359 } 360 361 /* 362 TS Post Step Function. Copy the solution back into the swarm for migration. We may also need to reform 363 the solution vector in cases of particle migration, but we forgo that here since there is no velocity space grid 364 to migrate between. 365 */ 366 static PetscErrorCode UpdateSwarm(TS ts) { 367 PetscInt idx, n; 368 const PetscScalar *u; 369 PetscScalar *velocity; 370 DM sw; 371 Vec sol; 372 373 PetscFunctionBeginUser; 374 PetscCall(TSGetDM(ts, &sw)); 375 PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&velocity)); 376 PetscCall(TSGetSolution(ts, &sol)); 377 PetscCall(VecGetArrayRead(sol, &u)); 378 PetscCall(VecGetLocalSize(sol, &n)); 379 for (idx = 0; idx < n; ++idx) velocity[idx] = u[idx]; 380 PetscCall(VecRestoreArrayRead(sol, &u)); 381 PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&velocity)); 382 PetscFunctionReturn(0); 383 } 384 385 static PetscErrorCode InitializeSolve(TS ts, Vec u) { 386 DM dm; 387 AppCtx *user; 388 389 PetscFunctionBeginUser; 390 PetscCall(TSGetDM(ts, &dm)); 391 PetscCall(DMGetApplicationContext(dm, &user)); 392 PetscCall(SetInitialCoordinates(dm)); 393 PetscCall(SetInitialConditions(dm, u)); 394 PetscFunctionReturn(0); 395 } 396 397 int main(int argc, char **argv) { 398 TS ts; /* nonlinear solver */ 399 DM dm, sw; /* Velocity space mesh and Particle Swarm */ 400 Vec u, v; /* problem vector */ 401 MPI_Comm comm; 402 AppCtx user; 403 404 PetscFunctionBeginUser; 405 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 406 comm = PETSC_COMM_WORLD; 407 PetscCall(ProcessOptions(comm, &user)); 408 /* Initialize objects and set initial conditions */ 409 PetscCall(CreateMesh(comm, &dm, &user)); 410 PetscCall(CreateParticles(dm, &sw, &user)); 411 PetscCall(DMSetApplicationContext(sw, &user)); 412 PetscCall(DMSwarmVectorDefineField(sw, "velocity")); 413 PetscCall(TSCreate(comm, &ts)); 414 PetscCall(TSSetDM(ts, sw)); 415 PetscCall(TSSetMaxTime(ts, 10.0)); 416 PetscCall(TSSetTimeStep(ts, 0.1)); 417 PetscCall(TSSetMaxSteps(ts, 1)); 418 PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 419 PetscCall(TSSetRHSFunction(ts, NULL, RHSFunctionParticles, &user)); 420 PetscCall(TSSetFromOptions(ts)); 421 PetscCall(TSSetComputeInitialCondition(ts, InitializeSolve)); 422 PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "velocity", &v)); 423 PetscCall(VecDuplicate(v, &u)); 424 PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "velocity", &v)); 425 PetscCall(TSComputeInitialCondition(ts, u)); 426 PetscCall(TSSetPostStep(ts, UpdateSwarm)); 427 PetscCall(TSSolve(ts, u)); 428 PetscCall(VecDestroy(&u)); 429 PetscCall(TSDestroy(&ts)); 430 PetscCall(DMDestroy(&sw)); 431 PetscCall(DMDestroy(&dm)); 432 PetscCall(PetscFinalize()); 433 return 0; 434 } 435 436 /*TEST 437 build: 438 requires: triangle !single !complex 439 test: 440 suffix: midpoint 441 args: -N 3 -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1,-1 -dm_plex_box_upper 1,1 -dm_view \ 442 -ts_type theta -ts_theta_theta 0.5 -ts_dmswarm_monitor_moments -ts_monitor_frequency 1 -snes_fd 443 TEST*/ 444