1c4762a1bSJed Brown static char help[] = "Vlasov example of particles orbiting around a central massive point.\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown #include <petsc/private/dmpleximpl.h> /* For norm */ 5c4762a1bSJed Brown #include <petsc/private/petscfeimpl.h> /* For CoordinatesRefToReal() */ 6c4762a1bSJed Brown #include <petscdmswarm.h> 7c4762a1bSJed Brown #include <petscts.h> 8c4762a1bSJed Brown 9c4762a1bSJed Brown typedef struct { 1030602db0SMatthew G. Knepley PetscInt dim; 11c4762a1bSJed Brown PetscInt particlesPerCell; /* The number of partices per cell */ 12c4762a1bSJed Brown PetscReal momentTol; /* Tolerance for checking moment conservation */ 13c4762a1bSJed Brown PetscBool monitor; /* Flag for using the TS monitor */ 14c4762a1bSJed Brown PetscBool error; /* Flag for printing the error */ 15c4762a1bSJed Brown PetscInt ostep; /* print the energy at each ostep time steps */ 16c4762a1bSJed Brown } AppCtx; 17c4762a1bSJed Brown 18c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 19c4762a1bSJed Brown { 20c4762a1bSJed Brown PetscErrorCode ierr; 21c4762a1bSJed Brown 22c4762a1bSJed Brown PetscFunctionBeginUser; 23c4762a1bSJed Brown options->monitor = PETSC_FALSE; 24c4762a1bSJed Brown options->error = PETSC_FALSE; 25c4762a1bSJed Brown options->particlesPerCell = 1; 26c4762a1bSJed Brown options->momentTol = 100.0*PETSC_MACHINE_EPSILON; 27c4762a1bSJed Brown options->ostep = 100; 28c4762a1bSJed Brown 29c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Vlasov Options", "DMPLEX");CHKERRQ(ierr); 305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-output_step", "Number of time steps between output", "ex4.c", options->ostep, &options->ostep, PETSC_NULL)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-monitor", "Flag to use the TS monitor", "ex4.c", options->monitor, &options->monitor, NULL)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-error", "Flag to print the error", "ex4.c", options->error, &options->error, NULL)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-particles_per_cell", "Number of particles per cell", "ex4.c", options->particlesPerCell, &options->particlesPerCell, NULL)); 34c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 35c4762a1bSJed Brown 36c4762a1bSJed Brown PetscFunctionReturn(0); 37c4762a1bSJed Brown } 38c4762a1bSJed Brown 39c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user) 40c4762a1bSJed Brown { 41c4762a1bSJed Brown PetscFunctionBeginUser; 425f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 445f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 465f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(*dm, &user->dim)); 47c4762a1bSJed Brown PetscFunctionReturn(0); 48c4762a1bSJed Brown } 49c4762a1bSJed Brown 50c4762a1bSJed Brown static PetscErrorCode SetInitialCoordinates(DM dmSw) 51c4762a1bSJed Brown { 52c4762a1bSJed Brown DM dm; 53c4762a1bSJed Brown AppCtx *user; 54c4762a1bSJed Brown PetscRandom rnd; 55c4762a1bSJed Brown PetscBool simplex; 56c4762a1bSJed Brown PetscReal *centroid, *coords, *xi0, *v0, *J, *invJ, detJ; 57c4762a1bSJed Brown PetscInt dim, d, cStart, cEnd, c, Np, p; 58c4762a1bSJed Brown 59c4762a1bSJed Brown PetscFunctionBeginUser; 605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PetscObjectComm((PetscObject) dmSw), &rnd)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetInterval(rnd, -1.0, 1.0)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rnd)); 63c4762a1bSJed Brown 645f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetApplicationContext(dmSw, &user)); 65c4762a1bSJed Brown Np = user->particlesPerCell; 665f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetCellDM(dmSw, &dm)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexIsSimplex(dm, &simplex)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc5(dim, ¢roid, dim, &xi0, dim, &v0, dim*dim, &J, dim*dim, &invJ)); 71c4762a1bSJed Brown for (d = 0; d < dim; ++d) xi0[d] = -1.0; 725f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords)); 73c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 74c4762a1bSJed Brown if (Np == 1) { 755f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL)); 76c4762a1bSJed Brown for (d = 0; d < dim; ++d) coords[c*dim+d] = centroid[d]; 77c4762a1bSJed Brown } else { 785f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); /* affine */ 79c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 80c4762a1bSJed Brown const PetscInt n = c*Np + p; 81c4762a1bSJed Brown PetscReal sum = 0.0, refcoords[3]; 82c4762a1bSJed Brown 83c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValueReal(rnd, &refcoords[d])); 85c4762a1bSJed Brown sum += refcoords[d]; 86c4762a1bSJed Brown } 87c4762a1bSJed Brown if (simplex && sum > 0.0) for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim)*sum; 88c4762a1bSJed Brown CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n*dim]); 89c4762a1bSJed Brown } 90c4762a1bSJed Brown } 91c4762a1bSJed Brown } 925f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree5(centroid, xi0, v0, J, invJ)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rnd)); 95c4762a1bSJed Brown PetscFunctionReturn(0); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown 98c4762a1bSJed Brown static PetscErrorCode SetInitialConditions(DM dmSw, Vec u) 99c4762a1bSJed Brown { 100c4762a1bSJed Brown DM dm; 101c4762a1bSJed Brown AppCtx *user; 102c4762a1bSJed Brown PetscScalar *initialConditions; 103c4762a1bSJed Brown PetscInt dim, cStart, cEnd, c, Np, p; 104c4762a1bSJed Brown 105c4762a1bSJed Brown PetscFunctionBeginUser; 1065f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetApplicationContext(dmSw, &user)); 107c4762a1bSJed Brown Np = user->particlesPerCell; 1085f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetCellDM(dmSw, &dm)); 1095f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(u, &initialConditions)); 112c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 113c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 114c4762a1bSJed Brown const PetscInt n = c*Np + p; 115c4762a1bSJed Brown 116c4762a1bSJed Brown initialConditions[(n*2 + 0)*dim + 0] = n+1; 117c4762a1bSJed Brown initialConditions[(n*2 + 0)*dim + 1] = 0; 118c4762a1bSJed Brown initialConditions[(n*2 + 1)*dim + 0] = 0; 119c4762a1bSJed Brown initialConditions[(n*2 + 1)*dim + 1] = PetscSqrtReal(1000./(n+1.)); 120c4762a1bSJed Brown } 121c4762a1bSJed Brown } 1225f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(u, &initialConditions)); 123c4762a1bSJed Brown PetscFunctionReturn(0); 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126c4762a1bSJed Brown static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user) 127c4762a1bSJed Brown { 128c4762a1bSJed Brown PetscInt *cellid; 129c4762a1bSJed Brown PetscInt dim, cStart, cEnd, c, Np = user->particlesPerCell, p; 130c4762a1bSJed Brown 131c4762a1bSJed Brown PetscFunctionBeginUser; 1325f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 1335f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(PetscObjectComm((PetscObject) dm), sw)); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*sw, DMSWARM)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetDimension(*sw, dim)); 136c4762a1bSJed Brown 1375f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetType(*sw, DMSWARM_PIC)); 1385f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetCellDM(*sw, dm)); 1395f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRegisterPetscDatatypeField(*sw, "kinematics", 2*dim, PETSC_REAL)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmFinalizeFieldRegister(*sw)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0)); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*sw)); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid)); 145c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 146c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 147c4762a1bSJed Brown const PetscInt n = c*Np + p; 148c4762a1bSJed Brown 149c4762a1bSJed Brown cellid[n] = c; 150c4762a1bSJed Brown } 151c4762a1bSJed Brown } 1525f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid)); 1535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) *sw, "Particles")); 1545f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*sw, NULL, "-sw_view")); 155c4762a1bSJed Brown PetscFunctionReturn(0); 156c4762a1bSJed Brown } 157c4762a1bSJed Brown 158c4762a1bSJed Brown /* Create particle RHS Functions for gravity with G = 1 for simplification */ 159c4762a1bSJed Brown static PetscErrorCode RHSFunction1(TS ts, PetscReal t, Vec V, Vec Xres, void *ctx) 160c4762a1bSJed Brown { 161c4762a1bSJed Brown const PetscScalar *v; 162c4762a1bSJed Brown PetscScalar *xres; 163c4762a1bSJed Brown PetscInt Np, p, dim, d; 164c4762a1bSJed Brown 165c4762a1bSJed Brown PetscFunctionBeginUser; 16630602db0SMatthew G. Knepley /* The DM is not currently pushed down to the splits */ 16730602db0SMatthew G. Knepley dim = ((AppCtx *) ctx)->dim; 1685f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(Xres, &Np)); 169c4762a1bSJed Brown Np /= dim; 1705f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(Xres, &xres)); 1715f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(V, &v)); 172c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 173c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 174c4762a1bSJed Brown xres[p*dim+d] = v[p*dim+d]; 175c4762a1bSJed Brown } 176c4762a1bSJed Brown } 1775f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(V,& v)); 1785f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(Xres, &xres)); 179c4762a1bSJed Brown PetscFunctionReturn(0); 180c4762a1bSJed Brown } 181c4762a1bSJed Brown 182c4762a1bSJed Brown static PetscErrorCode RHSFunction2(TS ts, PetscReal t, Vec X, Vec Vres, void *ctx) 183c4762a1bSJed Brown { 184c4762a1bSJed Brown const PetscScalar *x; 185c4762a1bSJed Brown PetscScalar *vres; 186c4762a1bSJed Brown PetscInt Np, p, dim, d; 187c4762a1bSJed Brown 188c4762a1bSJed Brown PetscFunctionBeginUser; 18930602db0SMatthew G. Knepley /* The DM is not currently pushed down to the splits */ 19030602db0SMatthew G. Knepley dim = ((AppCtx *) ctx)->dim; 1915f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(Vres, &Np)); 192c4762a1bSJed Brown Np /= dim; 1935f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(Vres, &vres)); 1945f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(X, &x)); 195c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 196c4762a1bSJed Brown const PetscScalar rsqr = DMPlex_NormD_Internal(dim, &x[p*dim]); 197c4762a1bSJed Brown 198c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 199c4762a1bSJed Brown vres[p*dim+d] = -(1000./(p+1.)) * x[p*dim+d]/PetscSqr(rsqr); 200c4762a1bSJed Brown } 201c4762a1bSJed Brown } 2025f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(Vres, &vres)); 2035f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(X, &x)); 204c4762a1bSJed Brown PetscFunctionReturn(0); 205c4762a1bSJed Brown } 206c4762a1bSJed Brown 207c4762a1bSJed Brown static PetscErrorCode RHSFunctionParticles(TS ts, PetscReal t , Vec U, Vec R, void *ctx) 208c4762a1bSJed Brown { 209c4762a1bSJed Brown DM dm; 210c4762a1bSJed Brown const PetscScalar *u; 211c4762a1bSJed Brown PetscScalar *r; 212c4762a1bSJed Brown PetscInt Np, p, dim, d; 213c4762a1bSJed Brown 214c4762a1bSJed Brown PetscFunctionBeginUser; 2155f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts, &dm)); 2165f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 2175f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(U, &Np)); 218c4762a1bSJed Brown Np /= 2*dim; 2195f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(U, &u)); 2205f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(R, &r)); 221c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 222c4762a1bSJed Brown const PetscScalar rsqr = DMPlex_NormD_Internal(dim, &u[p*2*dim]); 223c4762a1bSJed Brown 224c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 225c4762a1bSJed Brown r[p*2*dim+d] = u[p*2*dim+d+2]; 226c4762a1bSJed Brown r[p*2*dim+d+2] = -(1000./(1.+p)) * u[p*2*dim+d]/PetscSqr(rsqr); 227c4762a1bSJed Brown } 228c4762a1bSJed Brown } 2295f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(U, &u)); 2305f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(R, &r)); 231c4762a1bSJed Brown PetscFunctionReturn(0); 232c4762a1bSJed Brown } 233c4762a1bSJed Brown 234c4762a1bSJed Brown static PetscErrorCode InitializeSolve(TS ts, Vec u) 235c4762a1bSJed Brown { 236c4762a1bSJed Brown DM dm; 237c4762a1bSJed Brown AppCtx *user; 238c4762a1bSJed Brown 239c4762a1bSJed Brown PetscFunctionBeginUser; 2405f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts, &dm)); 2415f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetApplicationContext(dm, &user)); 2425f80ce2aSJacob Faibussowitsch CHKERRQ(SetInitialCoordinates(dm)); 2435f80ce2aSJacob Faibussowitsch CHKERRQ(SetInitialConditions(dm, u)); 244c4762a1bSJed Brown PetscFunctionReturn(0); 245c4762a1bSJed Brown } 246c4762a1bSJed Brown 247c4762a1bSJed Brown static PetscErrorCode ComputeError(TS ts, Vec U, Vec E) 248c4762a1bSJed Brown { 249c4762a1bSJed Brown MPI_Comm comm; 250c4762a1bSJed Brown DM sdm; 251c4762a1bSJed Brown AppCtx *user; 252c4762a1bSJed Brown const PetscScalar *u, *coords; 253c4762a1bSJed Brown PetscScalar *e; 254c4762a1bSJed Brown PetscReal t; 255c4762a1bSJed Brown PetscInt dim, Np, p; 256c4762a1bSJed Brown 257c4762a1bSJed Brown PetscFunctionBeginUser; 2585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject) ts, &comm)); 2595f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts, &sdm)); 2605f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetApplicationContext(sdm, &user)); 2615f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(sdm, &dim)); 2625f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSolveTime(ts, &t)); 2635f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(E, &e)); 2645f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(U, &u)); 2655f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(U, &Np)); 266c4762a1bSJed Brown Np /= 2*dim; 2675f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords)); 268c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 269c4762a1bSJed Brown const PetscScalar *x = &u[(p*2+0)*dim]; 270c4762a1bSJed Brown const PetscScalar *v = &u[(p*2+1)*dim]; 271c4762a1bSJed Brown const PetscReal x0 = p+1.; 272c4762a1bSJed Brown const PetscReal omega = PetscSqrtReal(1000./(p+1.))/x0; 273c4762a1bSJed Brown const PetscReal xe[3] = { x0*PetscCosReal(omega*t), x0*PetscSinReal(omega*t), 0.0}; 274c4762a1bSJed Brown const PetscReal ve[3] = {-x0*omega*PetscSinReal(omega*t), x0*omega*PetscCosReal(omega*t), 0.0}; 275c4762a1bSJed Brown PetscInt d; 276c4762a1bSJed Brown 277c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 278c4762a1bSJed Brown e[(p*2+0)*dim+d] = x[d] - xe[d]; 279c4762a1bSJed Brown e[(p*2+1)*dim+d] = v[d] - ve[d]; 280c4762a1bSJed Brown } 2815f80ce2aSJacob Faibussowitsch if (user->error) CHKERRQ(PetscPrintf(comm, "t %.4g: p%D error [%.2g %.2g] sol [(%.6lf %.6lf) (%.6lf %.6lf)] exact [(%.6lf %.6lf) (%.6lf %.6lf)] energy/exact energy %g / %g\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], 0.5*DMPlex_NormD_Internal(dim, v), (double) (0.5*(1000./(p+1))))); 282c4762a1bSJed Brown } 2835f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords)); 2845f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(U, &u)); 2855f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(E, &e)); 286c4762a1bSJed Brown PetscFunctionReturn(0); 287c4762a1bSJed Brown } 288c4762a1bSJed Brown 289c4762a1bSJed Brown int main(int argc, char **argv) 290c4762a1bSJed Brown { 291c4762a1bSJed Brown TS ts; 292c4762a1bSJed Brown TSConvergedReason reason; 293c4762a1bSJed Brown DM dm, sw; 294c4762a1bSJed Brown Vec u; 295c4762a1bSJed Brown IS is1, is2; 296c4762a1bSJed Brown PetscInt *idx1, *idx2; 297c4762a1bSJed Brown MPI_Comm comm; 298c4762a1bSJed Brown AppCtx user; 299c4762a1bSJed Brown const PetscScalar *endVals; 300c4762a1bSJed Brown PetscReal ftime = .1; 301c4762a1bSJed Brown PetscInt locSize, dim, d, Np, p, steps; 302c4762a1bSJed Brown 303*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, NULL, help)); 304c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 3055f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(comm, &user)); 306c4762a1bSJed Brown 3075f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(comm, &dm, &user)); 3085f80ce2aSJacob Faibussowitsch CHKERRQ(CreateParticles(dm, &sw, &user)); 3095f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetApplicationContext(sw, &user)); 310c4762a1bSJed Brown 3115f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(comm, &ts)); 3125f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts, TSBASICSYMPLECTIC)); 3135f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetDM(ts, sw)); 3145f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts, ftime)); 3155f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts, 0.0001)); 3165f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxSteps(ts, 10)); 3175f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 3185f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTime(ts, 0.0)); 3195f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts, NULL, RHSFunctionParticles, &user)); 320c4762a1bSJed Brown 3215f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmCreateGlobalVectorFromField(sw, "kinematics", &u)); 3225f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(sw, &dim)); 3235f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(u, &locSize)); 324c4762a1bSJed Brown Np = locSize/(2*dim); 3255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(locSize/2, &idx1)); 3265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(locSize/2, &idx2)); 327c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 328c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 329c4762a1bSJed Brown idx1[p*dim+d] = (p*2+0)*dim + d; 330c4762a1bSJed Brown idx2[p*dim+d] = (p*2+1)*dim + d; 331c4762a1bSJed Brown } 332c4762a1bSJed Brown } 3335f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(comm, locSize/2, idx1, PETSC_OWN_POINTER, &is1)); 3345f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(comm, locSize/2, idx2, PETSC_OWN_POINTER, &is2)); 3355f80ce2aSJacob Faibussowitsch CHKERRQ(TSRHSSplitSetIS(ts, "position", is1)); 3365f80ce2aSJacob Faibussowitsch CHKERRQ(TSRHSSplitSetIS(ts, "momentum", is2)); 3375f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is1)); 3385f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is2)); 3395f80ce2aSJacob Faibussowitsch CHKERRQ(TSRHSSplitSetRHSFunction(ts,"position",NULL,RHSFunction1,&user)); 3405f80ce2aSJacob Faibussowitsch CHKERRQ(TSRHSSplitSetRHSFunction(ts,"momentum",NULL,RHSFunction2,&user)); 341c4762a1bSJed Brown 3425f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 3435f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetComputeInitialCondition(ts, InitializeSolve)); 3445f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetComputeExactError(ts, ComputeError)); 3455f80ce2aSJacob Faibussowitsch CHKERRQ(TSComputeInitialCondition(ts, u)); 3465f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(u, NULL, "-init_view")); 3475f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts, u)); 3485f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSolveTime(ts, &ftime)); 3495f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetConvergedReason(ts, &reason)); 3505f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetStepNumber(ts, &steps)); 3515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(comm,"%s at time %g after %D steps\n", TSConvergedReasons[reason], (double) ftime, steps)); 352c4762a1bSJed Brown 3535f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(u, &endVals)); 354c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 355c4762a1bSJed Brown const PetscReal norm = DMPlex_NormD_Internal(dim, &endVals[(p*2 + 1)*dim]); 3565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(comm, "Particle %D initial Energy: %g Final Energy: %g\n", p, (double) (0.5*(1000./(p+1))), (double) 0.5*PetscSqr(norm))); 357c4762a1bSJed Brown } 3585f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(u, &endVals)); 3595f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmDestroyGlobalVectorFromField(sw, "kinematics", &u)); 3605f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 3615f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&sw)); 3625f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 363*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 364*b122ec5aSJacob Faibussowitsch return 0; 365c4762a1bSJed Brown } 366c4762a1bSJed Brown 367c4762a1bSJed Brown /*TEST 368c4762a1bSJed Brown 369c4762a1bSJed Brown build: 370c4762a1bSJed Brown requires: triangle !single !complex 371c4762a1bSJed Brown test: 372c4762a1bSJed Brown suffix: bsi1 373d7462660SMatthew Knepley args: -dm_plex_box_faces 1,1 -dm_plex_box_lower -5,-5 -dm_plex_box_upper 5,5 \ 374d7462660SMatthew Knepley -ts_basicsymplectic_type 1 -ts_max_time 0.1 -ts_convergence_estimate -convest_num_refine 2 \ 375d7462660SMatthew Knepley -dm_view -sw_view -ts_monitor_sp_swarm -ts_monitor_sp_swarm_retain -ts_monitor_sp_swarm_phase 0 376c4762a1bSJed Brown test: 377c4762a1bSJed Brown suffix: bsi2 378d7462660SMatthew Knepley args: -dm_plex_box_faces 1,1 -dm_plex_box_lower -5,-5 -dm_plex_box_upper 5,5 \ 379d7462660SMatthew Knepley -ts_basicsymplectic_type 2 -ts_max_time 0.1 -ts_convergence_estimate -convest_num_refine 2 \ 380d7462660SMatthew Knepley -dm_view -sw_view -ts_monitor_sp_swarm -ts_monitor_sp_swarm_retain -ts_monitor_sp_swarm_phase 0 381c4762a1bSJed Brown test: 382c4762a1bSJed Brown suffix: euler 383d7462660SMatthew Knepley args: -dm_plex_box_faces 1,1 -dm_plex_box_lower -5,-5 -dm_plex_box_upper 5,5 \ 384d7462660SMatthew Knepley -ts_type euler -ts_max_time 0.1 -ts_convergence_estimate -convest_num_refine 2 \ 385d7462660SMatthew Knepley -dm_view -sw_view -ts_monitor_sp_swarm -ts_monitor_sp_swarm_retain -ts_monitor_sp_swarm_phase 0 386c4762a1bSJed Brown 387c4762a1bSJed Brown TEST*/ 388