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