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