1*db4653c2SDaniel Finn static char help[] = "Comparing basic symplectic, theta and discrete gradients interators on a simple hamiltonian system (harmonic oscillator) with particles\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 { 10*db4653c2SDaniel Finn PetscInt dim; /* The topological mesh dimension */ 11*db4653c2SDaniel Finn PetscBool simplex; /* Flag for simplices or tensor cells */ 12*db4653c2SDaniel Finn char filename[PETSC_MAX_PATH_LEN]; /* Name of the mesh filename if any */ 13c4762a1bSJed Brown PetscReal omega; /* Oscillation frequency omega */ 14c4762a1bSJed Brown PetscInt particlesPerCell; /* The number of partices per cell */ 15*db4653c2SDaniel Finn PetscInt numberOfCells; /* Number of cells in mesh */ 16c4762a1bSJed Brown PetscReal momentTol; /* Tolerance for checking moment conservation */ 17c4762a1bSJed Brown PetscBool monitor; /* Flag for using the TS monitor */ 18c4762a1bSJed Brown PetscBool error; /* Flag for printing the error */ 19c4762a1bSJed Brown PetscInt ostep; /* print the energy at each ostep time steps */ 20c4762a1bSJed Brown } AppCtx; 21c4762a1bSJed Brown 22c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 23c4762a1bSJed Brown { 24c4762a1bSJed Brown PetscErrorCode ierr; 25c4762a1bSJed Brown 26c4762a1bSJed Brown PetscFunctionBeginUser; 27*db4653c2SDaniel Finn options->dim = 2; 28*db4653c2SDaniel Finn options->simplex = PETSC_TRUE; 29c4762a1bSJed Brown options->monitor = PETSC_FALSE; 30c4762a1bSJed Brown options->error = PETSC_FALSE; 31c4762a1bSJed Brown options->particlesPerCell = 1; 32*db4653c2SDaniel Finn options->numberOfCells = 2; 33c4762a1bSJed Brown options->momentTol = 100.0*PETSC_MACHINE_EPSILON; 34c4762a1bSJed Brown options->omega = 64.0; 35c4762a1bSJed Brown options->ostep = 100; 36c4762a1bSJed Brown 37*db4653c2SDaniel Finn ierr = PetscStrcpy(options->filename, "");CHKERRQ(ierr); 38*db4653c2SDaniel Finn 39c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Harmonic Oscillator Options", "DMPLEX");CHKERRQ(ierr); 40c4762a1bSJed Brown ierr = PetscOptionsInt("-output_step", "Number of time steps between output", "ex4.c", options->ostep, &options->ostep, PETSC_NULL);CHKERRQ(ierr); 41*db4653c2SDaniel Finn ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex4.c", options->dim, &options->dim, NULL);CHKERRQ(ierr); 42c4762a1bSJed Brown ierr = PetscOptionsBool("-monitor", "Flag to use the TS monitor", "ex4.c", options->monitor, &options->monitor, NULL);CHKERRQ(ierr); 43c4762a1bSJed Brown ierr = PetscOptionsBool("-error", "Flag to print the error", "ex4.c", options->error, &options->error, NULL);CHKERRQ(ierr); 44*db4653c2SDaniel Finn ierr = PetscOptionsBool("-simplex", "The flag for simplices or tensor cells", "ex4.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr); 45*db4653c2SDaniel Finn ierr = PetscOptionsString("-mesh", "Name of the mesh filename if any", "ex4.c", options->filename, options->filename, sizeof(options->filename), NULL);CHKERRQ(ierr); 46c4762a1bSJed Brown ierr = PetscOptionsInt("-particles_per_cell", "Number of particles per cell", "ex4.c", options->particlesPerCell, &options->particlesPerCell, NULL);CHKERRQ(ierr); 47c4762a1bSJed Brown ierr = PetscOptionsReal("-omega", "Oscillator frequency", "ex4.c", options->omega, &options->omega, PETSC_NULL);CHKERRQ(ierr); 48c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 49c4762a1bSJed Brown 50c4762a1bSJed Brown PetscFunctionReturn(0); 51*db4653c2SDaniel Finn 52c4762a1bSJed Brown } 53c4762a1bSJed Brown 54c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user) 55c4762a1bSJed Brown { 56c4762a1bSJed Brown PetscErrorCode ierr; 57c4762a1bSJed Brown 58c4762a1bSJed Brown PetscFunctionBeginUser; 5930602db0SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 6030602db0SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 61c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 62c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 63*db4653c2SDaniel Finn 64c4762a1bSJed Brown PetscFunctionReturn(0); 65c4762a1bSJed Brown } 66c4762a1bSJed Brown 67c4762a1bSJed Brown static PetscErrorCode SetInitialCoordinates(DM dmSw) 68c4762a1bSJed Brown { 69c4762a1bSJed Brown DM dm; 70c4762a1bSJed Brown AppCtx *user; 71c4762a1bSJed Brown PetscRandom rnd; 72c4762a1bSJed Brown PetscBool simplex; 73c4762a1bSJed Brown PetscReal *centroid, *coords, *xi0, *v0, *J, *invJ, detJ; 74c4762a1bSJed Brown PetscInt dim, d, cStart, cEnd, c, Np, p; 75c4762a1bSJed Brown PetscErrorCode ierr; 76c4762a1bSJed Brown 77c4762a1bSJed Brown PetscFunctionBeginUser; 78c4762a1bSJed Brown ierr = PetscRandomCreate(PetscObjectComm((PetscObject) dmSw), &rnd);CHKERRQ(ierr); 79c4762a1bSJed Brown ierr = PetscRandomSetInterval(rnd, -1.0, 1.0);CHKERRQ(ierr); 80c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rnd);CHKERRQ(ierr); 81c4762a1bSJed Brown 823ec1f749SStefano Zampini ierr = DMGetApplicationContext(dmSw, &user);CHKERRQ(ierr); 83c4762a1bSJed Brown Np = user->particlesPerCell; 84c4762a1bSJed Brown ierr = DMSwarmGetCellDM(dmSw, &dm);CHKERRQ(ierr); 85c4762a1bSJed Brown ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 8630602db0SMatthew G. Knepley ierr = DMPlexIsSimplex(dm, &simplex);CHKERRQ(ierr); 87c4762a1bSJed Brown ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 88*db4653c2SDaniel Finn user->numberOfCells = cEnd - cStart; 89c4762a1bSJed Brown ierr = PetscMalloc5(dim, ¢roid, dim, &xi0, dim, &v0, dim*dim, &J, dim*dim, &invJ);CHKERRQ(ierr); 90c4762a1bSJed Brown for (d = 0; d < dim; ++d) xi0[d] = -1.0; 91c4762a1bSJed Brown ierr = DMSwarmGetField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 92c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 93c4762a1bSJed Brown if (Np == 1) { 94c4762a1bSJed Brown ierr = DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL);CHKERRQ(ierr); 95c4762a1bSJed Brown for (d = 0; d < dim; ++d) coords[c*dim+d] = centroid[d]; 96c4762a1bSJed Brown } else { 97c4762a1bSJed Brown ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); /* affine */ 98c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 99c4762a1bSJed Brown const PetscInt n = c*Np + p; 100c4762a1bSJed Brown PetscReal sum = 0.0, refcoords[3]; 101c4762a1bSJed Brown 102c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 103c4762a1bSJed Brown ierr = PetscRandomGetValueReal(rnd, &refcoords[d]);CHKERRQ(ierr); 104c4762a1bSJed Brown sum += refcoords[d]; 105c4762a1bSJed Brown } 106c4762a1bSJed Brown if (simplex && sum > 0.0) for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim)*sum; 107c4762a1bSJed Brown CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n*dim]); 108c4762a1bSJed Brown } 109c4762a1bSJed Brown } 110c4762a1bSJed Brown } 111*db4653c2SDaniel Finn 112c4762a1bSJed Brown ierr = DMSwarmRestoreField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 113c4762a1bSJed Brown ierr = PetscFree5(centroid, xi0, v0, J, invJ);CHKERRQ(ierr); 114c4762a1bSJed Brown ierr = PetscRandomDestroy(&rnd);CHKERRQ(ierr); 115c4762a1bSJed Brown PetscFunctionReturn(0); 116c4762a1bSJed Brown } 117c4762a1bSJed Brown 118c4762a1bSJed Brown static PetscErrorCode SetInitialConditions(DM dmSw, Vec u) 119c4762a1bSJed Brown { 120c4762a1bSJed Brown DM dm; 121c4762a1bSJed Brown AppCtx *user; 122c4762a1bSJed Brown PetscReal *coords; 123c4762a1bSJed Brown PetscScalar *initialConditions; 124c4762a1bSJed Brown PetscInt dim, cStart, cEnd, c, Np, p; 125c4762a1bSJed Brown PetscErrorCode ierr; 126c4762a1bSJed Brown 127c4762a1bSJed Brown PetscFunctionBeginUser; 1283ec1f749SStefano Zampini ierr = DMGetApplicationContext(dmSw, &user);CHKERRQ(ierr); 129c4762a1bSJed Brown Np = user->particlesPerCell; 130c4762a1bSJed Brown ierr = DMSwarmGetCellDM(dmSw, &dm);CHKERRQ(ierr); 131c4762a1bSJed Brown ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 132c4762a1bSJed Brown ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 133c4762a1bSJed Brown ierr = DMSwarmGetField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 134c4762a1bSJed Brown ierr = VecGetArray(u, &initialConditions);CHKERRQ(ierr); 135c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 136c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 137c4762a1bSJed Brown const PetscInt n = c*Np + p; 138c4762a1bSJed Brown initialConditions[n*2+0] = DMPlex_NormD_Internal(dim, &coords[n*dim]); 139c4762a1bSJed Brown initialConditions[n*2+1] = 0.0; 140c4762a1bSJed Brown } 141c4762a1bSJed Brown } 142c4762a1bSJed Brown ierr = VecRestoreArray(u, &initialConditions);CHKERRQ(ierr); 143c4762a1bSJed Brown ierr = DMSwarmRestoreField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 144c4762a1bSJed Brown PetscFunctionReturn(0); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 147c4762a1bSJed Brown static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user) 148c4762a1bSJed Brown { 149c4762a1bSJed Brown PetscInt *cellid; 150*db4653c2SDaniel Finn PetscInt dim, cStart, cEnd, c, Np, p; 151c4762a1bSJed Brown PetscErrorCode ierr; 152c4762a1bSJed Brown 153c4762a1bSJed Brown PetscFunctionBeginUser; 154c4762a1bSJed Brown ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 155c4762a1bSJed Brown ierr = DMCreate(PetscObjectComm((PetscObject) dm), sw);CHKERRQ(ierr); 156c4762a1bSJed Brown ierr = DMSetType(*sw, DMSWARM);CHKERRQ(ierr); 157c4762a1bSJed Brown ierr = DMSetDimension(*sw, dim);CHKERRQ(ierr); 158*db4653c2SDaniel Finn Np = user->particlesPerCell; 159c4762a1bSJed Brown ierr = DMSwarmSetType(*sw, DMSWARM_PIC);CHKERRQ(ierr); 160c4762a1bSJed Brown ierr = DMSwarmSetCellDM(*sw, dm);CHKERRQ(ierr); 161c4762a1bSJed Brown ierr = DMSwarmRegisterPetscDatatypeField(*sw, "kinematics", 2, PETSC_REAL);CHKERRQ(ierr); 162c4762a1bSJed Brown ierr = DMSwarmFinalizeFieldRegister(*sw);CHKERRQ(ierr); 163c4762a1bSJed Brown ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 164c4762a1bSJed Brown ierr = DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0);CHKERRQ(ierr); 165c4762a1bSJed Brown ierr = DMSetFromOptions(*sw);CHKERRQ(ierr); 166c4762a1bSJed Brown ierr = DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr); 167c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 168c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 169c4762a1bSJed Brown const PetscInt n = c*Np + p; 170c4762a1bSJed Brown 171c4762a1bSJed Brown cellid[n] = c; 172c4762a1bSJed Brown } 173c4762a1bSJed Brown } 174c4762a1bSJed Brown ierr = DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr); 175c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) *sw, "Particles");CHKERRQ(ierr); 176c4762a1bSJed Brown ierr = DMViewFromOptions(*sw, NULL, "-sw_view");CHKERRQ(ierr); 177c4762a1bSJed Brown PetscFunctionReturn(0); 178c4762a1bSJed Brown } 179c4762a1bSJed Brown 180c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) 181c4762a1bSJed Brown { 182c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 183c4762a1bSJed Brown const PetscReal omega = user->omega; 184c4762a1bSJed Brown const PetscScalar *u; 185c4762a1bSJed Brown MPI_Comm comm; 186c4762a1bSJed Brown PetscReal dt; 187c4762a1bSJed Brown PetscInt Np, p; 188c4762a1bSJed Brown PetscErrorCode ierr; 189c4762a1bSJed Brown 190c4762a1bSJed Brown PetscFunctionBeginUser; 191c4762a1bSJed Brown if (step%user->ostep == 0) { 192c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr); 193c4762a1bSJed Brown if (!step) {ierr = PetscPrintf(comm, "Time Step Part Energy Mod Energy\n");CHKERRQ(ierr);} 194c4762a1bSJed Brown ierr = TSGetTimeStep(ts, &dt);CHKERRQ(ierr); 195c4762a1bSJed Brown ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr); 196c4762a1bSJed Brown ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr); 197c4762a1bSJed Brown Np /= 2; 198c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 199c4762a1bSJed Brown const PetscReal x = PetscRealPart(u[p*2+0]); 200c4762a1bSJed Brown const PetscReal v = PetscRealPart(u[p*2+1]); 201c4762a1bSJed Brown const PetscReal E = 0.5*(v*v + PetscSqr(omega)*x*x); 202c4762a1bSJed Brown const PetscReal mE = 0.5*(v*v + PetscSqr(omega)*x*x - PetscSqr(omega)*dt*x*v); 203c4762a1bSJed Brown ierr = PetscPrintf(comm, "%.6lf %4D %4D %10.4lf %10.4lf\n", t, step, p, (double) E, (double) mE);CHKERRQ(ierr); 204c4762a1bSJed Brown } 205c4762a1bSJed Brown ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr); 206c4762a1bSJed Brown } 207c4762a1bSJed Brown PetscFunctionReturn(0); 208c4762a1bSJed Brown } 209c4762a1bSJed Brown 210c4762a1bSJed Brown static PetscErrorCode InitializeSolve(TS ts, Vec u) 211c4762a1bSJed Brown { 212c4762a1bSJed Brown DM dm; 213c4762a1bSJed Brown AppCtx *user; 214c4762a1bSJed Brown PetscErrorCode ierr; 215c4762a1bSJed Brown 216c4762a1bSJed Brown PetscFunctionBeginUser; 217c4762a1bSJed Brown ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 2183ec1f749SStefano Zampini ierr = DMGetApplicationContext(dm, &user);CHKERRQ(ierr); 219c4762a1bSJed Brown ierr = SetInitialCoordinates(dm);CHKERRQ(ierr); 220c4762a1bSJed Brown ierr = SetInitialConditions(dm, u);CHKERRQ(ierr); 221c4762a1bSJed Brown PetscFunctionReturn(0); 222c4762a1bSJed Brown } 223c4762a1bSJed Brown 224c4762a1bSJed Brown static PetscErrorCode ComputeError(TS ts, Vec U, Vec E) 225c4762a1bSJed Brown { 226c4762a1bSJed Brown MPI_Comm comm; 227c4762a1bSJed Brown DM sdm; 228c4762a1bSJed Brown AppCtx *user; 229c4762a1bSJed Brown const PetscScalar *u, *coords; 230c4762a1bSJed Brown PetscScalar *e; 231c4762a1bSJed Brown PetscReal t, omega; 232c4762a1bSJed Brown PetscInt dim, Np, p; 233c4762a1bSJed Brown PetscErrorCode ierr; 234c4762a1bSJed Brown 235c4762a1bSJed Brown PetscFunctionBeginUser; 236c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr); 237c4762a1bSJed Brown ierr = TSGetDM(ts, &sdm);CHKERRQ(ierr); 2383ec1f749SStefano Zampini ierr = DMGetApplicationContext(sdm, &user);CHKERRQ(ierr); 239c4762a1bSJed Brown omega = user->omega; 240c4762a1bSJed Brown ierr = DMGetDimension(sdm, &dim);CHKERRQ(ierr); 241c4762a1bSJed Brown ierr = TSGetSolveTime(ts, &t);CHKERRQ(ierr); 242c4762a1bSJed Brown ierr = VecGetArray(E, &e);CHKERRQ(ierr); 243c4762a1bSJed Brown ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr); 244c4762a1bSJed Brown ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr); 245c4762a1bSJed Brown Np /= 2; 246c4762a1bSJed Brown ierr = DMSwarmGetField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 247c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 248c4762a1bSJed Brown const PetscReal x = PetscRealPart(u[p*2+0]); 249c4762a1bSJed Brown const PetscReal v = PetscRealPart(u[p*2+1]); 250c4762a1bSJed Brown const PetscReal x0 = DMPlex_NormD_Internal(dim, &coords[p*dim]); 251c4762a1bSJed Brown const PetscReal ex = x0*PetscCosReal(omega*t); 252c4762a1bSJed Brown const PetscReal ev = -x0*omega*PetscSinReal(omega*t); 253c4762a1bSJed Brown 254c4762a1bSJed Brown if (user->error) {ierr = PetscPrintf(comm, "p%D error [%.2g %.2g] sol [%.6lf %.6lf] exact [%.6lf %.6lf] energy/exact energy %g / %g\n", p, (double) PetscAbsReal(x-ex), (double) PetscAbsReal(v-ev), (double) x, (double) v, (double) ex, (double) ev, 0.5*(v*v + PetscSqr(omega)*x*x), (double) 0.5*PetscSqr(omega*x0));} 255c4762a1bSJed Brown e[p*2+0] = x - ex; 256c4762a1bSJed Brown e[p*2+1] = v - ev; 257c4762a1bSJed Brown } 258c4762a1bSJed Brown ierr = DMSwarmRestoreField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 25940be0ff1SMatthew G. Knepley 260c4762a1bSJed Brown ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr); 261c4762a1bSJed Brown ierr = VecRestoreArray(E, &e);CHKERRQ(ierr); 262c4762a1bSJed Brown PetscFunctionReturn(0); 263c4762a1bSJed Brown } 264c4762a1bSJed Brown 26540be0ff1SMatthew G. Knepley /*---------------------Create particle RHS Functions--------------------------*/ 26640be0ff1SMatthew G. Knepley static PetscErrorCode RHSFunction1(TS ts, PetscReal t, Vec V, Vec Xres, void *ctx) 26740be0ff1SMatthew G. Knepley { 26840be0ff1SMatthew G. Knepley const PetscScalar *v; 26940be0ff1SMatthew G. Knepley PetscScalar *xres; 27040be0ff1SMatthew G. Knepley PetscInt Np, p; 27140be0ff1SMatthew G. Knepley PetscErrorCode ierr; 27240be0ff1SMatthew G. Knepley 27340be0ff1SMatthew G. Knepley PetscFunctionBeginUser; 27440be0ff1SMatthew G. Knepley ierr = VecGetArray(Xres, &xres);CHKERRQ(ierr); 27540be0ff1SMatthew G. Knepley ierr = VecGetArrayRead(V, &v);CHKERRQ(ierr); 27640be0ff1SMatthew G. Knepley ierr = VecGetLocalSize(Xres, &Np);CHKERRQ(ierr); 27740be0ff1SMatthew G. Knepley for (p = 0; p < Np; ++p) { 27840be0ff1SMatthew G. Knepley xres[p] = v[p]; 27940be0ff1SMatthew G. Knepley } 28040be0ff1SMatthew G. Knepley ierr = VecRestoreArrayRead(V, &v);CHKERRQ(ierr); 28140be0ff1SMatthew G. Knepley ierr = VecRestoreArray(Xres, &xres);CHKERRQ(ierr); 28240be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 28340be0ff1SMatthew G. Knepley } 28440be0ff1SMatthew G. Knepley 28540be0ff1SMatthew G. Knepley static PetscErrorCode RHSFunction2(TS ts, PetscReal t, Vec X, Vec Vres, void *ctx) 28640be0ff1SMatthew G. Knepley { 28740be0ff1SMatthew G. Knepley AppCtx *user = (AppCtx *)ctx; 28840be0ff1SMatthew G. Knepley const PetscScalar *x; 28940be0ff1SMatthew G. Knepley PetscScalar *vres; 29040be0ff1SMatthew G. Knepley PetscInt Np, p; 29140be0ff1SMatthew G. Knepley PetscErrorCode ierr; 29240be0ff1SMatthew G. Knepley 29340be0ff1SMatthew G. Knepley PetscFunctionBeginUser; 29440be0ff1SMatthew G. Knepley ierr = VecGetArray(Vres, &vres);CHKERRQ(ierr); 29540be0ff1SMatthew G. Knepley ierr = VecGetArrayRead(X, &x);CHKERRQ(ierr); 29640be0ff1SMatthew G. Knepley ierr = VecGetLocalSize(Vres, &Np);CHKERRQ(ierr); 29740be0ff1SMatthew G. Knepley for (p = 0; p < Np; ++p) { 29840be0ff1SMatthew G. Knepley vres[p] = -PetscSqr(user->omega)*x[p]; 29940be0ff1SMatthew G. Knepley } 30040be0ff1SMatthew G. Knepley ierr = VecRestoreArrayRead(X, &x);CHKERRQ(ierr); 30140be0ff1SMatthew G. Knepley ierr = VecRestoreArray(Vres, &vres);CHKERRQ(ierr); 30240be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 30340be0ff1SMatthew G. Knepley } 30440be0ff1SMatthew G. Knepley 30540be0ff1SMatthew G. Knepley /*----------------------------------------------------------------------------*/ 30640be0ff1SMatthew G. Knepley 30740be0ff1SMatthew G. Knepley /*--------------------Define RHSFunction, RHSJacobian (Theta)-----------------*/ 30840be0ff1SMatthew G. Knepley static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec G, void *ctx) 30940be0ff1SMatthew G. Knepley { 31040be0ff1SMatthew G. Knepley AppCtx *user = (AppCtx *) ctx; 31140be0ff1SMatthew G. Knepley DM dm; 31240be0ff1SMatthew G. Knepley const PetscScalar *u; 31340be0ff1SMatthew G. Knepley PetscScalar *g; 31440be0ff1SMatthew G. Knepley PetscInt Np, p; 31540be0ff1SMatthew G. Knepley PetscErrorCode ierr; 31640be0ff1SMatthew G. Knepley 31740be0ff1SMatthew G. Knepley PetscFunctionBeginUser; 31840be0ff1SMatthew G. Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 31940be0ff1SMatthew G. Knepley ierr = VecGetArray(G, &g);CHKERRQ(ierr); 32040be0ff1SMatthew G. Knepley ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr); 32140be0ff1SMatthew G. Knepley ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr); 32240be0ff1SMatthew G. Knepley Np /= 2; 32340be0ff1SMatthew G. Knepley for (p = 0; p < Np; ++p) { 32440be0ff1SMatthew G. Knepley g[p*2+0] = u[p*2+1]; 32540be0ff1SMatthew G. Knepley g[p*2+1] = -PetscSqr(user->omega)*u[p*2+0]; 32640be0ff1SMatthew G. Knepley } 32740be0ff1SMatthew G. Knepley ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr); 32840be0ff1SMatthew G. Knepley ierr = VecRestoreArray(G, &g);CHKERRQ(ierr); 32940be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 33040be0ff1SMatthew G. Knepley } 33140be0ff1SMatthew G. Knepley 33240be0ff1SMatthew G. Knepley /*Ji = dFi/dxj 33340be0ff1SMatthew G. Knepley J= (0 1) 33440be0ff1SMatthew G. Knepley (-w^2 0) 33540be0ff1SMatthew G. Knepley */ 33640be0ff1SMatthew G. Knepley static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U , Mat J, Mat P, void *ctx) 33740be0ff1SMatthew G. Knepley { 33840be0ff1SMatthew G. Knepley AppCtx *user = (AppCtx *) ctx; 339*db4653c2SDaniel Finn 340*db4653c2SDaniel Finn PetscInt Np; 34140be0ff1SMatthew G. Knepley PetscInt i, m, n; 34240be0ff1SMatthew G. Knepley const PetscScalar *u; 34340be0ff1SMatthew G. Knepley PetscScalar vals[4] = {0., 1., -PetscSqr(user->omega), 0.}; 34440be0ff1SMatthew G. Knepley PetscErrorCode ierr; 34540be0ff1SMatthew G. Knepley 346362febeeSStefano Zampini PetscFunctionBeginUser; 34740be0ff1SMatthew G. Knepley ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr); 348*db4653c2SDaniel Finn ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr); 349*db4653c2SDaniel Finn Np /= 2; 35040be0ff1SMatthew G. Knepley ierr = MatGetOwnershipRange(J, &m, &n);CHKERRQ(ierr); 35140be0ff1SMatthew G. Knepley for (i = 0; i < Np; ++i) { 35240be0ff1SMatthew G. Knepley const PetscInt rows[2] = {2*i, 2*i+1}; 35340be0ff1SMatthew G. Knepley ierr = MatSetValues(J, 2, rows, 2, rows, vals, INSERT_VALUES);CHKERRQ(ierr); 35440be0ff1SMatthew G. Knepley } 35540be0ff1SMatthew G. Knepley ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 35640be0ff1SMatthew G. Knepley ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 35740be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 35840be0ff1SMatthew G. Knepley 35940be0ff1SMatthew G. Knepley } 360*db4653c2SDaniel Finn 36140be0ff1SMatthew G. Knepley /*----------------------------------------------------------------------------*/ 36240be0ff1SMatthew G. Knepley 36340be0ff1SMatthew G. Knepley /*----------------Define S, F, G Functions (Discrete Gradients)---------------*/ 36440be0ff1SMatthew G. Knepley /* 36540be0ff1SMatthew G. Knepley u_t = S * gradF 36640be0ff1SMatthew G. Knepley --or-- 36740be0ff1SMatthew G. Knepley u_t = S * G 36840be0ff1SMatthew G. Knepley 36940be0ff1SMatthew G. Knepley + Sfunc - constructor for the S matrix from the formulation 37040be0ff1SMatthew G. Knepley . Ffunc - functional F from the formulation 37140be0ff1SMatthew G. Knepley - Gfunc - constructor for the gradient of F from the formulation 37240be0ff1SMatthew G. Knepley */ 37340be0ff1SMatthew G. Knepley 37440be0ff1SMatthew G. Knepley PetscErrorCode Sfunc(TS ts, PetscReal t, Vec U, Mat S, void *ctx) 37540be0ff1SMatthew G. Knepley { 376*db4653c2SDaniel Finn PetscInt Np; 37740be0ff1SMatthew G. Knepley PetscInt i, m, n; 378*db4653c2SDaniel Finn 37940be0ff1SMatthew G. Knepley const PetscScalar *u; 38040be0ff1SMatthew G. Knepley PetscScalar vals[4] = {0., 1., -1, 0.}; 38140be0ff1SMatthew G. Knepley PetscErrorCode ierr; 38240be0ff1SMatthew G. Knepley 38340be0ff1SMatthew G. Knepley PetscFunctionBeginUser; 384*db4653c2SDaniel Finn 38540be0ff1SMatthew G. Knepley ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr); 386*db4653c2SDaniel Finn ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr); 387*db4653c2SDaniel Finn Np /= 2; 38840be0ff1SMatthew G. Knepley ierr = MatGetOwnershipRange(S, &m, &n);CHKERRQ(ierr); 38940be0ff1SMatthew G. Knepley for (i = 0; i < Np; ++i) { 39040be0ff1SMatthew G. Knepley const PetscInt rows[2] = {2*i, 2*i+1}; 39140be0ff1SMatthew G. Knepley ierr = MatSetValues(S, 2, rows, 2, rows, vals, INSERT_VALUES);CHKERRQ(ierr); 39240be0ff1SMatthew G. Knepley } 39340be0ff1SMatthew G. Knepley ierr = MatAssemblyBegin(S,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 39440be0ff1SMatthew G. Knepley ierr = MatAssemblyEnd(S,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 39540be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 396*db4653c2SDaniel Finn 39740be0ff1SMatthew G. Knepley } 39840be0ff1SMatthew G. Knepley 39940be0ff1SMatthew G. Knepley PetscErrorCode Ffunc(TS ts, PetscReal t, Vec U, PetscScalar *F, void *ctx) 40040be0ff1SMatthew G. Knepley { 40140be0ff1SMatthew G. Knepley AppCtx *user = (AppCtx *) ctx; 40240be0ff1SMatthew G. Knepley DM dm; 40340be0ff1SMatthew G. Knepley const PetscScalar *u; 404*db4653c2SDaniel Finn PetscInt Np; 40540be0ff1SMatthew G. Knepley PetscInt p; 40640be0ff1SMatthew G. Knepley PetscErrorCode ierr; 40740be0ff1SMatthew G. Knepley 40840be0ff1SMatthew G. Knepley PetscFunctionBeginUser; 40940be0ff1SMatthew G. Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 410*db4653c2SDaniel Finn 411*db4653c2SDaniel Finn /* Define F */ 41240be0ff1SMatthew G. Knepley ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr); 413*db4653c2SDaniel Finn ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr); 414*db4653c2SDaniel Finn Np /= 2; 41540be0ff1SMatthew G. Knepley for (p = 0; p < Np; ++p) { 416*db4653c2SDaniel Finn *F += 0.5*PetscSqr(user->omega)*PetscSqr(u[p*2+0]) + 0.5*PetscSqr(u[p*2+1]); 41740be0ff1SMatthew G. Knepley } 41840be0ff1SMatthew G. Knepley ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr); 41940be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 42040be0ff1SMatthew G. Knepley } 42140be0ff1SMatthew G. Knepley 42240be0ff1SMatthew G. Knepley PetscErrorCode gradFfunc(TS ts, PetscReal t, Vec U, Vec gradF, void *ctx) 42340be0ff1SMatthew G. Knepley { 42440be0ff1SMatthew G. Knepley AppCtx *user = (AppCtx *) ctx; 42540be0ff1SMatthew G. Knepley DM dm; 42640be0ff1SMatthew G. Knepley const PetscScalar *u; 42740be0ff1SMatthew G. Knepley PetscScalar *g; 428*db4653c2SDaniel Finn PetscInt Np; 42940be0ff1SMatthew G. Knepley PetscInt p; 43040be0ff1SMatthew G. Knepley PetscErrorCode ierr; 43140be0ff1SMatthew G. Knepley 43240be0ff1SMatthew G. Knepley PetscFunctionBeginUser; 43340be0ff1SMatthew G. Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 43440be0ff1SMatthew G. Knepley ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr); 435*db4653c2SDaniel Finn ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr); 436*db4653c2SDaniel Finn Np /= 2; 437*db4653c2SDaniel Finn /*Define gradF*/ 43840be0ff1SMatthew G. Knepley ierr = VecGetArray(gradF, &g);CHKERRQ(ierr); 43940be0ff1SMatthew G. Knepley for (p = 0; p < Np; ++p) { 44040be0ff1SMatthew G. Knepley g[p*2+0] = PetscSqr(user->omega)*u[p*2+0]; /*dF/dx*/ 44140be0ff1SMatthew G. Knepley g[p*2+1] = u[p*2+1]; /*dF/dv*/ 44240be0ff1SMatthew G. Knepley } 44340be0ff1SMatthew G. Knepley ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr); 44440be0ff1SMatthew G. Knepley ierr = VecRestoreArray(gradF, &g);CHKERRQ(ierr); 44540be0ff1SMatthew G. Knepley PetscFunctionReturn(0); 44640be0ff1SMatthew G. Knepley } 447*db4653c2SDaniel Finn 44840be0ff1SMatthew G. Knepley /*----------------------------------------------------------------------------*/ 44940be0ff1SMatthew G. Knepley 450c4762a1bSJed Brown int main(int argc,char **argv) 451c4762a1bSJed Brown { 452c4762a1bSJed Brown TS ts; /* nonlinear solver */ 453c4762a1bSJed Brown DM dm, sw; /* Mesh and particle managers */ 454c4762a1bSJed Brown Vec u; /* swarm vector */ 45540be0ff1SMatthew G. Knepley PetscInt n; /* swarm vector size */ 456c4762a1bSJed Brown IS is1, is2; 457c4762a1bSJed Brown MPI_Comm comm; 45840be0ff1SMatthew G. Knepley Mat J; /* Jacobian matrix */ 459c4762a1bSJed Brown AppCtx user; 460c4762a1bSJed Brown PetscErrorCode ierr; 461c4762a1bSJed Brown 46240be0ff1SMatthew G. Knepley /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 46340be0ff1SMatthew G. Knepley Initialize program 46440be0ff1SMatthew G. Knepley - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 465c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 466c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 467c4762a1bSJed Brown ierr = ProcessOptions(comm, &user);CHKERRQ(ierr); 468c4762a1bSJed Brown 46940be0ff1SMatthew G. Knepley /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 47040be0ff1SMatthew G. Knepley Create Particle-Mesh 47140be0ff1SMatthew G. Knepley - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 472c4762a1bSJed Brown ierr = CreateMesh(comm, &dm, &user);CHKERRQ(ierr); 473c4762a1bSJed Brown ierr = CreateParticles(dm, &sw, &user);CHKERRQ(ierr); 474c4762a1bSJed Brown ierr = DMSetApplicationContext(sw, &user);CHKERRQ(ierr); 475c4762a1bSJed Brown 47640be0ff1SMatthew G. Knepley /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 47740be0ff1SMatthew G. Knepley Setup timestepping solver 47840be0ff1SMatthew G. Knepley - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 479c4762a1bSJed Brown ierr = TSCreate(comm, &ts);CHKERRQ(ierr); 48040be0ff1SMatthew G. Knepley ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 48140be0ff1SMatthew G. Knepley 482c4762a1bSJed Brown ierr = TSSetDM(ts, sw);CHKERRQ(ierr); 483c4762a1bSJed Brown ierr = TSSetMaxTime(ts, 0.1);CHKERRQ(ierr); 484c4762a1bSJed Brown ierr = TSSetTimeStep(ts, 0.00001);CHKERRQ(ierr); 485c4762a1bSJed Brown ierr = TSSetMaxSteps(ts, 100);CHKERRQ(ierr); 486c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 487c4762a1bSJed Brown if (user.monitor) {ierr = TSMonitorSet(ts, Monitor, &user, NULL);CHKERRQ(ierr);} 488c4762a1bSJed Brown 48940be0ff1SMatthew G. Knepley /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 49040be0ff1SMatthew G. Knepley Prepare to solve 49140be0ff1SMatthew G. Knepley - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 492c4762a1bSJed Brown ierr = DMSwarmCreateGlobalVectorFromField(sw, "kinematics", &u);CHKERRQ(ierr); 493c4762a1bSJed Brown ierr = VecGetLocalSize(u, &n);CHKERRQ(ierr); 49440be0ff1SMatthew G. Knepley ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 49540be0ff1SMatthew G. Knepley ierr = TSSetComputeInitialCondition(ts, InitializeSolve);CHKERRQ(ierr); 49640be0ff1SMatthew G. Knepley ierr = TSSetComputeExactError(ts, ComputeError);CHKERRQ(ierr); 49740be0ff1SMatthew G. Knepley 49840be0ff1SMatthew G. Knepley /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 49940be0ff1SMatthew G. Knepley Define function F(U, Udot , x , t) = G(U, x, t) 50040be0ff1SMatthew G. Knepley - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 50140be0ff1SMatthew G. Knepley 50240be0ff1SMatthew G. Knepley /* - - - - - - - Basic Symplectic - - - - - - - - - - - - - - - - - - - - - -*/ 503c4762a1bSJed Brown ierr = ISCreateStride(comm, n/2, 0, 2, &is1);CHKERRQ(ierr); 504c4762a1bSJed Brown ierr = ISCreateStride(comm, n/2, 1, 2, &is2);CHKERRQ(ierr); 505c4762a1bSJed Brown ierr = TSRHSSplitSetIS(ts, "position", is1);CHKERRQ(ierr); 506c4762a1bSJed Brown ierr = TSRHSSplitSetIS(ts, "momentum", is2);CHKERRQ(ierr); 507c4762a1bSJed Brown ierr = ISDestroy(&is1);CHKERRQ(ierr); 508c4762a1bSJed Brown ierr = ISDestroy(&is2);CHKERRQ(ierr); 509c4762a1bSJed Brown ierr = TSRHSSplitSetRHSFunction(ts, "position", NULL, RHSFunction1, &user);CHKERRQ(ierr); 510c4762a1bSJed Brown ierr = TSRHSSplitSetRHSFunction(ts, "momentum", NULL, RHSFunction2, &user);CHKERRQ(ierr); 511c4762a1bSJed Brown 51240be0ff1SMatthew G. Knepley /* - - - - - - - Theta (Implicit Midpoint) - - - - - - - - - - - - - - - - -*/ 513*db4653c2SDaniel Finn 51440be0ff1SMatthew G. Knepley ierr = TSSetRHSFunction(ts, NULL, RHSFunction, &user);CHKERRQ(ierr); 51540be0ff1SMatthew G. Knepley 51640be0ff1SMatthew G. Knepley ierr = MatCreate(PETSC_COMM_WORLD,&J);CHKERRQ(ierr); 51740be0ff1SMatthew G. Knepley ierr = MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n);CHKERRQ(ierr); 51840be0ff1SMatthew G. Knepley ierr = MatSetFromOptions(J);CHKERRQ(ierr); 51940be0ff1SMatthew G. Knepley ierr = MatSetUp(J);CHKERRQ(ierr); 52040be0ff1SMatthew G. Knepley ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 52140be0ff1SMatthew G. Knepley ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 52240be0ff1SMatthew G. Knepley ierr = TSSetRHSJacobian(ts,J,J,RHSJacobian,&user);CHKERRQ(ierr); 52340be0ff1SMatthew G. Knepley 52440be0ff1SMatthew G. Knepley /* - - - - - - - Discrete Gradients - - - - - - - - - - - - - - - - - - - - */ 52540be0ff1SMatthew G. Knepley ierr = TSDiscGradSetFormulation(ts, Sfunc, Ffunc, gradFfunc, &user);CHKERRQ(ierr); 52640be0ff1SMatthew G. Knepley 52740be0ff1SMatthew G. Knepley /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 52840be0ff1SMatthew G. Knepley Solve 52940be0ff1SMatthew G. Knepley - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 530c4762a1bSJed Brown ierr = TSComputeInitialCondition(ts, u);CHKERRQ(ierr); 531c4762a1bSJed Brown ierr = TSSolve(ts, u);CHKERRQ(ierr); 532c4762a1bSJed Brown 53340be0ff1SMatthew G. Knepley /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 53440be0ff1SMatthew G. Knepley Clean up workspace 53540be0ff1SMatthew G. Knepley - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 53640be0ff1SMatthew G. Knepley ierr = MatDestroy(&J);CHKERRQ(ierr); 53740be0ff1SMatthew G. Knepley ierr = DMSwarmDestroyGlobalVectorFromField(sw, "kinematics", &u);CHKERRQ(ierr); 538c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 539c4762a1bSJed Brown ierr = DMDestroy(&sw);CHKERRQ(ierr); 540c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 541c4762a1bSJed Brown ierr = PetscFinalize(); 542c4762a1bSJed Brown return ierr; 543c4762a1bSJed Brown } 544c4762a1bSJed Brown 545c4762a1bSJed Brown /*TEST 546c4762a1bSJed Brown 547c4762a1bSJed Brown build: 548c4762a1bSJed Brown requires: triangle !single !complex 549c4762a1bSJed Brown test: 550c4762a1bSJed Brown suffix: 1 551*db4653c2SDaniel Finn args: -dm_plex_box_faces 1,1 -ts_type basicsymplectic -ts_basicsymplectic_type 1 -ts_convergence_estimate -convest_num_refine 2 -dm_view -monitor -output_step 50 -error 55240be0ff1SMatthew G. Knepley 553c4762a1bSJed Brown test: 554c4762a1bSJed Brown suffix: 2 555*db4653c2SDaniel Finn args: -dm_plex_box_faces 1,1 -ts_type basicsymplectic -ts_basicsymplectic_type 2 -ts_convergence_estimate -convest_num_refine 2 -dm_view -monitor -output_step 50 -error 55640be0ff1SMatthew G. Knepley 557c4762a1bSJed Brown test: 558c4762a1bSJed Brown suffix: 3 559*db4653c2SDaniel Finn args: -dm_plex_box_faces 1,1 -ts_type basicsymplectic -ts_basicsymplectic_type 3 -ts_convergence_estimate -convest_num_refine 2 -ts_dt 0.0001 -dm_view -monitor -output_step 50 -error 56040be0ff1SMatthew G. Knepley 561c4762a1bSJed Brown test: 562c4762a1bSJed Brown suffix: 4 563*db4653c2SDaniel Finn args: -dm_plex_box_faces 1,1 -ts_type basicsymplectic -ts_basicsymplectic_type 4 -ts_convergence_estimate -convest_num_refine 2 -ts_dt 0.0001 -dm_view -monitor -output_step 50 -error 56440be0ff1SMatthew G. Knepley 56540be0ff1SMatthew G. Knepley test: 56640be0ff1SMatthew G. Knepley suffix: 5 567*db4653c2SDaniel Finn args: -dm_plex_box_faces 1,1 -ts_type theta -ts_theta_theta 0.5 -monitor -output_step 50 -error -ts_convergence_estimate -convest_num_refine 2 -dm_view 56840be0ff1SMatthew G. Knepley 56940be0ff1SMatthew G. Knepley test: 57040be0ff1SMatthew G. Knepley suffix: 6 571*db4653c2SDaniel Finn args: -dm_plex_box_faces 1,1 -ts_type discgrad -monitor -output_step 50 -error -ts_convergence_estimate -convest_num_refine 2 -dm_view 572*db4653c2SDaniel Finn 573*db4653c2SDaniel Finn test: 574*db4653c2SDaniel Finn suffix: 7 575*db4653c2SDaniel Finn args: -dm_plex_box_faces 1,1 -ts_type discgrad -ts_discgrad_gonzalez -monitor -output_step 50 -error -ts_convergence_estimate -convest_num_refine 2 -dm_view 576c4762a1bSJed Brown 577c4762a1bSJed Brown TEST*/ 578