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