1 static char help[] = "Vlasov example of central orbits\n"; 2 3 /* 4 To visualize the orbit, we can used 5 6 -ts_monitor_sp_swarm -ts_monitor_sp_swarm_retain -1 -ts_monitor_sp_swarm_phase 0 -draw_size 500,500 7 8 and we probably want it to run fast and not check convergence 9 10 -convest_num_refine 0 -ts_dt 0.01 -ts_max_steps 100 -ts_max_time 100 -output_step 10 11 */ 12 13 #include <petscts.h> 14 #include <petscdmplex.h> 15 #include <petscdmswarm.h> 16 #include <petsc/private/dmpleximpl.h> /* For norm and dot */ 17 18 PETSC_EXTERN PetscErrorCode circleSingleX(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *); 19 PETSC_EXTERN PetscErrorCode circleSingleV(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *); 20 PETSC_EXTERN PetscErrorCode circleMultipleX(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *); 21 PETSC_EXTERN PetscErrorCode circleMultipleV(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *); 22 23 typedef struct { 24 PetscBool error; /* Flag for printing the error */ 25 PetscInt ostep; /* print the energy at each ostep time steps */ 26 } AppCtx; 27 28 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 29 { 30 PetscErrorCode ierr; 31 32 PetscFunctionBeginUser; 33 options->error = PETSC_FALSE; 34 options->ostep = 100; 35 36 ierr = PetscOptionsBegin(comm, "", "Central Orbit Options", "DMSWARM");PetscCall(ierr); 37 PetscCall(PetscOptionsBool("-error", "Flag to print the error", "ex5.c", options->error, &options->error, NULL)); 38 PetscCall(PetscOptionsInt("-output_step", "Number of time steps between output", "ex5.c", options->ostep, &options->ostep, PETSC_NULL)); 39 ierr = PetscOptionsEnd();PetscCall(ierr); 40 PetscFunctionReturn(0); 41 } 42 43 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 44 { 45 PetscFunctionBeginUser; 46 PetscCall(DMCreate(comm, dm)); 47 PetscCall(DMSetType(*dm, DMPLEX)); 48 PetscCall(DMSetFromOptions(*dm)); 49 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 50 PetscFunctionReturn(0); 51 } 52 53 static PetscErrorCode CreateSwarm(DM dm, AppCtx *user, DM *sw) 54 { 55 PetscReal v0[1] = {1.}; 56 PetscInt dim; 57 58 PetscFunctionBeginUser; 59 PetscCall(DMGetDimension(dm, &dim)); 60 PetscCall(DMCreate(PetscObjectComm((PetscObject) dm), sw)); 61 PetscCall(DMSetType(*sw, DMSWARM)); 62 PetscCall(DMSetDimension(*sw, dim)); 63 PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC)); 64 PetscCall(DMSwarmSetCellDM(*sw, dm)); 65 PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR)); 66 PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "velocity", dim, PETSC_REAL)); 67 PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "species", 1, PETSC_INT)); 68 PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "initCoordinates", dim, PETSC_REAL)); 69 PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "initVelocity", dim, PETSC_REAL)); 70 PetscCall(DMSwarmFinalizeFieldRegister(*sw)); 71 PetscCall(DMSwarmComputeLocalSizeFromOptions(*sw)); 72 PetscCall(DMSwarmInitializeCoordinates(*sw)); 73 PetscCall(DMSwarmInitializeVelocitiesFromOptions(*sw, v0)); 74 PetscCall(DMSetFromOptions(*sw)); 75 PetscCall(DMSetApplicationContext(*sw, user)); 76 PetscCall(PetscObjectSetName((PetscObject) *sw, "Particles")); 77 PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view")); 78 { 79 Vec gc, gc0, gv, gv0; 80 81 PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, DMSwarmPICField_coor, &gc)); 82 PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, "initCoordinates", &gc0)); 83 PetscCall(VecCopy(gc, gc0)); 84 PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, DMSwarmPICField_coor, &gc)); 85 PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, "initCoordinates", &gc0)); 86 PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, "velocity", &gv)); 87 PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, "initVelocity", &gv0)); 88 PetscCall(VecCopy(gv, gv0)); 89 PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, "velocity", &gv)); 90 PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, "initVelocity", &gv0)); 91 } 92 PetscFunctionReturn(0); 93 } 94 95 static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec G, void *ctx) 96 { 97 DM sw; 98 const PetscReal *coords, *vel; 99 const PetscScalar *u; 100 PetscScalar *g; 101 PetscInt dim, d, Np, p; 102 103 PetscFunctionBeginUser; 104 PetscCall(TSGetDM(ts, &sw)); 105 PetscCall(DMGetDimension(sw, &dim)); 106 PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **) &coords)); 107 PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **) &vel)); 108 PetscCall(VecGetLocalSize(U, &Np)); 109 PetscCall(VecGetArrayRead(U, &u)); 110 PetscCall(VecGetArray(G, &g)); 111 Np /= 2*dim; 112 for (p = 0; p < Np; ++p) { 113 const PetscReal x0 = coords[p*dim+0]; 114 const PetscReal vy0 = vel[p*dim+1]; 115 const PetscReal omega = vy0/x0; 116 117 for (d = 0; d < dim; ++d) { 118 g[(p*2+0)*dim + d] = u[(p*2+1)*dim + d]; 119 g[(p*2+1)*dim + d] = -PetscSqr(omega) * u[(p*2+0)*dim + d]; 120 } 121 } 122 PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **) &coords)); 123 PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **) &vel)); 124 PetscCall(VecRestoreArrayRead(U, &u)); 125 PetscCall(VecRestoreArray(G, &g)); 126 PetscFunctionReturn(0); 127 } 128 129 /* J_{ij} = dF_i/dx_j 130 J_p = ( 0 1) 131 (-w^2 0) 132 */ 133 static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U , Mat J, Mat P, void *ctx) 134 { 135 DM sw; 136 const PetscReal *coords, *vel; 137 PetscInt dim, d, Np, p, rStart; 138 139 PetscFunctionBeginUser; 140 PetscCall(TSGetDM(ts, &sw)); 141 PetscCall(DMGetDimension(sw, &dim)); 142 PetscCall(VecGetLocalSize(U, &Np)); 143 PetscCall(MatGetOwnershipRange(J, &rStart, NULL)); 144 PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **) &coords)); 145 PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **) &vel)); 146 Np /= 2*dim; 147 for (p = 0; p < Np; ++p) { 148 const PetscReal x0 = coords[p*dim+0]; 149 const PetscReal vy0 = vel[p*dim+1]; 150 const PetscReal omega = vy0/x0; 151 PetscScalar vals[4] = {0., 1., -PetscSqr(omega), 0.}; 152 153 for (d = 0; d < dim; ++d) { 154 const PetscInt rows[2] = {(p*2+0)*dim+d + rStart, (p*2+1)*dim+d + rStart}; 155 PetscCall(MatSetValues(J, 2, rows, 2, rows, vals, INSERT_VALUES)); 156 } 157 } 158 PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **) &coords)); 159 PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **) &vel)); 160 PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 161 PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 162 PetscFunctionReturn(0); 163 } 164 165 static PetscErrorCode RHSFunctionX(TS ts, PetscReal t, Vec V, Vec Xres, void *ctx) 166 { 167 const PetscScalar *v; 168 PetscScalar *xres; 169 PetscInt Np, p; 170 171 PetscFunctionBeginUser; 172 PetscCall(VecGetLocalSize(Xres, &Np)); 173 PetscCall(VecGetArrayRead(V, &v)); 174 PetscCall(VecGetArray(Xres, &xres)); 175 for (p = 0; p < Np; ++p) xres[p] = v[p]; 176 PetscCall(VecRestoreArrayRead(V, &v)); 177 PetscCall(VecRestoreArray(Xres, &xres)); 178 PetscFunctionReturn(0); 179 } 180 181 static PetscErrorCode RHSFunctionV(TS ts, PetscReal t, Vec X, Vec Vres, void *ctx) 182 { 183 DM sw; 184 const PetscScalar *x; 185 const PetscReal *coords, *vel; 186 PetscScalar *vres; 187 PetscInt Np, p, dim, d; 188 189 PetscFunctionBeginUser; 190 PetscCall(TSGetDM(ts, &sw)); 191 PetscCall(DMGetDimension(sw, &dim)); 192 PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **) &coords)); 193 PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **) &vel)); 194 PetscCall(VecGetLocalSize(Vres, &Np)); 195 PetscCall(VecGetArrayRead(X, &x)); 196 PetscCall(VecGetArray(Vres, &vres)); 197 PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Dimension must be 2"); 198 Np /= dim; 199 for (p = 0; p < Np; ++p) { 200 const PetscReal x0 = coords[p*dim+0]; 201 const PetscReal vy0 = vel[p*dim+1]; 202 const PetscReal omega = vy0/x0; 203 204 for (d = 0; d < dim; ++d) vres[p*dim + d] = -PetscSqr(omega)*x[p*dim + d]; 205 } 206 PetscCall(VecRestoreArrayRead(X, &x)); 207 PetscCall(VecRestoreArray(Vres, &vres)); 208 PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **) &coords)); 209 PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **) &vel)); 210 PetscFunctionReturn(0); 211 } 212 213 static PetscErrorCode CreateSolution(TS ts) 214 { 215 DM sw; 216 Vec u; 217 PetscInt dim, Np; 218 219 PetscFunctionBegin; 220 PetscCall(TSGetDM(ts, &sw)); 221 PetscCall(DMGetDimension(sw, &dim)); 222 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 223 PetscCall(VecCreate(PETSC_COMM_WORLD, &u)); 224 PetscCall(VecSetBlockSize(u, dim)); 225 PetscCall(VecSetSizes(u, 2*Np*dim, PETSC_DECIDE)); 226 PetscCall(VecSetUp(u)); 227 PetscCall(TSSetSolution(ts, u)); 228 PetscCall(VecDestroy(&u)); 229 PetscFunctionReturn(0); 230 } 231 232 static PetscErrorCode SetProblem(TS ts) 233 { 234 AppCtx *user; 235 DM sw; 236 237 PetscFunctionBegin; 238 PetscCall(TSGetDM(ts, &sw)); 239 PetscCall(DMGetApplicationContext(sw, (void **) &user)); 240 // Define unified system for (X, V) 241 { 242 Mat J; 243 PetscInt dim, Np; 244 245 PetscCall(DMGetDimension(sw, &dim)); 246 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 247 PetscCall(MatCreate(PETSC_COMM_WORLD, &J)); 248 PetscCall(MatSetSizes(J, 2*Np*dim, 2*Np*dim, PETSC_DECIDE, PETSC_DECIDE)); 249 PetscCall(MatSetBlockSize(J, 2*dim)); 250 PetscCall(MatSetFromOptions(J)); 251 PetscCall(MatSetUp(J)); 252 PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &user)); 253 PetscCall(TSSetRHSJacobian(ts, J, J, RHSJacobian, &user)); 254 PetscCall(MatDestroy(&J)); 255 } 256 // Define split system for X and V 257 { 258 Vec u; 259 IS isx, isv, istmp; 260 const PetscInt *idx; 261 PetscInt dim, Np, rstart; 262 263 PetscCall(TSGetSolution(ts, &u)); 264 PetscCall(DMGetDimension(sw, &dim)); 265 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 266 PetscCall(VecGetOwnershipRange(u, &rstart, NULL)); 267 PetscCall(ISCreateStride(PETSC_COMM_WORLD, Np, (rstart/dim)+0, 2, &istmp)); 268 PetscCall(ISGetIndices(istmp, &idx)); 269 PetscCall(ISCreateBlock(PETSC_COMM_WORLD, dim, Np, idx, PETSC_COPY_VALUES, &isx)); 270 PetscCall(ISRestoreIndices(istmp, &idx)); 271 PetscCall(ISDestroy(&istmp)); 272 PetscCall(ISCreateStride(PETSC_COMM_WORLD, Np, (rstart/dim)+1, 2, &istmp)); 273 PetscCall(ISGetIndices(istmp, &idx)); 274 PetscCall(ISCreateBlock(PETSC_COMM_WORLD, dim, Np, idx, PETSC_COPY_VALUES, &isv)); 275 PetscCall(ISRestoreIndices(istmp, &idx)); 276 PetscCall(ISDestroy(&istmp)); 277 PetscCall(TSRHSSplitSetIS(ts, "position", isx)); 278 PetscCall(TSRHSSplitSetIS(ts, "momentum", isv)); 279 PetscCall(ISDestroy(&isx)); 280 PetscCall(ISDestroy(&isv)); 281 PetscCall(TSRHSSplitSetRHSFunction(ts, "position", NULL, RHSFunctionX, &user)); 282 PetscCall(TSRHSSplitSetRHSFunction(ts, "momentum", NULL, RHSFunctionV, &user)); 283 } 284 // Define symplectic formulation U_t = S . G, where G = grad F 285 { 286 //PetscCall(TSDiscGradSetFormulation(ts, RHSJacobianS, RHSObjectiveF, RHSFunctionG, &user)); 287 } 288 PetscFunctionReturn(0); 289 } 290 291 static PetscErrorCode DMSwarmTSRedistribute(TS ts) 292 { 293 DM sw; 294 Vec u; 295 PetscReal t, maxt, dt; 296 PetscInt n, maxn; 297 298 PetscFunctionBegin; 299 PetscCall(TSGetDM(ts, &sw)); 300 PetscCall(TSGetTime(ts, &t)); 301 PetscCall(TSGetMaxTime(ts, &maxt)); 302 PetscCall(TSGetTimeStep(ts, &dt)); 303 PetscCall(TSGetStepNumber(ts, &n)); 304 PetscCall(TSGetMaxSteps(ts, &maxn)); 305 306 PetscCall(TSReset(ts)); 307 PetscCall(TSSetDM(ts, sw)); 308 /* TODO Check whether TS was set from options */ 309 PetscCall(TSSetFromOptions(ts)); 310 PetscCall(TSSetTime(ts, t)); 311 PetscCall(TSSetMaxTime(ts, maxt)); 312 PetscCall(TSSetTimeStep(ts, dt)); 313 PetscCall(TSSetStepNumber(ts, n)); 314 PetscCall(TSSetMaxSteps(ts, maxn)); 315 316 PetscCall(CreateSolution(ts)); 317 PetscCall(SetProblem(ts)); 318 PetscCall(TSGetSolution(ts, &u)); 319 PetscFunctionReturn(0); 320 } 321 322 PetscErrorCode circleSingleX(PetscInt dim, PetscReal time, const PetscReal dummy[], PetscInt p, PetscScalar x[], void *ctx) 323 { 324 x[0] = p+1.; 325 x[1] = 0.; 326 return 0; 327 } 328 329 PetscErrorCode circleSingleV(PetscInt dim, PetscReal time, const PetscReal dummy[], PetscInt p, PetscScalar v[], void *ctx) 330 { 331 v[0] = 0.; 332 v[1] = PetscSqrtReal(1000. / (p+1.)); 333 return 0; 334 } 335 336 /* Put 5 particles into each circle */ 337 PetscErrorCode circleMultipleX(PetscInt dim, PetscReal time, const PetscReal dummy[], PetscInt p, PetscScalar x[], void *ctx) 338 { 339 const PetscInt n = 5; 340 const PetscReal r0 = (p / n) + 1.; 341 const PetscReal th0 = (2.*PETSC_PI * (p%n))/n; 342 343 x[0] = r0 * PetscCosReal(th0); 344 x[1] = r0 * PetscSinReal(th0); 345 return 0; 346 } 347 348 /* Put 5 particles into each circle */ 349 PetscErrorCode circleMultipleV(PetscInt dim, PetscReal time, const PetscReal dummy[], PetscInt p, PetscScalar v[], void *ctx) 350 { 351 const PetscInt n = 5; 352 const PetscReal r0 = (p / n) + 1.; 353 const PetscReal th0 = (2.*PETSC_PI * (p%n))/n; 354 const PetscReal omega = PetscSqrtReal(1000. / r0) / r0; 355 356 v[0] = -r0 * omega * PetscSinReal(th0); 357 v[1] = r0 * omega * PetscCosReal(th0); 358 return 0; 359 } 360 361 /* 362 InitializeSolveAndSwarm - Set the solution values to the swarm coordinates and velocities, and also possibly set the initial values. 363 364 Input Parameters: 365 + ts - The TS 366 - useInitial - Flag to also set the initial conditions to the current coodinates and velocities and setup the problem 367 368 Output Parameters: 369 . u - The initialized solution vector 370 371 Level: advanced 372 373 .seealso: InitializeSolve() 374 */ 375 static PetscErrorCode InitializeSolveAndSwarm(TS ts, PetscBool useInitial) 376 { 377 DM sw; 378 Vec u, gc, gv, gc0, gv0; 379 IS isx, isv; 380 AppCtx *user; 381 382 PetscFunctionBeginUser; 383 PetscCall(TSGetDM(ts, &sw)); 384 PetscCall(DMGetApplicationContext(sw, &user)); 385 if (useInitial) { 386 PetscReal v0[1] = {1.}; 387 388 PetscCall(DMSwarmInitializeCoordinates(sw)); 389 PetscCall(DMSwarmInitializeVelocitiesFromOptions(sw, v0)); 390 PetscCall(DMSwarmMigrate(sw, PETSC_TRUE)); 391 PetscCall(DMSwarmTSRedistribute(ts)); 392 } 393 PetscCall(TSGetSolution(ts, &u)); 394 PetscCall(TSRHSSplitGetIS(ts, "position", &isx)); 395 PetscCall(TSRHSSplitGetIS(ts, "momentum", &isv)); 396 PetscCall(DMSwarmCreateGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc)); 397 PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "initCoordinates", &gc0)); 398 if (useInitial) PetscCall(VecCopy(gc, gc0)); 399 PetscCall(VecISCopy(u, isx, SCATTER_FORWARD, gc)); 400 PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc)); 401 PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "initCoordinates", &gc0)); 402 PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "velocity", &gv)); 403 PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "initVelocity", &gv0)); 404 if (useInitial) PetscCall(VecCopy(gv, gv0)); 405 PetscCall(VecISCopy(u, isv, SCATTER_FORWARD, gv)); 406 PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "velocity", &gv)); 407 PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "initVelocity", &gv0)); 408 PetscFunctionReturn(0); 409 } 410 411 static PetscErrorCode InitializeSolve(TS ts, Vec u) 412 { 413 PetscFunctionBegin; 414 PetscCall(TSSetSolution(ts, u)); 415 PetscCall(InitializeSolveAndSwarm(ts, PETSC_TRUE)); 416 PetscFunctionReturn(0); 417 } 418 419 static PetscErrorCode ComputeError(TS ts, Vec U, Vec E) 420 { 421 MPI_Comm comm; 422 DM sw; 423 AppCtx *user; 424 const PetscScalar *u; 425 const PetscReal *coords, *vel; 426 PetscScalar *e; 427 PetscReal t; 428 PetscInt dim, Np, p; 429 430 PetscFunctionBeginUser; 431 PetscCall(PetscObjectGetComm((PetscObject) ts, &comm)); 432 PetscCall(TSGetDM(ts, &sw)); 433 PetscCall(DMGetApplicationContext(sw, &user)); 434 PetscCall(DMGetDimension(sw, &dim)); 435 PetscCall(TSGetSolveTime(ts, &t)); 436 PetscCall(VecGetArray(E, &e)); 437 PetscCall(VecGetArrayRead(U, &u)); 438 PetscCall(VecGetLocalSize(U, &Np)); 439 PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **) &coords)); 440 PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **) &vel)); 441 Np /= 2*dim; 442 for (p = 0; p < Np; ++p) { 443 /* TODO generalize initial conditions and project into plane instead of assuming x-y */ 444 const PetscReal r0 = DMPlex_NormD_Internal(dim, &coords[p*dim]); 445 const PetscReal th0 = PetscAtan2Real(coords[p*dim+1], coords[p*dim+0]); 446 const PetscReal v0 = DMPlex_NormD_Internal(dim, &vel[p*dim]); 447 const PetscReal omega = v0/r0; 448 const PetscReal ct = PetscCosReal(omega*t + th0); 449 const PetscReal st = PetscSinReal(omega*t + th0); 450 const PetscScalar *x = &u[(p*2+0)*dim]; 451 const PetscScalar *v = &u[(p*2+1)*dim]; 452 const PetscReal xe[3] = { r0*ct, r0*st, 0.0}; 453 const PetscReal ve[3] = {-v0*st, v0*ct, 0.0}; 454 PetscInt d; 455 456 for (d = 0; d < dim; ++d) { 457 e[(p*2+0)*dim+d] = x[d] - xe[d]; 458 e[(p*2+1)*dim+d] = v[d] - ve[d]; 459 } 460 if (user->error) { 461 const PetscReal en = 0.5*DMPlex_DotRealD_Internal(dim, v, v); 462 const PetscReal exen = 0.5*PetscSqr(v0); 463 PetscCall(PetscPrintf(comm, "t %.4g: p%D error [%.2g %.2g] sol [(%.6lf %.6lf) (%.6lf %.6lf)] exact [(%.6lf %.6lf) (%.6lf %.6lf)] energy/exact energy %g / %g (%.10lf%%)\n", t, p, (double) DMPlex_NormD_Internal(dim, &e[(p*2+0)*dim]), (double) DMPlex_NormD_Internal(dim, &e[(p*2+1)*dim]), (double) x[0], (double) x[1], (double) v[0], (double) v[1], (double) xe[0], (double) xe[1], (double) ve[0], (double) ve[1], (double) en, (double) exen, (double) PetscAbsReal(exen - en)*100./exen)); 464 } 465 } 466 PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **) &coords)); 467 PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **) &vel)); 468 PetscCall(VecRestoreArrayRead(U, &u)); 469 PetscCall(VecRestoreArray(E, &e)); 470 PetscFunctionReturn(0); 471 } 472 473 static PetscErrorCode EnergyMonitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) 474 { 475 const PetscInt ostep = ((AppCtx *) ctx)->ostep; 476 DM sw; 477 const PetscScalar *u; 478 PetscInt dim, Np, p; 479 480 PetscFunctionBeginUser; 481 if (step%ostep == 0) { 482 PetscCall(TSGetDM(ts, &sw)); 483 PetscCall(DMGetDimension(sw, &dim)); 484 PetscCall(VecGetArrayRead(U, &u)); 485 PetscCall(VecGetLocalSize(U, &Np)); 486 Np /= 2*dim; 487 if (!step) PetscCall(PetscPrintf(PetscObjectComm((PetscObject) ts), "Time Step Part Energy\n")); 488 for (p = 0; p < Np; ++p) { 489 const PetscReal v2 = DMPlex_DotRealD_Internal(dim, &u[(p*2+1)*dim], &u[(p*2+1)*dim]); 490 491 PetscCall(PetscSynchronizedPrintf(PetscObjectComm((PetscObject) ts), "%.6lf %4D %4D %10.4lf\n", t, step, p, (double) 0.5*v2)); 492 } 493 PetscCall(PetscSynchronizedFlush(PetscObjectComm((PetscObject) ts), NULL)); 494 PetscCall(VecRestoreArrayRead(U, &u)); 495 } 496 PetscFunctionReturn(0); 497 } 498 499 static PetscErrorCode MigrateParticles(TS ts) 500 { 501 DM sw; 502 503 PetscFunctionBeginUser; 504 PetscCall(TSGetDM(ts, &sw)); 505 PetscCall(DMViewFromOptions(sw, NULL, "-migrate_view_pre")); 506 { 507 Vec u, gc, gv; 508 IS isx, isv; 509 510 PetscCall(TSGetSolution(ts, &u)); 511 PetscCall(TSRHSSplitGetIS(ts, "position", &isx)); 512 PetscCall(TSRHSSplitGetIS(ts, "momentum", &isv)); 513 PetscCall(DMSwarmCreateGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc)); 514 PetscCall(VecISCopy(u, isx, SCATTER_REVERSE, gc)); 515 PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc)); 516 PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "velocity", &gv)); 517 PetscCall(VecISCopy(u, isv, SCATTER_REVERSE, gv)); 518 PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "velocity", &gv)); 519 } 520 PetscCall(DMSwarmMigrate(sw, PETSC_TRUE)); 521 PetscCall(DMSwarmTSRedistribute(ts)); 522 PetscCall(InitializeSolveAndSwarm(ts, PETSC_FALSE)); 523 //PetscCall(VecViewFromOptions(u, NULL, "-sol_migrate_view")); 524 PetscFunctionReturn(0); 525 } 526 527 int main(int argc, char **argv) 528 { 529 DM dm, sw; 530 TS ts; 531 Vec u; 532 AppCtx user; 533 534 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 535 PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 536 PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 537 PetscCall(CreateSwarm(dm, &user, &sw)); 538 PetscCall(DMSetApplicationContext(sw, &user)); 539 540 PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 541 PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 542 PetscCall(TSSetDM(ts, sw)); 543 PetscCall(TSSetMaxTime(ts, 0.1)); 544 PetscCall(TSSetTimeStep(ts, 0.00001)); 545 PetscCall(TSSetMaxSteps(ts, 100)); 546 PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 547 PetscCall(TSMonitorSet(ts, EnergyMonitor, &user, NULL)); 548 PetscCall(TSSetFromOptions(ts)); 549 PetscCall(TSSetComputeInitialCondition(ts, InitializeSolve)); 550 PetscCall(TSSetComputeExactError(ts, ComputeError)); 551 PetscCall(TSSetPostStep(ts, MigrateParticles)); 552 553 PetscCall(CreateSolution(ts)); 554 PetscCall(TSGetSolution(ts, &u)); 555 PetscCall(TSComputeInitialCondition(ts, u)); 556 PetscCall(TSSolve(ts, NULL)); 557 558 PetscCall(TSDestroy(&ts)); 559 PetscCall(DMDestroy(&sw)); 560 PetscCall(DMDestroy(&dm)); 561 PetscCall(PetscFinalize()); 562 return 0; 563 } 564 565 /*TEST 566 567 build: 568 requires: double !complex 569 570 testset: 571 requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT) 572 args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -5,-5 -dm_plex_box_upper 5,5 \ 573 -dm_swarm_num_particles 2 -dm_swarm_coordinate_function circleSingleX -dm_swarm_velocity_function circleSingleV \ 574 -ts_type basicsymplectic -ts_convergence_estimate -convest_num_refine 2 \ 575 -dm_view -output_step 50 -error 576 test: 577 suffix: bsi_2d_1 578 args: -ts_basicsymplectic_type 1 579 test: 580 suffix: bsi_2d_2 581 args: -ts_basicsymplectic_type 2 582 test: 583 suffix: bsi_2d_3 584 args: -ts_basicsymplectic_type 3 585 test: 586 suffix: bsi_2d_4 587 args: -ts_basicsymplectic_type 4 -ts_dt 0.0001 588 589 testset: 590 requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT) 591 args: -dm_swarm_num_particles 2 -dm_swarm_coordinate_function circleSingleX -dm_swarm_velocity_function circleSingleV \ 592 -ts_type theta -ts_theta_theta 0.5 -ts_convergence_estimate -convest_num_refine 2 \ 593 -mat_type baij -ksp_error_if_not_converged -pc_type lu \ 594 -dm_view -output_step 50 -error 595 test: 596 suffix: im_2d_0 597 args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -5,-5 -dm_plex_box_upper 5,5 598 599 testset: 600 requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT) 601 args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 10,10 -dm_plex_box_lower -5,-5 -dm_plex_box_upper 5,5 -petscpartitioner_type simple \ 602 -dm_swarm_num_particles 2 -dm_swarm_coordinate_function circleSingleX -dm_swarm_velocity_function circleSingleV \ 603 -ts_type basicsymplectic -ts_convergence_estimate -convest_num_refine 2 \ 604 -dm_view -output_step 50 605 test: 606 suffix: bsi_2d_mesh_1 607 args: -ts_basicsymplectic_type 1 -error 608 test: 609 suffix: bsi_2d_mesh_1_par_2 610 nsize: 2 611 args: -ts_basicsymplectic_type 1 612 test: 613 suffix: bsi_2d_mesh_1_par_3 614 nsize: 3 615 args: -ts_basicsymplectic_type 1 616 test: 617 suffix: bsi_2d_mesh_1_par_4 618 nsize: 4 619 args: -ts_basicsymplectic_type 1 -dm_swarm_num_particles 0,0,2,0 620 621 testset: 622 requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT) 623 args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 10,10 -dm_plex_box_lower -5,-5 -dm_plex_box_upper 5,5 \ 624 -dm_swarm_num_particles 10 -dm_swarm_coordinate_function circleMultipleX -dm_swarm_velocity_function circleMultipleV \ 625 -ts_convergence_estimate -convest_num_refine 2 \ 626 -dm_view -output_step 50 -error 627 test: 628 suffix: bsi_2d_multiple_1 629 args: -ts_type basicsymplectic -ts_basicsymplectic_type 1 630 test: 631 suffix: bsi_2d_multiple_2 632 args: -ts_type basicsymplectic -ts_basicsymplectic_type 2 633 test: 634 suffix: bsi_2d_multiple_3 635 args: -ts_type basicsymplectic -ts_basicsymplectic_type 3 -ts_dt 0.001 636 test: 637 suffix: im_2d_multiple_0 638 args: -ts_type theta -ts_theta_theta 0.5 \ 639 -mat_type baij -ksp_error_if_not_converged -pc_type lu 640 641 TEST*/ 642