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