1*c4762a1bSJed Brown static char help[] = "Vlasov example of particles orbiting around a central massive point.\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscdmplex.h> 4*c4762a1bSJed Brown #include <petsc/private/dmpleximpl.h> /* For norm */ 5*c4762a1bSJed Brown #include <petsc/private/petscfeimpl.h> /* For CoordinatesRefToReal() */ 6*c4762a1bSJed Brown #include <petscdmswarm.h> 7*c4762a1bSJed Brown #include <petscts.h> 8*c4762a1bSJed Brown 9*c4762a1bSJed Brown typedef struct { 10*c4762a1bSJed Brown PetscInt dim; /* The topological mesh dimension */ 11*c4762a1bSJed Brown PetscBool simplex; /* Flag for simplices or tensor cells */ 12*c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN]; /* Name of the mesh filename if any */ 13*c4762a1bSJed Brown PetscInt particlesPerCell; /* The number of partices per cell */ 14*c4762a1bSJed Brown PetscReal momentTol; /* Tolerance for checking moment conservation */ 15*c4762a1bSJed Brown PetscBool monitor; /* Flag for using the TS monitor */ 16*c4762a1bSJed Brown PetscBool error; /* Flag for printing the error */ 17*c4762a1bSJed Brown PetscInt ostep; /* print the energy at each ostep time steps */ 18*c4762a1bSJed Brown } AppCtx; 19*c4762a1bSJed Brown 20*c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 21*c4762a1bSJed Brown { 22*c4762a1bSJed Brown PetscErrorCode ierr; 23*c4762a1bSJed Brown 24*c4762a1bSJed Brown PetscFunctionBeginUser; 25*c4762a1bSJed Brown options->dim = 2; 26*c4762a1bSJed Brown options->simplex = PETSC_TRUE; 27*c4762a1bSJed Brown options->monitor = PETSC_FALSE; 28*c4762a1bSJed Brown options->error = PETSC_FALSE; 29*c4762a1bSJed Brown options->particlesPerCell = 1; 30*c4762a1bSJed Brown options->momentTol = 100.0*PETSC_MACHINE_EPSILON; 31*c4762a1bSJed Brown options->ostep = 100; 32*c4762a1bSJed Brown 33*c4762a1bSJed Brown ierr = PetscStrcpy(options->filename, "");CHKERRQ(ierr); 34*c4762a1bSJed Brown 35*c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Vlasov Options", "DMPLEX");CHKERRQ(ierr); 36*c4762a1bSJed Brown ierr = PetscOptionsInt("-output_step", "Number of time steps between output", "ex4.c", options->ostep, &options->ostep, PETSC_NULL);CHKERRQ(ierr); 37*c4762a1bSJed Brown ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex4.c", options->dim, &options->dim, NULL);CHKERRQ(ierr); 38*c4762a1bSJed Brown ierr = PetscOptionsBool("-monitor", "Flag to use the TS monitor", "ex4.c", options->monitor, &options->monitor, NULL);CHKERRQ(ierr); 39*c4762a1bSJed Brown ierr = PetscOptionsBool("-error", "Flag to print the error", "ex4.c", options->error, &options->error, NULL);CHKERRQ(ierr); 40*c4762a1bSJed Brown ierr = PetscOptionsBool("-simplex", "The flag for simplices or tensor cells", "ex4.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr); 41*c4762a1bSJed Brown ierr = PetscOptionsString("-mesh", "Name of the mesh filename if any", "ex4.c", options->filename, options->filename, PETSC_MAX_PATH_LEN, NULL);CHKERRQ(ierr); 42*c4762a1bSJed Brown ierr = PetscOptionsInt("-particles_per_cell", "Number of particles per cell", "ex4.c", options->particlesPerCell, &options->particlesPerCell, NULL);CHKERRQ(ierr); 43*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 44*c4762a1bSJed Brown 45*c4762a1bSJed Brown PetscFunctionReturn(0); 46*c4762a1bSJed Brown } 47*c4762a1bSJed Brown 48*c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user) 49*c4762a1bSJed Brown { 50*c4762a1bSJed Brown PetscBool flg; 51*c4762a1bSJed Brown PetscErrorCode ierr; 52*c4762a1bSJed Brown 53*c4762a1bSJed Brown PetscFunctionBeginUser; 54*c4762a1bSJed Brown ierr = PetscStrcmp(user->filename, "", &flg);CHKERRQ(ierr); 55*c4762a1bSJed Brown if (flg) { 56*c4762a1bSJed Brown ierr = DMPlexCreateBoxMesh(comm, user->dim, user->simplex, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr); 57*c4762a1bSJed Brown } else { 58*c4762a1bSJed Brown ierr = DMPlexCreateFromFile(comm, user->filename, PETSC_TRUE, dm);CHKERRQ(ierr); 59*c4762a1bSJed Brown ierr = DMGetDimension(*dm, &user->dim);CHKERRQ(ierr); 60*c4762a1bSJed Brown } 61*c4762a1bSJed Brown { 62*c4762a1bSJed Brown DM distributedMesh = NULL; 63*c4762a1bSJed Brown 64*c4762a1bSJed Brown ierr = DMPlexDistribute(*dm, 0, NULL, &distributedMesh);CHKERRQ(ierr); 65*c4762a1bSJed Brown if (distributedMesh) { 66*c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 67*c4762a1bSJed Brown *dm = distributedMesh; 68*c4762a1bSJed Brown } 69*c4762a1bSJed Brown } 70*c4762a1bSJed Brown ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr); /* needed for periodic */ 71*c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 72*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr); 73*c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 74*c4762a1bSJed Brown PetscFunctionReturn(0); 75*c4762a1bSJed Brown } 76*c4762a1bSJed Brown 77*c4762a1bSJed Brown static PetscErrorCode SetInitialCoordinates(DM dmSw) 78*c4762a1bSJed Brown { 79*c4762a1bSJed Brown DM dm; 80*c4762a1bSJed Brown AppCtx *user; 81*c4762a1bSJed Brown PetscRandom rnd; 82*c4762a1bSJed Brown PetscBool simplex; 83*c4762a1bSJed Brown PetscReal *centroid, *coords, *xi0, *v0, *J, *invJ, detJ; 84*c4762a1bSJed Brown PetscInt dim, d, cStart, cEnd, c, Np, p; 85*c4762a1bSJed Brown PetscErrorCode ierr; 86*c4762a1bSJed Brown 87*c4762a1bSJed Brown PetscFunctionBeginUser; 88*c4762a1bSJed Brown ierr = PetscRandomCreate(PetscObjectComm((PetscObject) dmSw), &rnd);CHKERRQ(ierr); 89*c4762a1bSJed Brown ierr = PetscRandomSetInterval(rnd, -1.0, 1.0);CHKERRQ(ierr); 90*c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rnd);CHKERRQ(ierr); 91*c4762a1bSJed Brown 92*c4762a1bSJed Brown ierr = DMGetApplicationContext(dmSw, (void **) &user);CHKERRQ(ierr); 93*c4762a1bSJed Brown simplex = user->simplex; 94*c4762a1bSJed Brown Np = user->particlesPerCell; 95*c4762a1bSJed Brown ierr = DMSwarmGetCellDM(dmSw, &dm);CHKERRQ(ierr); 96*c4762a1bSJed Brown ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 97*c4762a1bSJed Brown ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 98*c4762a1bSJed Brown ierr = PetscMalloc5(dim, ¢roid, dim, &xi0, dim, &v0, dim*dim, &J, dim*dim, &invJ);CHKERRQ(ierr); 99*c4762a1bSJed Brown for (d = 0; d < dim; ++d) xi0[d] = -1.0; 100*c4762a1bSJed Brown ierr = DMSwarmGetField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 101*c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 102*c4762a1bSJed Brown if (Np == 1) { 103*c4762a1bSJed Brown ierr = DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL);CHKERRQ(ierr); 104*c4762a1bSJed Brown for (d = 0; d < dim; ++d) coords[c*dim+d] = centroid[d]; 105*c4762a1bSJed Brown } else { 106*c4762a1bSJed Brown ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); /* affine */ 107*c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 108*c4762a1bSJed Brown const PetscInt n = c*Np + p; 109*c4762a1bSJed Brown PetscReal sum = 0.0, refcoords[3]; 110*c4762a1bSJed Brown 111*c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 112*c4762a1bSJed Brown ierr = PetscRandomGetValueReal(rnd, &refcoords[d]);CHKERRQ(ierr); 113*c4762a1bSJed Brown sum += refcoords[d]; 114*c4762a1bSJed Brown } 115*c4762a1bSJed Brown if (simplex && sum > 0.0) for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim)*sum; 116*c4762a1bSJed Brown CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n*dim]); 117*c4762a1bSJed Brown } 118*c4762a1bSJed Brown } 119*c4762a1bSJed Brown } 120*c4762a1bSJed Brown ierr = DMSwarmRestoreField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 121*c4762a1bSJed Brown ierr = PetscFree5(centroid, xi0, v0, J, invJ);CHKERRQ(ierr); 122*c4762a1bSJed Brown ierr = PetscRandomDestroy(&rnd);CHKERRQ(ierr); 123*c4762a1bSJed Brown PetscFunctionReturn(0); 124*c4762a1bSJed Brown } 125*c4762a1bSJed Brown 126*c4762a1bSJed Brown static PetscErrorCode SetInitialConditions(DM dmSw, Vec u) 127*c4762a1bSJed Brown { 128*c4762a1bSJed Brown DM dm; 129*c4762a1bSJed Brown AppCtx *user; 130*c4762a1bSJed Brown PetscScalar *initialConditions; 131*c4762a1bSJed Brown PetscInt dim, cStart, cEnd, c, Np, p; 132*c4762a1bSJed Brown PetscErrorCode ierr; 133*c4762a1bSJed Brown 134*c4762a1bSJed Brown PetscFunctionBeginUser; 135*c4762a1bSJed Brown ierr = DMGetApplicationContext(dmSw, (void **) &user);CHKERRQ(ierr); 136*c4762a1bSJed Brown Np = user->particlesPerCell; 137*c4762a1bSJed Brown ierr = DMSwarmGetCellDM(dmSw, &dm);CHKERRQ(ierr); 138*c4762a1bSJed Brown ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 139*c4762a1bSJed Brown ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 140*c4762a1bSJed Brown ierr = VecGetArray(u, &initialConditions);CHKERRQ(ierr); 141*c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 142*c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 143*c4762a1bSJed Brown const PetscInt n = c*Np + p; 144*c4762a1bSJed Brown 145*c4762a1bSJed Brown initialConditions[(n*2 + 0)*dim + 0] = n+1; 146*c4762a1bSJed Brown initialConditions[(n*2 + 0)*dim + 1] = 0; 147*c4762a1bSJed Brown initialConditions[(n*2 + 1)*dim + 0] = 0; 148*c4762a1bSJed Brown initialConditions[(n*2 + 1)*dim + 1] = PetscSqrtReal(1000./(n+1.)); 149*c4762a1bSJed Brown } 150*c4762a1bSJed Brown } 151*c4762a1bSJed Brown ierr = VecRestoreArray(u, &initialConditions);CHKERRQ(ierr); 152*c4762a1bSJed Brown PetscFunctionReturn(0); 153*c4762a1bSJed Brown } 154*c4762a1bSJed Brown 155*c4762a1bSJed Brown static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user) 156*c4762a1bSJed Brown { 157*c4762a1bSJed Brown PetscInt *cellid; 158*c4762a1bSJed Brown PetscInt dim, cStart, cEnd, c, Np = user->particlesPerCell, p; 159*c4762a1bSJed Brown PetscErrorCode ierr; 160*c4762a1bSJed Brown 161*c4762a1bSJed Brown PetscFunctionBeginUser; 162*c4762a1bSJed Brown ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 163*c4762a1bSJed Brown ierr = DMCreate(PetscObjectComm((PetscObject) dm), sw);CHKERRQ(ierr); 164*c4762a1bSJed Brown ierr = DMSetType(*sw, DMSWARM);CHKERRQ(ierr); 165*c4762a1bSJed Brown ierr = DMSetDimension(*sw, dim);CHKERRQ(ierr); 166*c4762a1bSJed Brown 167*c4762a1bSJed Brown ierr = DMSwarmSetType(*sw, DMSWARM_PIC);CHKERRQ(ierr); 168*c4762a1bSJed Brown ierr = DMSwarmSetCellDM(*sw, dm);CHKERRQ(ierr); 169*c4762a1bSJed Brown ierr = DMSwarmRegisterPetscDatatypeField(*sw, "kinematics", 2*dim, PETSC_REAL);CHKERRQ(ierr); 170*c4762a1bSJed Brown ierr = DMSwarmFinalizeFieldRegister(*sw);CHKERRQ(ierr); 171*c4762a1bSJed Brown ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 172*c4762a1bSJed Brown ierr = DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0);CHKERRQ(ierr); 173*c4762a1bSJed Brown ierr = DMSetFromOptions(*sw);CHKERRQ(ierr); 174*c4762a1bSJed Brown ierr = DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr); 175*c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 176*c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 177*c4762a1bSJed Brown const PetscInt n = c*Np + p; 178*c4762a1bSJed Brown 179*c4762a1bSJed Brown cellid[n] = c; 180*c4762a1bSJed Brown } 181*c4762a1bSJed Brown } 182*c4762a1bSJed Brown ierr = DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr); 183*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) *sw, "Particles");CHKERRQ(ierr); 184*c4762a1bSJed Brown ierr = DMViewFromOptions(*sw, NULL, "-sw_view");CHKERRQ(ierr); 185*c4762a1bSJed Brown PetscFunctionReturn(0); 186*c4762a1bSJed Brown } 187*c4762a1bSJed Brown 188*c4762a1bSJed Brown /* Create particle RHS Functions for gravity with G = 1 for simplification */ 189*c4762a1bSJed Brown static PetscErrorCode RHSFunction1(TS ts, PetscReal t, Vec V, Vec Xres, void *ctx) 190*c4762a1bSJed Brown { 191*c4762a1bSJed Brown DM dm; 192*c4762a1bSJed Brown const PetscScalar *v; 193*c4762a1bSJed Brown PetscScalar *xres; 194*c4762a1bSJed Brown PetscInt Np, p, dim, d; 195*c4762a1bSJed Brown PetscErrorCode ierr; 196*c4762a1bSJed Brown 197*c4762a1bSJed Brown PetscFunctionBeginUser; 198*c4762a1bSJed Brown ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 199*c4762a1bSJed Brown ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 200*c4762a1bSJed Brown ierr = VecGetLocalSize(Xres, &Np);CHKERRQ(ierr); 201*c4762a1bSJed Brown Np /= dim; 202*c4762a1bSJed Brown ierr = VecGetArray(Xres, &xres);CHKERRQ(ierr); 203*c4762a1bSJed Brown ierr = VecGetArrayRead(V, &v);CHKERRQ(ierr); 204*c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 205*c4762a1bSJed Brown for(d = 0; d < dim; ++d) { 206*c4762a1bSJed Brown xres[p*dim+d] = v[p*dim+d]; 207*c4762a1bSJed Brown } 208*c4762a1bSJed Brown } 209*c4762a1bSJed Brown ierr = VecRestoreArrayRead(V,& v);CHKERRQ(ierr); 210*c4762a1bSJed Brown ierr = VecRestoreArray(Xres, &xres);CHKERRQ(ierr); 211*c4762a1bSJed Brown PetscFunctionReturn(0); 212*c4762a1bSJed Brown } 213*c4762a1bSJed Brown 214*c4762a1bSJed Brown static PetscErrorCode RHSFunction2(TS ts, PetscReal t, Vec X, Vec Vres, void *ctx) 215*c4762a1bSJed Brown { 216*c4762a1bSJed Brown DM dm; 217*c4762a1bSJed Brown const PetscScalar *x; 218*c4762a1bSJed Brown PetscScalar *vres; 219*c4762a1bSJed Brown PetscInt Np, p, dim, d; 220*c4762a1bSJed Brown PetscErrorCode ierr; 221*c4762a1bSJed Brown 222*c4762a1bSJed Brown 223*c4762a1bSJed Brown PetscFunctionBeginUser; 224*c4762a1bSJed Brown ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 225*c4762a1bSJed Brown ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 226*c4762a1bSJed Brown ierr = VecGetLocalSize(Vres, &Np);CHKERRQ(ierr); 227*c4762a1bSJed Brown Np /= dim; 228*c4762a1bSJed Brown ierr = VecGetArray(Vres, &vres);CHKERRQ(ierr); 229*c4762a1bSJed Brown ierr = VecGetArrayRead(X, &x);CHKERRQ(ierr); 230*c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 231*c4762a1bSJed Brown const PetscScalar rsqr = DMPlex_NormD_Internal(dim, &x[p*dim]); 232*c4762a1bSJed Brown 233*c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 234*c4762a1bSJed Brown vres[p*dim+d] = -(1000./(p+1.)) * x[p*dim+d]/PetscSqr(rsqr); 235*c4762a1bSJed Brown } 236*c4762a1bSJed Brown } 237*c4762a1bSJed Brown ierr = VecRestoreArray(Vres, &vres);CHKERRQ(ierr); 238*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X, &x);CHKERRQ(ierr); 239*c4762a1bSJed Brown PetscFunctionReturn(0); 240*c4762a1bSJed Brown } 241*c4762a1bSJed Brown 242*c4762a1bSJed Brown static PetscErrorCode RHSFunctionParticles(TS ts, PetscReal t , Vec U, Vec R, void *ctx) 243*c4762a1bSJed Brown { 244*c4762a1bSJed Brown DM dm; 245*c4762a1bSJed Brown const PetscScalar *u; 246*c4762a1bSJed Brown PetscScalar *r; 247*c4762a1bSJed Brown PetscInt Np, p, dim, d; 248*c4762a1bSJed Brown PetscErrorCode ierr; 249*c4762a1bSJed Brown 250*c4762a1bSJed Brown PetscFunctionBeginUser; 251*c4762a1bSJed Brown ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 252*c4762a1bSJed Brown ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 253*c4762a1bSJed Brown ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr); 254*c4762a1bSJed Brown Np /= 2*dim; 255*c4762a1bSJed Brown ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr); 256*c4762a1bSJed Brown ierr = VecGetArray(R, &r);CHKERRQ(ierr); 257*c4762a1bSJed Brown for( p = 0; p < Np; ++p) { 258*c4762a1bSJed Brown const PetscScalar rsqr = DMPlex_NormD_Internal(dim, &u[p*2*dim]); 259*c4762a1bSJed Brown 260*c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 261*c4762a1bSJed Brown r[p*2*dim+d] = u[p*2*dim+d+2]; 262*c4762a1bSJed Brown r[p*2*dim+d+2] = -(1000./(1.+p)) * u[p*2*dim+d]/PetscSqr(rsqr); 263*c4762a1bSJed Brown } 264*c4762a1bSJed Brown } 265*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr); 266*c4762a1bSJed Brown ierr = VecRestoreArray(R, &r);CHKERRQ(ierr); 267*c4762a1bSJed Brown PetscFunctionReturn(0); 268*c4762a1bSJed Brown } 269*c4762a1bSJed Brown 270*c4762a1bSJed Brown static PetscErrorCode InitializeSolve(TS ts, Vec u) 271*c4762a1bSJed Brown { 272*c4762a1bSJed Brown DM dm; 273*c4762a1bSJed Brown AppCtx *user; 274*c4762a1bSJed Brown PetscErrorCode ierr; 275*c4762a1bSJed Brown 276*c4762a1bSJed Brown PetscFunctionBeginUser; 277*c4762a1bSJed Brown ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 278*c4762a1bSJed Brown ierr = DMGetApplicationContext(dm, (void **) &user);CHKERRQ(ierr); 279*c4762a1bSJed Brown ierr = SetInitialCoordinates(dm);CHKERRQ(ierr); 280*c4762a1bSJed Brown ierr = SetInitialConditions(dm, u);CHKERRQ(ierr); 281*c4762a1bSJed Brown PetscFunctionReturn(0); 282*c4762a1bSJed Brown } 283*c4762a1bSJed Brown 284*c4762a1bSJed Brown static PetscErrorCode ComputeError(TS ts, Vec U, Vec E) 285*c4762a1bSJed Brown { 286*c4762a1bSJed Brown MPI_Comm comm; 287*c4762a1bSJed Brown DM sdm; 288*c4762a1bSJed Brown AppCtx *user; 289*c4762a1bSJed Brown const PetscScalar *u, *coords; 290*c4762a1bSJed Brown PetscScalar *e; 291*c4762a1bSJed Brown PetscReal t; 292*c4762a1bSJed Brown PetscInt dim, Np, p; 293*c4762a1bSJed Brown PetscErrorCode ierr; 294*c4762a1bSJed Brown 295*c4762a1bSJed Brown PetscFunctionBeginUser; 296*c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr); 297*c4762a1bSJed Brown ierr = TSGetDM(ts, &sdm);CHKERRQ(ierr); 298*c4762a1bSJed Brown ierr = DMGetApplicationContext(sdm, (void **) &user);CHKERRQ(ierr); 299*c4762a1bSJed Brown ierr = DMGetDimension(sdm, &dim);CHKERRQ(ierr); 300*c4762a1bSJed Brown ierr = TSGetSolveTime(ts, &t);CHKERRQ(ierr); 301*c4762a1bSJed Brown ierr = VecGetArray(E, &e);CHKERRQ(ierr); 302*c4762a1bSJed Brown ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr); 303*c4762a1bSJed Brown ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr); 304*c4762a1bSJed Brown Np /= 2*dim; 305*c4762a1bSJed Brown ierr = DMSwarmGetField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 306*c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 307*c4762a1bSJed Brown const PetscScalar *x = &u[(p*2+0)*dim]; 308*c4762a1bSJed Brown const PetscScalar *v = &u[(p*2+1)*dim]; 309*c4762a1bSJed Brown const PetscReal x0 = p+1.; 310*c4762a1bSJed Brown const PetscReal omega = PetscSqrtReal(1000./(p+1.))/x0; 311*c4762a1bSJed Brown const PetscReal xe[3] = { x0*PetscCosReal(omega*t), x0*PetscSinReal(omega*t), 0.0}; 312*c4762a1bSJed Brown const PetscReal ve[3] = {-x0*omega*PetscSinReal(omega*t), x0*omega*PetscCosReal(omega*t), 0.0}; 313*c4762a1bSJed Brown PetscInt d; 314*c4762a1bSJed Brown 315*c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 316*c4762a1bSJed Brown e[(p*2+0)*dim+d] = x[d] - xe[d]; 317*c4762a1bSJed Brown e[(p*2+1)*dim+d] = v[d] - ve[d]; 318*c4762a1bSJed Brown } 319*c4762a1bSJed Brown if (user->error) {ierr = PetscPrintf(comm, "p%D error [%.2g %.2g] sol [(%.6lf %.6lf) (%.6lf %.6lf)] exact [(%.6lf %.6lf) (%.6lf %.6lf)] energy/exact energy %g / %g\n", 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))));} 320*c4762a1bSJed Brown } 321*c4762a1bSJed Brown ierr = DMSwarmRestoreField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 322*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr); 323*c4762a1bSJed Brown ierr = VecRestoreArray(E, &e);CHKERRQ(ierr); 324*c4762a1bSJed Brown PetscFunctionReturn(0); 325*c4762a1bSJed Brown } 326*c4762a1bSJed Brown 327*c4762a1bSJed Brown int main(int argc, char **argv) 328*c4762a1bSJed Brown { 329*c4762a1bSJed Brown TS ts; 330*c4762a1bSJed Brown TSConvergedReason reason; 331*c4762a1bSJed Brown DM dm, sw; 332*c4762a1bSJed Brown Vec u; 333*c4762a1bSJed Brown IS is1, is2; 334*c4762a1bSJed Brown PetscInt *idx1, *idx2; 335*c4762a1bSJed Brown MPI_Comm comm; 336*c4762a1bSJed Brown AppCtx user; 337*c4762a1bSJed Brown const PetscScalar *endVals; 338*c4762a1bSJed Brown PetscReal ftime = .1; 339*c4762a1bSJed Brown PetscInt locSize, dim, d, Np, p, steps; 340*c4762a1bSJed Brown PetscErrorCode ierr; 341*c4762a1bSJed Brown 342*c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 343*c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 344*c4762a1bSJed Brown ierr = ProcessOptions(comm, &user);CHKERRQ(ierr); 345*c4762a1bSJed Brown 346*c4762a1bSJed Brown ierr = CreateMesh(comm, &dm, &user);CHKERRQ(ierr); 347*c4762a1bSJed Brown ierr = CreateParticles(dm, &sw, &user);CHKERRQ(ierr); 348*c4762a1bSJed Brown ierr = DMSetApplicationContext(sw, &user);CHKERRQ(ierr); 349*c4762a1bSJed Brown 350*c4762a1bSJed Brown ierr = TSCreate(comm, &ts);CHKERRQ(ierr); 351*c4762a1bSJed Brown ierr = TSSetType(ts, TSBASICSYMPLECTIC);CHKERRQ(ierr); 352*c4762a1bSJed Brown ierr = TSSetDM(ts, sw);CHKERRQ(ierr); 353*c4762a1bSJed Brown ierr = TSSetMaxTime(ts, ftime);CHKERRQ(ierr); 354*c4762a1bSJed Brown ierr = TSSetTimeStep(ts, 0.0001);CHKERRQ(ierr); 355*c4762a1bSJed Brown ierr = TSSetMaxSteps(ts, 10);CHKERRQ(ierr); 356*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 357*c4762a1bSJed Brown ierr = TSSetTime(ts, 0.0);CHKERRQ(ierr); 358*c4762a1bSJed Brown ierr = TSSetRHSFunction(ts, NULL, RHSFunctionParticles, &user);CHKERRQ(ierr); 359*c4762a1bSJed Brown 360*c4762a1bSJed Brown ierr = DMSwarmCreateGlobalVectorFromField(sw, "kinematics", &u);CHKERRQ(ierr); 361*c4762a1bSJed Brown ierr = DMGetDimension(sw, &dim);CHKERRQ(ierr); 362*c4762a1bSJed Brown ierr = VecGetLocalSize(u, &locSize);CHKERRQ(ierr); 363*c4762a1bSJed Brown Np = locSize/(2*dim); 364*c4762a1bSJed Brown ierr = PetscMalloc1(locSize/2, &idx1);CHKERRQ(ierr); 365*c4762a1bSJed Brown ierr = PetscMalloc1(locSize/2, &idx2);CHKERRQ(ierr); 366*c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 367*c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 368*c4762a1bSJed Brown idx1[p*dim+d] = (p*2+0)*dim + d; 369*c4762a1bSJed Brown idx2[p*dim+d] = (p*2+1)*dim + d; 370*c4762a1bSJed Brown } 371*c4762a1bSJed Brown } 372*c4762a1bSJed Brown ierr = ISCreateGeneral(comm, locSize/2, idx1, PETSC_OWN_POINTER, &is1);CHKERRQ(ierr); 373*c4762a1bSJed Brown ierr = ISCreateGeneral(comm, locSize/2, idx2, PETSC_OWN_POINTER, &is2);CHKERRQ(ierr); 374*c4762a1bSJed Brown ierr = TSRHSSplitSetIS(ts, "position", is1);CHKERRQ(ierr); 375*c4762a1bSJed Brown ierr = TSRHSSplitSetIS(ts, "momentum", is2);CHKERRQ(ierr); 376*c4762a1bSJed Brown ierr = ISDestroy(&is1);CHKERRQ(ierr); 377*c4762a1bSJed Brown ierr = ISDestroy(&is2);CHKERRQ(ierr); 378*c4762a1bSJed Brown ierr = TSRHSSplitSetRHSFunction(ts,"position",NULL,RHSFunction1,&user);CHKERRQ(ierr); 379*c4762a1bSJed Brown ierr = TSRHSSplitSetRHSFunction(ts,"momentum",NULL,RHSFunction2,&user);CHKERRQ(ierr); 380*c4762a1bSJed Brown 381*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 382*c4762a1bSJed Brown ierr = TSSetComputeInitialCondition(ts, InitializeSolve);CHKERRQ(ierr); 383*c4762a1bSJed Brown ierr = TSSetComputeExactError(ts, ComputeError);CHKERRQ(ierr); 384*c4762a1bSJed Brown ierr = TSComputeInitialCondition(ts, u);CHKERRQ(ierr); 385*c4762a1bSJed Brown ierr = VecViewFromOptions(u, NULL, "-init_view");CHKERRQ(ierr); 386*c4762a1bSJed Brown ierr = TSSolve(ts, u);CHKERRQ(ierr); 387*c4762a1bSJed Brown ierr = TSGetSolveTime(ts, &ftime);CHKERRQ(ierr); 388*c4762a1bSJed Brown ierr = TSGetConvergedReason(ts, &reason);CHKERRQ(ierr); 389*c4762a1bSJed Brown ierr = TSGetStepNumber(ts, &steps);CHKERRQ(ierr); 390*c4762a1bSJed Brown ierr = PetscPrintf(comm,"%s at time %g after %D steps\n", TSConvergedReasons[reason], (double) ftime, steps);CHKERRQ(ierr); 391*c4762a1bSJed Brown 392*c4762a1bSJed Brown ierr = VecGetArrayRead(u, &endVals);CHKERRQ(ierr); 393*c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 394*c4762a1bSJed Brown const PetscReal norm = DMPlex_NormD_Internal(dim, &endVals[(p*2 + 1)*dim]); 395*c4762a1bSJed Brown ierr = PetscPrintf(comm, "Particle %D initial Energy: %g Final Energy: %g\n", p, (double) (0.5*(1000./(p+1))), (double) 0.5*PetscSqr(norm));CHKERRQ(ierr); 396*c4762a1bSJed Brown } 397*c4762a1bSJed Brown ierr = VecRestoreArrayRead(u, &endVals);CHKERRQ(ierr); 398*c4762a1bSJed Brown ierr = DMSwarmDestroyGlobalVectorFromField(sw, "kinematics", &u);CHKERRQ(ierr); 399*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 400*c4762a1bSJed Brown ierr = DMDestroy(&sw);CHKERRQ(ierr); 401*c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 402*c4762a1bSJed Brown ierr = PetscFinalize(); 403*c4762a1bSJed Brown return ierr; 404*c4762a1bSJed Brown } 405*c4762a1bSJed Brown 406*c4762a1bSJed Brown 407*c4762a1bSJed Brown /*TEST 408*c4762a1bSJed Brown 409*c4762a1bSJed Brown build: 410*c4762a1bSJed Brown requires: triangle !single !complex 411*c4762a1bSJed Brown test: 412*c4762a1bSJed Brown suffix: bsi1 413*c4762a1bSJed Brown args: -dim 2 -dm_plex_box_faces 1,1 -dm_view -sw_view -ts_basicsymplectic_type 1 -ts_max_time 0.1 -ts_monitor_sp_swarm -ts_convergence_estimate -convest_num_refine 2 414*c4762a1bSJed Brown test: 415*c4762a1bSJed Brown suffix: bsi2 416*c4762a1bSJed Brown args: -dim 2 -dm_plex_box_faces 1,1 -dm_view -sw_view -ts_basicsymplectic_type 2 -ts_max_time 0.1 -ts_monitor_sp_swarm -ts_convergence_estimate -convest_num_refine 2 417*c4762a1bSJed Brown test: 418*c4762a1bSJed Brown suffix: euler 419*c4762a1bSJed Brown args: -dim 2 -dm_plex_box_faces 1,1 -dm_view -sw_view -ts_type euler -ts_max_time 0.1 -ts_monitor_sp_swarm -ts_convergence_estimate -convest_num_refine 2 420*c4762a1bSJed Brown 421*c4762a1bSJed Brown TEST*/ 422