1*fa0c1127SMark Adams static char help[] = "Landau Damping test using Vlasov-Poisson equations\n"; 2*fa0c1127SMark Adams 3*fa0c1127SMark Adams /* 4*fa0c1127SMark Adams This example solves the Vlasov-Poisson system for Landau damping (6X + 6V). 5*fa0c1127SMark Adams The system is solved using a Particle-In-Cell (PIC) method with DMSwarm for particles and DMPlex/PetscFE for the field. 6*fa0c1127SMark Adams This particular example uses the velocity mesh from DMPlexLandauCreateVelocitySpace for 3D velocity space. 7*fa0c1127SMark Adams 8*fa0c1127SMark Adams Options: 9*fa0c1127SMark Adams -particle_monitor [prefix] : Output particle data (x, v, w, E) to binary files at each output step. 10*fa0c1127SMark Adams Optional prefix for filenames (default: "particles"). 11*fa0c1127SMark Adams 12*fa0c1127SMark Adams */ 13*fa0c1127SMark Adams #include <petscts.h> 14*fa0c1127SMark Adams #include <petscdmplex.h> 15*fa0c1127SMark Adams #include <petscdmswarm.h> 16*fa0c1127SMark Adams #include <petscfe.h> 17*fa0c1127SMark Adams #include <petscds.h> 18*fa0c1127SMark Adams #include <petscbag.h> 19*fa0c1127SMark Adams #include <petscdraw.h> 20*fa0c1127SMark Adams #include <petscviewer.h> 21*fa0c1127SMark Adams #include <petsclandau.h> 22*fa0c1127SMark Adams #include <petscdmcomposite.h> 23*fa0c1127SMark Adams #include <petsc/private/dmpleximpl.h> /* For norm and dot */ 24*fa0c1127SMark Adams #include <petsc/private/petscfeimpl.h> /* For interpolation */ 25*fa0c1127SMark Adams #include <petsc/private/dmswarmimpl.h> /* For swarm debugging */ 26*fa0c1127SMark Adams #include "petscdm.h" 27*fa0c1127SMark Adams #include "petscdmlabel.h" 28*fa0c1127SMark Adams 29*fa0c1127SMark Adams PETSC_EXTERN PetscErrorCode stream(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *); 30*fa0c1127SMark Adams PETSC_EXTERN PetscErrorCode line(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *); 31*fa0c1127SMark Adams 32*fa0c1127SMark Adams typedef enum { 33*fa0c1127SMark Adams V0, 34*fa0c1127SMark Adams X0, 35*fa0c1127SMark Adams T0, 36*fa0c1127SMark Adams M0, 37*fa0c1127SMark Adams Q0, 38*fa0c1127SMark Adams PHI0, 39*fa0c1127SMark Adams POISSON, 40*fa0c1127SMark Adams VLASOV, 41*fa0c1127SMark Adams SIGMA, 42*fa0c1127SMark Adams NUM_CONSTANTS 43*fa0c1127SMark Adams } ConstantType; 44*fa0c1127SMark Adams typedef struct { 45*fa0c1127SMark Adams PetscScalar v0; /* Velocity scale, often the thermal velocity */ 46*fa0c1127SMark Adams PetscScalar t0; /* Time scale */ 47*fa0c1127SMark Adams PetscScalar x0; /* Space scale */ 48*fa0c1127SMark Adams PetscScalar m0; /* Mass scale */ 49*fa0c1127SMark Adams PetscScalar q0; /* Charge scale */ 50*fa0c1127SMark Adams PetscScalar kb; 51*fa0c1127SMark Adams PetscScalar epsi0; 52*fa0c1127SMark Adams PetscScalar phi0; /* Potential scale */ 53*fa0c1127SMark Adams PetscScalar poissonNumber; /* Non-Dimensional Poisson Number */ 54*fa0c1127SMark Adams PetscScalar vlasovNumber; /* Non-Dimensional Vlasov Number */ 55*fa0c1127SMark Adams PetscReal sigma; /* Nondimensional charge per length in x */ 56*fa0c1127SMark Adams } Parameter; 57*fa0c1127SMark Adams 58*fa0c1127SMark Adams typedef struct { 59*fa0c1127SMark Adams PetscBag bag; // Problem parameters 60*fa0c1127SMark Adams PetscBool error; // Flag for printing the error 61*fa0c1127SMark Adams PetscBool efield_monitor; // Flag to show electric field monitor 62*fa0c1127SMark Adams PetscBool moment_monitor; // Flag to show distribution moment monitor 63*fa0c1127SMark Adams char particle_monitor_prefix[PETSC_MAX_PATH_LEN]; 64*fa0c1127SMark Adams PetscBool particle_monitor; // Flag to output particle data 65*fa0c1127SMark Adams PetscInt ostep; // Print the energy at each ostep time steps 66*fa0c1127SMark Adams PetscInt numParticles; 67*fa0c1127SMark Adams PetscReal timeScale; /* Nondimensionalizing time scale */ 68*fa0c1127SMark Adams PetscReal charges[2]; /* The charges of each species */ 69*fa0c1127SMark Adams PetscReal masses[2]; /* The masses of each species */ 70*fa0c1127SMark Adams PetscReal thermal_energy[2]; /* Thermal Energy (used to get other constants)*/ 71*fa0c1127SMark Adams PetscReal cosine_coefficients[2]; /*(alpha, k)*/ 72*fa0c1127SMark Adams PetscReal totalWeight; 73*fa0c1127SMark Adams PetscReal stepSize; 74*fa0c1127SMark Adams PetscInt steps; 75*fa0c1127SMark Adams PetscReal initVel; 76*fa0c1127SMark Adams SNES snes; // EM solver 77*fa0c1127SMark Adams DM dmPot; // The DM for potential 78*fa0c1127SMark Adams Mat M; // The finite element mass matrix for potential 79*fa0c1127SMark Adams PetscFEGeom *fegeom; // Geometric information for the DM cells 80*fa0c1127SMark Adams PetscBool validE; // Flag to indicate E-field in swarm is valid 81*fa0c1127SMark Adams PetscReal drawlgEmin; // The minimum lg(E) to plot 82*fa0c1127SMark Adams PetscDrawLG drawlgE; // Logarithm of maximum electric field 83*fa0c1127SMark Adams DM swarm; 84*fa0c1127SMark Adams PetscBool checkweights; 85*fa0c1127SMark Adams PetscInt checkVRes; // Flag to check/output velocity residuals for nightly tests 86*fa0c1127SMark Adams DM landau_pack; 87*fa0c1127SMark Adams PetscBool use_landau_velocity_space; 88*fa0c1127SMark Adams PetscLogEvent RhsXEvent, RhsVEvent, ESolveEvent, ETabEvent; 89*fa0c1127SMark Adams } AppCtx; 90*fa0c1127SMark Adams 91*fa0c1127SMark Adams static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 92*fa0c1127SMark Adams { 93*fa0c1127SMark Adams PetscFunctionBeginUser; 94*fa0c1127SMark Adams PetscInt d = 2; 95*fa0c1127SMark Adams PetscInt maxSpecies = 2; 96*fa0c1127SMark Adams options->error = PETSC_FALSE; 97*fa0c1127SMark Adams options->efield_monitor = PETSC_FALSE; 98*fa0c1127SMark Adams options->moment_monitor = PETSC_FALSE; 99*fa0c1127SMark Adams options->particle_monitor = PETSC_FALSE; 100*fa0c1127SMark Adams options->ostep = 100; 101*fa0c1127SMark Adams options->timeScale = 2.0e-14; 102*fa0c1127SMark Adams options->charges[0] = -1.0; 103*fa0c1127SMark Adams options->charges[1] = 1.0; 104*fa0c1127SMark Adams options->masses[0] = 1.0; 105*fa0c1127SMark Adams options->masses[1] = 1000.0; 106*fa0c1127SMark Adams options->thermal_energy[0] = 1.0; 107*fa0c1127SMark Adams options->thermal_energy[1] = 1.0; 108*fa0c1127SMark Adams options->cosine_coefficients[0] = 0.01; 109*fa0c1127SMark Adams options->cosine_coefficients[1] = 0.5; 110*fa0c1127SMark Adams options->initVel = 1; 111*fa0c1127SMark Adams options->totalWeight = 1.0; 112*fa0c1127SMark Adams options->validE = PETSC_FALSE; 113*fa0c1127SMark Adams options->drawlgEmin = -6; 114*fa0c1127SMark Adams options->drawlgE = NULL; 115*fa0c1127SMark Adams options->snes = NULL; 116*fa0c1127SMark Adams options->dmPot = NULL; 117*fa0c1127SMark Adams options->M = NULL; 118*fa0c1127SMark Adams options->numParticles = 32768; 119*fa0c1127SMark Adams options->checkweights = PETSC_FALSE; 120*fa0c1127SMark Adams options->checkVRes = 0; 121*fa0c1127SMark Adams options->landau_pack = NULL; 122*fa0c1127SMark Adams options->use_landau_velocity_space = PETSC_FALSE; 123*fa0c1127SMark Adams 124*fa0c1127SMark Adams PetscOptionsBegin(comm, "", "Landau Damping options", "DMSWARM"); 125*fa0c1127SMark Adams PetscCall(PetscOptionsBool("-error", "Flag to print the error", "ex3.c", options->error, &options->error, NULL)); 126*fa0c1127SMark Adams PetscCall(PetscOptionsBool("-efield_monitor", "Flag to plot log(max E) over time", "ex3.c", options->efield_monitor, &options->efield_monitor, NULL)); 127*fa0c1127SMark Adams PetscCall(PetscOptionsReal("-efield_min_monitor", "Minimum E field to plot", "ex3.c", options->drawlgEmin, &options->drawlgEmin, NULL)); 128*fa0c1127SMark Adams PetscCall(PetscOptionsBool("-moments_monitor", "Flag to show moments table", "ex3.c", options->moment_monitor, &options->moment_monitor, NULL)); 129*fa0c1127SMark Adams PetscCall(PetscOptionsString("-particle_monitor", "Prefix for particle data files", "ex3.c", options->particle_monitor_prefix, options->particle_monitor_prefix, sizeof(options->particle_monitor_prefix), &options->particle_monitor)); 130*fa0c1127SMark Adams PetscCall(PetscOptionsBool("-check_weights", "Ensure all particle weights are positive", "ex3.c", options->checkweights, &options->checkweights, NULL)); 131*fa0c1127SMark Adams PetscCall(PetscOptionsInt("-check_vel_res", "Check particle velocity residuals for nightly tests", "ex3.c", options->checkVRes, &options->checkVRes, NULL)); 132*fa0c1127SMark Adams PetscCall(PetscOptionsInt("-output_step", "Number of time steps between output", "ex3.c", options->ostep, &options->ostep, NULL)); 133*fa0c1127SMark Adams PetscCall(PetscOptionsReal("-timeScale", "Nondimensionalizing time scale", "ex3.c", options->timeScale, &options->timeScale, NULL)); 134*fa0c1127SMark Adams PetscCall(PetscOptionsReal("-initial_velocity", "Initial velocity of perturbed particle", "ex3.c", options->initVel, &options->initVel, NULL)); 135*fa0c1127SMark Adams PetscCall(PetscOptionsReal("-total_weight", "Total weight of all particles", "ex3.c", options->totalWeight, &options->totalWeight, NULL)); 136*fa0c1127SMark Adams PetscCall(PetscOptionsRealArray("-cosine_coefficients", "Amplitude and frequency of cosine equation used in initialization", "ex3.c", options->cosine_coefficients, &d, NULL)); 137*fa0c1127SMark Adams PetscCall(PetscOptionsRealArray("-charges", "Species charges", "ex3.c", options->charges, &maxSpecies, NULL)); 138*fa0c1127SMark Adams PetscCall(PetscOptionsBool("-use_landau_velocity_space", "Use Landau velocity space", "ex3.c", options->use_landau_velocity_space, &options->use_landau_velocity_space, NULL)); 139*fa0c1127SMark Adams PetscOptionsEnd(); 140*fa0c1127SMark Adams 141*fa0c1127SMark Adams PetscCall(PetscLogEventRegister("RhsX", TS_CLASSID, &options->RhsXEvent)); 142*fa0c1127SMark Adams PetscCall(PetscLogEventRegister("RhsV", TS_CLASSID, &options->RhsVEvent)); 143*fa0c1127SMark Adams PetscCall(PetscLogEventRegister("ESolve", TS_CLASSID, &options->ESolveEvent)); 144*fa0c1127SMark Adams PetscCall(PetscLogEventRegister("ETab", TS_CLASSID, &options->ETabEvent)); 145*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 146*fa0c1127SMark Adams } 147*fa0c1127SMark Adams 148*fa0c1127SMark Adams static PetscErrorCode SetupContext(DM dm, DM sw, AppCtx *user) 149*fa0c1127SMark Adams { 150*fa0c1127SMark Adams MPI_Comm comm; 151*fa0c1127SMark Adams 152*fa0c1127SMark Adams PetscFunctionBeginUser; 153*fa0c1127SMark Adams PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 154*fa0c1127SMark Adams if (user->efield_monitor) { 155*fa0c1127SMark Adams PetscDraw draw; 156*fa0c1127SMark Adams PetscDrawAxis axis; 157*fa0c1127SMark Adams 158*fa0c1127SMark Adams PetscCall(PetscDrawCreate(comm, NULL, "Max Electric Field", 0, 300, 400, 300, &draw)); 159*fa0c1127SMark Adams PetscCall(PetscDrawSetSave(draw, "ex3_Efield")); 160*fa0c1127SMark Adams PetscCall(PetscDrawSetFromOptions(draw)); 161*fa0c1127SMark Adams PetscCall(PetscDrawLGCreate(draw, 1, &user->drawlgE)); 162*fa0c1127SMark Adams PetscCall(PetscDrawDestroy(&draw)); 163*fa0c1127SMark Adams PetscCall(PetscDrawLGGetAxis(user->drawlgE, &axis)); 164*fa0c1127SMark Adams PetscCall(PetscDrawAxisSetLabels(axis, "Electron Electric Field", "time", "E_max")); 165*fa0c1127SMark Adams PetscCall(PetscDrawLGSetLimits(user->drawlgE, 0., user->steps * user->stepSize, user->drawlgEmin, 0.)); 166*fa0c1127SMark Adams PetscCall(PetscDrawLGSetUseMarkers(user->drawlgE, PETSC_FALSE)); 167*fa0c1127SMark Adams } 168*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 169*fa0c1127SMark Adams } 170*fa0c1127SMark Adams 171*fa0c1127SMark Adams static PetscErrorCode DestroyContext(AppCtx *user) 172*fa0c1127SMark Adams { 173*fa0c1127SMark Adams PetscFunctionBeginUser; 174*fa0c1127SMark Adams if (user->landau_pack) PetscCall(DMPlexLandauDestroyVelocitySpace(&user->landau_pack)); 175*fa0c1127SMark Adams PetscCall(PetscDrawLGDestroy(&user->drawlgE)); 176*fa0c1127SMark Adams PetscCall(PetscBagDestroy(&user->bag)); 177*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 178*fa0c1127SMark Adams } 179*fa0c1127SMark Adams 180*fa0c1127SMark Adams static PetscErrorCode CheckNonNegativeWeights(DM sw, AppCtx *user) 181*fa0c1127SMark Adams { 182*fa0c1127SMark Adams const PetscScalar *w; 183*fa0c1127SMark Adams PetscInt Np; 184*fa0c1127SMark Adams 185*fa0c1127SMark Adams PetscFunctionBeginUser; 186*fa0c1127SMark Adams if (!user->checkweights) PetscFunctionReturn(PETSC_SUCCESS); 187*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&w)); 188*fa0c1127SMark Adams PetscCall(DMSwarmGetLocalSize(sw, &Np)); 189*fa0c1127SMark Adams for (PetscInt p = 0; p < Np; ++p) PetscCheck(w[p] >= 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Particle %" PetscInt_FMT " has negative weight %g", p, w[p]); 190*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&w)); 191*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 192*fa0c1127SMark Adams } 193*fa0c1127SMark Adams 194*fa0c1127SMark Adams static void f0_Dirichlet(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 195*fa0c1127SMark Adams { 196*fa0c1127SMark Adams for (PetscInt d = 0; d < dim; ++d) f0[0] += 0.5 * PetscSqr(u_x[d]); 197*fa0c1127SMark Adams } 198*fa0c1127SMark Adams 199*fa0c1127SMark Adams static PetscErrorCode computeFieldEnergy(DM dm, Vec u, PetscReal *En) 200*fa0c1127SMark Adams { 201*fa0c1127SMark Adams PetscDS ds; 202*fa0c1127SMark Adams const PetscInt field = 0; 203*fa0c1127SMark Adams PetscInt Nf; 204*fa0c1127SMark Adams void *ctx; 205*fa0c1127SMark Adams 206*fa0c1127SMark Adams PetscFunctionBegin; 207*fa0c1127SMark Adams PetscCall(DMGetApplicationContext(dm, &ctx)); 208*fa0c1127SMark Adams PetscCall(DMGetDS(dm, &ds)); 209*fa0c1127SMark Adams PetscCall(PetscDSGetNumFields(ds, &Nf)); 210*fa0c1127SMark Adams PetscCheck(Nf == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "We currently only support 1 field, not %" PetscInt_FMT, Nf); 211*fa0c1127SMark Adams PetscCall(PetscDSSetObjective(ds, field, &f0_Dirichlet)); 212*fa0c1127SMark Adams PetscCall(DMPlexComputeIntegralFEM(dm, u, En, ctx)); 213*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 214*fa0c1127SMark Adams } 215*fa0c1127SMark Adams 216*fa0c1127SMark Adams static void f0_grad_phi2(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 217*fa0c1127SMark Adams { 218*fa0c1127SMark Adams f0[0] = 0.; 219*fa0c1127SMark Adams for (PetscInt d = 0; d < dim; ++d) f0[0] += PetscSqr(u_x[uOff_x[0] + d]); // + d * dim cause segfault 220*fa0c1127SMark Adams } 221*fa0c1127SMark Adams 222*fa0c1127SMark Adams static PetscErrorCode MonitorEField(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) 223*fa0c1127SMark Adams { 224*fa0c1127SMark Adams AppCtx *user = (AppCtx *)ctx; 225*fa0c1127SMark Adams DM sw; 226*fa0c1127SMark Adams PetscScalar intESq; 227*fa0c1127SMark Adams PetscReal *E, *x, *weight; 228*fa0c1127SMark Adams PetscReal Enorm = 0., lgEnorm, lgEmax, sum = 0., Emax = 0., chargesum = 0.; 229*fa0c1127SMark Adams PetscReal pmoments[16]; /* \int f, \int v f, \int v^2 f */ 230*fa0c1127SMark Adams PetscInt *species, dim, Np, gNp; 231*fa0c1127SMark Adams MPI_Comm comm; 232*fa0c1127SMark Adams 233*fa0c1127SMark Adams PetscFunctionBeginUser; 234*fa0c1127SMark Adams if (step < 0 || !user->validE) PetscFunctionReturn(PETSC_SUCCESS); 235*fa0c1127SMark Adams PetscCall(PetscObjectGetComm((PetscObject)ts, &comm)); 236*fa0c1127SMark Adams PetscCall(TSGetDM(ts, &sw)); 237*fa0c1127SMark Adams PetscCall(DMGetDimension(sw, &dim)); 238*fa0c1127SMark Adams PetscCall(DMSwarmGetLocalSize(sw, &Np)); 239*fa0c1127SMark Adams PetscCall(DMSwarmGetSize(sw, &gNp)); 240*fa0c1127SMark Adams PetscCall(DMSwarmSortGetAccess(sw)); 241*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x)); 242*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, "E_field", NULL, NULL, (void **)&E)); 243*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 244*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&weight)); 245*fa0c1127SMark Adams 246*fa0c1127SMark Adams for (PetscInt p = 0; p < Np; ++p) { 247*fa0c1127SMark Adams PetscReal Emag = 0.; 248*fa0c1127SMark Adams for (PetscInt d = 0; d < dim; ++d) { 249*fa0c1127SMark Adams PetscReal temp = PetscAbsReal(E[p * dim + d]); 250*fa0c1127SMark Adams if (temp > Emax) Emax = temp; 251*fa0c1127SMark Adams Emag += PetscSqr(E[p * dim + d]); 252*fa0c1127SMark Adams } 253*fa0c1127SMark Adams Enorm += PetscSqrtReal(Emag); 254*fa0c1127SMark Adams sum += E[p * dim]; 255*fa0c1127SMark Adams chargesum += user->charges[0] * weight[p]; 256*fa0c1127SMark Adams } 257*fa0c1127SMark Adams PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &Emax, 1, MPIU_REAL, MPIU_MAX, comm)); 258*fa0c1127SMark Adams PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &Enorm, 1, MPIU_REAL, MPIU_SUM, comm)); 259*fa0c1127SMark Adams PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &sum, 1, MPIU_REAL, MPIU_SUM, comm)); 260*fa0c1127SMark Adams PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &chargesum, 1, MPIU_REAL, MPIU_SUM, comm)); 261*fa0c1127SMark Adams lgEnorm = Enorm != 0 ? PetscLog10Real(Enorm) : -16.; 262*fa0c1127SMark Adams lgEmax = Emax != 0 ? PetscLog10Real(Emax) : user->drawlgEmin; 263*fa0c1127SMark Adams if (lgEmax < user->drawlgEmin) lgEmax = user->drawlgEmin; 264*fa0c1127SMark Adams 265*fa0c1127SMark Adams PetscDS ds; 266*fa0c1127SMark Adams Vec phi; 267*fa0c1127SMark Adams 268*fa0c1127SMark Adams PetscCall(DMGetNamedGlobalVector(user->dmPot, "phi", &phi)); 269*fa0c1127SMark Adams PetscCall(DMGetDS(user->dmPot, &ds)); 270*fa0c1127SMark Adams PetscCall(PetscDSSetObjective(ds, 0, &f0_grad_phi2)); 271*fa0c1127SMark Adams PetscCall(DMPlexComputeIntegralFEM(user->dmPot, phi, &intESq, user)); 272*fa0c1127SMark Adams PetscCall(DMRestoreNamedGlobalVector(user->dmPot, "phi", &phi)); 273*fa0c1127SMark Adams 274*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x)); 275*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight)); 276*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, "E_field", NULL, NULL, (void **)&E)); 277*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 278*fa0c1127SMark Adams PetscCall(PetscDrawLGAddPoint(user->drawlgE, &t, &lgEmax)); 279*fa0c1127SMark Adams PetscCall(PetscDrawLGDraw(user->drawlgE)); 280*fa0c1127SMark Adams PetscDraw draw; 281*fa0c1127SMark Adams PetscCall(PetscDrawLGGetDraw(user->drawlgE, &draw)); 282*fa0c1127SMark Adams PetscCall(PetscDrawSave(draw)); 283*fa0c1127SMark Adams 284*fa0c1127SMark Adams PetscCall(DMSwarmComputeMoments(sw, "velocity", "w_q", pmoments)); 285*fa0c1127SMark Adams PetscCall(PetscPrintf(comm, "E: %f\t%+e\t%e\t%f\t%20.15e\t%f\t%f\t%f\t%20.15e\t%20.15e\t%20.15e\t%" PetscInt_FMT "\t(%" PetscInt_FMT ")\n", (double)t, (double)sum, (double)Enorm, (double)lgEnorm, (double)Emax, (double)lgEmax, (double)chargesum, (double)pmoments[0], (double)pmoments[1], (double)pmoments[1 + dim], (double)PetscSqrtReal(intESq), gNp, step)); 286*fa0c1127SMark Adams PetscCall(DMViewFromOptions(sw, NULL, "-sw_efield_view")); 287*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 288*fa0c1127SMark Adams } 289*fa0c1127SMark Adams 290*fa0c1127SMark Adams static PetscErrorCode MonitorMoments(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) 291*fa0c1127SMark Adams { 292*fa0c1127SMark Adams DM sw; 293*fa0c1127SMark Adams PetscReal pmoments[16]; /* \int f, \int v f, \int v^2 f */ 294*fa0c1127SMark Adams 295*fa0c1127SMark Adams PetscFunctionBeginUser; 296*fa0c1127SMark Adams if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); 297*fa0c1127SMark Adams PetscCall(TSGetDM(ts, &sw)); 298*fa0c1127SMark Adams 299*fa0c1127SMark Adams PetscCall(DMSwarmComputeMoments(sw, "velocity", "w_q", pmoments)); 300*fa0c1127SMark Adams 301*fa0c1127SMark Adams PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%f\t%f\t%f\t%f\n", (double)t, (double)pmoments[0], (double)pmoments[1], (double)pmoments[3])); 302*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 303*fa0c1127SMark Adams } 304*fa0c1127SMark Adams 305*fa0c1127SMark Adams static PetscErrorCode MonitorParticles(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) 306*fa0c1127SMark Adams { 307*fa0c1127SMark Adams AppCtx *user = (AppCtx *)ctx; 308*fa0c1127SMark Adams DM sw; 309*fa0c1127SMark Adams PetscInt Np, dim; 310*fa0c1127SMark Adams const PetscReal *x, *v, *E; 311*fa0c1127SMark Adams const PetscScalar *w; 312*fa0c1127SMark Adams PetscViewer viewer; 313*fa0c1127SMark Adams char filename[64]; 314*fa0c1127SMark Adams MPI_Comm comm; 315*fa0c1127SMark Adams 316*fa0c1127SMark Adams PetscFunctionBeginUser; 317*fa0c1127SMark Adams if (step < 0 || step % user->ostep != 0) PetscFunctionReturn(PETSC_SUCCESS); 318*fa0c1127SMark Adams PetscCall(TSGetDM(ts, &sw)); 319*fa0c1127SMark Adams PetscCall(PetscObjectGetComm((PetscObject)ts, &comm)); 320*fa0c1127SMark Adams PetscCall(DMGetDimension(sw, &dim)); 321*fa0c1127SMark Adams PetscCall(DMSwarmGetLocalSize(sw, &Np)); 322*fa0c1127SMark Adams 323*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x)); 324*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&v)); 325*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&w)); 326*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, "E_field", NULL, NULL, (void **)&E)); 327*fa0c1127SMark Adams 328*fa0c1127SMark Adams if (user->particle_monitor_prefix[0]) { 329*fa0c1127SMark Adams PetscCall(PetscSNPrintf(filename, sizeof(filename), "%s_step_%d.bin", user->particle_monitor_prefix, (int)step)); 330*fa0c1127SMark Adams } else { 331*fa0c1127SMark Adams PetscCall(PetscSNPrintf(filename, sizeof(filename), "particles_step_%d.bin", (int)step)); 332*fa0c1127SMark Adams } 333*fa0c1127SMark Adams PetscCall(PetscViewerBinaryOpen(comm, filename, FILE_MODE_WRITE, &viewer)); 334*fa0c1127SMark Adams 335*fa0c1127SMark Adams { 336*fa0c1127SMark Adams Vec vx, vv, vw, vE; 337*fa0c1127SMark Adams PetscCall(VecCreateMPIWithArray(comm, dim, Np * dim, PETSC_DECIDE, x, &vx)); 338*fa0c1127SMark Adams PetscCall(VecCreateMPIWithArray(comm, dim, Np * dim, PETSC_DECIDE, v, &vv)); 339*fa0c1127SMark Adams PetscCall(VecCreateMPIWithArray(comm, 1, Np, PETSC_DECIDE, w, &vw)); 340*fa0c1127SMark Adams PetscCall(VecCreateMPIWithArray(comm, dim, Np * dim, PETSC_DECIDE, E, &vE)); 341*fa0c1127SMark Adams 342*fa0c1127SMark Adams PetscCall(VecView(vx, viewer)); 343*fa0c1127SMark Adams PetscCall(VecView(vv, viewer)); 344*fa0c1127SMark Adams PetscCall(VecView(vw, viewer)); 345*fa0c1127SMark Adams PetscCall(VecView(vE, viewer)); 346*fa0c1127SMark Adams 347*fa0c1127SMark Adams PetscCall(VecDestroy(&vx)); 348*fa0c1127SMark Adams PetscCall(VecDestroy(&vv)); 349*fa0c1127SMark Adams PetscCall(VecDestroy(&vw)); 350*fa0c1127SMark Adams PetscCall(VecDestroy(&vE)); 351*fa0c1127SMark Adams } 352*fa0c1127SMark Adams PetscCall(PetscViewerDestroy(&viewer)); 353*fa0c1127SMark Adams 354*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x)); 355*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&v)); 356*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&w)); 357*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, "E_field", NULL, NULL, (void **)&E)); 358*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 359*fa0c1127SMark Adams } 360*fa0c1127SMark Adams 361*fa0c1127SMark Adams PetscErrorCode PetscPDFPertubedConstant2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 362*fa0c1127SMark Adams { 363*fa0c1127SMark Adams p[0] = (1. + 0.01 * PetscCosReal(0.5 * x[0])) / (4. * PETSC_PI * 4. * PETSC_PI); 364*fa0c1127SMark Adams return PETSC_SUCCESS; 365*fa0c1127SMark Adams } 366*fa0c1127SMark Adams PetscErrorCode PetscPDFPertubedConstant1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 367*fa0c1127SMark Adams { 368*fa0c1127SMark Adams p[0] = (1. + 0.01 * PetscCosReal(0.5 * x[0])) / (2 * PETSC_PI); 369*fa0c1127SMark Adams return PETSC_SUCCESS; 370*fa0c1127SMark Adams } 371*fa0c1127SMark Adams 372*fa0c1127SMark Adams PetscErrorCode PetscPDFCosine1D(const PetscReal x[], const PetscReal scale[], PetscReal p[]) 373*fa0c1127SMark Adams { 374*fa0c1127SMark Adams const PetscReal alpha = scale ? scale[0] : 0.0; 375*fa0c1127SMark Adams const PetscReal k = scale ? scale[1] : 1.; 376*fa0c1127SMark Adams p[0] = (1 + alpha * PetscCosReal(k * x[0])); 377*fa0c1127SMark Adams return PETSC_SUCCESS; 378*fa0c1127SMark Adams } 379*fa0c1127SMark Adams 380*fa0c1127SMark Adams PetscErrorCode PetscPDFCosine2D(const PetscReal x[], const PetscReal scale[], PetscReal p[]) 381*fa0c1127SMark Adams { 382*fa0c1127SMark Adams const PetscReal alpha = scale ? scale[0] : 0.; 383*fa0c1127SMark Adams const PetscReal k = scale ? scale[1] : 1.; 384*fa0c1127SMark Adams p[0] = (1 + alpha * PetscCosReal(k * x[0])); 385*fa0c1127SMark Adams return PETSC_SUCCESS; 386*fa0c1127SMark Adams } 387*fa0c1127SMark Adams 388*fa0c1127SMark Adams PetscErrorCode PetscPDFCosine3D(const PetscReal x[], const PetscReal scale[], PetscReal p[]) 389*fa0c1127SMark Adams { 390*fa0c1127SMark Adams const PetscReal alpha = scale ? scale[0] : 0.; 391*fa0c1127SMark Adams const PetscReal k = scale ? scale[1] : 1.; 392*fa0c1127SMark Adams p[0] = (1 + alpha * PetscCosReal(k * x[0])); 393*fa0c1127SMark Adams return PETSC_SUCCESS; 394*fa0c1127SMark Adams } 395*fa0c1127SMark Adams 396*fa0c1127SMark Adams static PetscErrorCode CreateVelocityDM(DM sw, DM *vdm) 397*fa0c1127SMark Adams { 398*fa0c1127SMark Adams PetscFE fe; 399*fa0c1127SMark Adams DMPolytopeType ct; 400*fa0c1127SMark Adams PetscInt dim, cStart; 401*fa0c1127SMark Adams const char *prefix = "v"; 402*fa0c1127SMark Adams AppCtx *user; 403*fa0c1127SMark Adams 404*fa0c1127SMark Adams PetscFunctionBegin; 405*fa0c1127SMark Adams PetscCall(DMGetDimension(sw, &dim)); 406*fa0c1127SMark Adams PetscCall(DMGetApplicationContext(sw, &user)); 407*fa0c1127SMark Adams if (dim == 3 && user->use_landau_velocity_space) { 408*fa0c1127SMark Adams LandauCtx *ctx; 409*fa0c1127SMark Adams Vec X; 410*fa0c1127SMark Adams Mat J; 411*fa0c1127SMark Adams 412*fa0c1127SMark Adams PetscCall(DMPlexLandauCreateVelocitySpace(PETSC_COMM_SELF, dim, prefix, &X, &J, &user->landau_pack)); 413*fa0c1127SMark Adams PetscCall(DMGetApplicationContext(user->landau_pack, &ctx)); 414*fa0c1127SMark Adams *vdm = ctx->plex[0]; 415*fa0c1127SMark Adams PetscCall(PetscObjectReference((PetscObject)*vdm)); 416*fa0c1127SMark Adams PetscCall(VecDestroy(&X)); 417*fa0c1127SMark Adams PetscCall(PetscObjectSetName((PetscObject)*vdm, "velocity")); 418*fa0c1127SMark Adams } else { 419*fa0c1127SMark Adams PetscCall(DMCreate(PETSC_COMM_SELF, vdm)); 420*fa0c1127SMark Adams PetscCall(DMSetType(*vdm, DMPLEX)); 421*fa0c1127SMark Adams PetscCall(DMPlexSetOptionsPrefix(*vdm, prefix)); 422*fa0c1127SMark Adams PetscCall(DMSetFromOptions(*vdm)); 423*fa0c1127SMark Adams PetscCall(PetscObjectSetName((PetscObject)*vdm, "velocity")); 424*fa0c1127SMark Adams } 425*fa0c1127SMark Adams PetscCall(DMViewFromOptions(*vdm, NULL, "-dm_view")); 426*fa0c1127SMark Adams 427*fa0c1127SMark Adams PetscCall(DMGetDimension(*vdm, &dim)); 428*fa0c1127SMark Adams PetscCall(DMPlexGetHeightStratum(*vdm, 0, &cStart, NULL)); 429*fa0c1127SMark Adams PetscCall(DMPlexGetCellType(*vdm, cStart, &ct)); 430*fa0c1127SMark Adams PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, 1, ct, prefix, PETSC_DETERMINE, &fe)); 431*fa0c1127SMark Adams PetscCall(PetscObjectSetName((PetscObject)fe, "distribution")); 432*fa0c1127SMark Adams PetscCall(DMSetField(*vdm, 0, NULL, (PetscObject)fe)); 433*fa0c1127SMark Adams PetscCall(DMCreateDS(*vdm)); 434*fa0c1127SMark Adams PetscCall(PetscFEDestroy(&fe)); 435*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 436*fa0c1127SMark Adams } 437*fa0c1127SMark Adams 438*fa0c1127SMark Adams /* 439*fa0c1127SMark Adams InitializeParticles_Centroid - Initialize a regular grid of particles. 440*fa0c1127SMark Adams 441*fa0c1127SMark Adams Input Parameters: 442*fa0c1127SMark Adams + sw - The `DMSWARM` 443*fa0c1127SMark Adams - force1D - Treat the spatial domain as 1D 444*fa0c1127SMark Adams 445*fa0c1127SMark Adams Notes: 446*fa0c1127SMark Adams This functions sets the species, cellid, spatial coordinate, and velocity fields for all particles. 447*fa0c1127SMark Adams 448*fa0c1127SMark Adams It places one particle in the centroid of each cell in the implicit tensor product of the spatial 449*fa0c1127SMark Adams and velocity meshes. 450*fa0c1127SMark Adams */ 451*fa0c1127SMark Adams static PetscErrorCode InitializeParticles_Centroid(DM sw) 452*fa0c1127SMark Adams { 453*fa0c1127SMark Adams DM_Swarm *swarm = (DM_Swarm *)sw->data; 454*fa0c1127SMark Adams DMSwarmCellDM celldm; 455*fa0c1127SMark Adams DM xdm, vdm; 456*fa0c1127SMark Adams PetscReal vmin[3], vmax[3]; 457*fa0c1127SMark Adams PetscReal *x, *v; 458*fa0c1127SMark Adams PetscInt *species, *cellid; 459*fa0c1127SMark Adams PetscInt dim, xcStart, xcEnd, vcStart, vcEnd, Ns, Np, Npc, debug; 460*fa0c1127SMark Adams PetscBool flg; 461*fa0c1127SMark Adams MPI_Comm comm; 462*fa0c1127SMark Adams const char *cellidname; 463*fa0c1127SMark Adams 464*fa0c1127SMark Adams PetscFunctionBegin; 465*fa0c1127SMark Adams PetscCall(PetscObjectGetComm((PetscObject)sw, &comm)); 466*fa0c1127SMark Adams 467*fa0c1127SMark Adams PetscOptionsBegin(comm, "", "DMSwarm Options", "DMSWARM"); 468*fa0c1127SMark Adams PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 469*fa0c1127SMark Adams PetscCall(PetscOptionsInt("-dm_swarm_num_species", "The number of species", "DMSwarmSetNumSpecies", Ns, &Ns, &flg)); 470*fa0c1127SMark Adams if (flg) PetscCall(DMSwarmSetNumSpecies(sw, Ns)); 471*fa0c1127SMark Adams PetscCall(PetscOptionsBoundedInt("-dm_swarm_print_coords", "Debug output level for particle coordinate computations", "InitializeParticles", 0, &swarm->printCoords, NULL, 0)); 472*fa0c1127SMark Adams PetscCall(PetscOptionsBoundedInt("-dm_swarm_print_weights", "Debug output level for particle weight computations", "InitializeWeights", 0, &swarm->printWeights, NULL, 0)); 473*fa0c1127SMark Adams PetscOptionsEnd(); 474*fa0c1127SMark Adams debug = swarm->printCoords; 475*fa0c1127SMark Adams 476*fa0c1127SMark Adams PetscCall(DMGetDimension(sw, &dim)); 477*fa0c1127SMark Adams PetscCall(DMSwarmGetCellDM(sw, &xdm)); 478*fa0c1127SMark Adams PetscCall(DMPlexGetHeightStratum(xdm, 0, &xcStart, &xcEnd)); 479*fa0c1127SMark Adams 480*fa0c1127SMark Adams PetscCall(DMSwarmGetCellDMByName(sw, "velocity", &celldm)); 481*fa0c1127SMark Adams PetscCall(DMSwarmCellDMGetDM(celldm, &vdm)); 482*fa0c1127SMark Adams PetscCall(DMPlexGetHeightStratum(vdm, 0, &vcStart, &vcEnd)); 483*fa0c1127SMark Adams 484*fa0c1127SMark Adams // One particle per centroid on the tensor product grid 485*fa0c1127SMark Adams Npc = (vcEnd - vcStart) * Ns; 486*fa0c1127SMark Adams Np = (xcEnd - xcStart) * Npc; 487*fa0c1127SMark Adams PetscCall(DMSwarmSetLocalSizes(sw, Np, 0)); 488*fa0c1127SMark Adams if (debug) { 489*fa0c1127SMark Adams PetscInt gNp, gNc, Nc = xcEnd - xcStart; 490*fa0c1127SMark Adams PetscCallMPI(MPIU_Allreduce(&Np, &gNp, 1, MPIU_INT, MPIU_SUM, comm)); 491*fa0c1127SMark Adams PetscCall(PetscPrintf(comm, "Global Np = %" PetscInt_FMT "\n", gNp)); 492*fa0c1127SMark Adams PetscCallMPI(MPIU_Allreduce(&Nc, &gNc, 1, MPIU_INT, MPIU_SUM, comm)); 493*fa0c1127SMark Adams PetscCall(PetscPrintf(comm, "Global X-cells = %" PetscInt_FMT "\n", gNc)); 494*fa0c1127SMark Adams PetscCall(PetscPrintf(comm, "Global V-cells = %" PetscInt_FMT "\n", vcEnd - vcStart)); 495*fa0c1127SMark Adams } 496*fa0c1127SMark Adams 497*fa0c1127SMark Adams // Set species and cellid 498*fa0c1127SMark Adams PetscCall(DMSwarmGetCellDMActive(sw, &celldm)); 499*fa0c1127SMark Adams PetscCall(DMSwarmCellDMGetCellID(celldm, &cellidname)); 500*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 501*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, cellidname, NULL, NULL, (void **)&cellid)); 502*fa0c1127SMark Adams for (PetscInt c = 0, p = 0; c < xcEnd - xcStart; ++c) { 503*fa0c1127SMark Adams for (PetscInt s = 0; s < Ns; ++s) { 504*fa0c1127SMark Adams for (PetscInt q = 0; q < Npc / Ns; ++q, ++p) { 505*fa0c1127SMark Adams species[p] = s; 506*fa0c1127SMark Adams cellid[p] = c; 507*fa0c1127SMark Adams } 508*fa0c1127SMark Adams } 509*fa0c1127SMark Adams } 510*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 511*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, cellidname, NULL, NULL, (void **)&cellid)); 512*fa0c1127SMark Adams 513*fa0c1127SMark Adams // Set particle coordinates 514*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x)); 515*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&v)); 516*fa0c1127SMark Adams PetscCall(DMSwarmSortGetAccess(sw)); 517*fa0c1127SMark Adams PetscCall(DMGetBoundingBox(vdm, vmin, vmax)); 518*fa0c1127SMark Adams PetscCall(DMGetCoordinatesLocalSetUp(xdm)); 519*fa0c1127SMark Adams PetscCall(DMGetCoordinatesLocalSetUp(vdm)); 520*fa0c1127SMark Adams for (PetscInt c = 0; c < xcEnd - xcStart; ++c) { 521*fa0c1127SMark Adams const PetscInt cell = c + xcStart; 522*fa0c1127SMark Adams PetscInt *pidx, Npc; 523*fa0c1127SMark Adams PetscReal centroid[3], volume; 524*fa0c1127SMark Adams 525*fa0c1127SMark Adams PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Npc, &pidx)); 526*fa0c1127SMark Adams PetscCall(DMPlexComputeCellGeometryFVM(xdm, cell, &volume, centroid, NULL)); 527*fa0c1127SMark Adams for (PetscInt s = 0; s < Ns; ++s) { 528*fa0c1127SMark Adams for (PetscInt q = 0; q < Npc / Ns; ++q) { 529*fa0c1127SMark Adams const PetscInt p = pidx[q * Ns + s]; 530*fa0c1127SMark Adams const PetscInt vc = vcStart + q; 531*fa0c1127SMark Adams PetscReal vcentroid[3], vvolume; 532*fa0c1127SMark Adams 533*fa0c1127SMark Adams PetscCall(DMPlexComputeCellGeometryFVM(vdm, vc, &vvolume, vcentroid, NULL)); 534*fa0c1127SMark Adams for (PetscInt d = 0; d < dim; ++d) { 535*fa0c1127SMark Adams x[p * dim + d] = centroid[d]; 536*fa0c1127SMark Adams v[p * dim + d] = vcentroid[d]; 537*fa0c1127SMark Adams } 538*fa0c1127SMark Adams 539*fa0c1127SMark Adams if (debug > 1) { 540*fa0c1127SMark Adams PetscCall(PetscPrintf(PETSC_COMM_SELF, "Particle %4" PetscInt_FMT " ", p)); 541*fa0c1127SMark Adams PetscCall(PetscPrintf(PETSC_COMM_SELF, " x: (")); 542*fa0c1127SMark Adams for (PetscInt d = 0; d < dim; ++d) { 543*fa0c1127SMark Adams if (d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 544*fa0c1127SMark Adams PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", x[p * dim + d])); 545*fa0c1127SMark Adams } 546*fa0c1127SMark Adams PetscCall(PetscPrintf(PETSC_COMM_SELF, ") v:(")); 547*fa0c1127SMark Adams for (PetscInt d = 0; d < dim; ++d) { 548*fa0c1127SMark Adams if (d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 549*fa0c1127SMark Adams PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", v[p * dim + d])); 550*fa0c1127SMark Adams } 551*fa0c1127SMark Adams PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n")); 552*fa0c1127SMark Adams } 553*fa0c1127SMark Adams } 554*fa0c1127SMark Adams } 555*fa0c1127SMark Adams PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Npc, &pidx)); 556*fa0c1127SMark Adams } 557*fa0c1127SMark Adams PetscCall(DMSwarmSortRestoreAccess(sw)); 558*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x)); 559*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&v)); 560*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 561*fa0c1127SMark Adams } 562*fa0c1127SMark Adams 563*fa0c1127SMark Adams /* 564*fa0c1127SMark Adams InitializeWeights - Compute weight for each local particle 565*fa0c1127SMark Adams 566*fa0c1127SMark Adams Input Parameters: 567*fa0c1127SMark Adams + sw - The `DMSwarm` 568*fa0c1127SMark Adams . totalWeight - The sum of all particle weights 569*fa0c1127SMark Adams . func - The PDF for the particle spatial distribution 570*fa0c1127SMark Adams - param - The PDF parameters 571*fa0c1127SMark Adams 572*fa0c1127SMark Adams Notes: 573*fa0c1127SMark Adams The PDF for velocity is assumed to be a Gaussian 574*fa0c1127SMark Adams 575*fa0c1127SMark Adams The particle weights are returned in the `w_q` field of `sw`. 576*fa0c1127SMark Adams */ 577*fa0c1127SMark Adams static PetscErrorCode InitializeWeights(DM sw, PetscReal totalWeight, PetscProbFn *func, const PetscReal param[]) 578*fa0c1127SMark Adams { 579*fa0c1127SMark Adams DM xdm, vdm; 580*fa0c1127SMark Adams DMSwarmCellDM celldm; 581*fa0c1127SMark Adams PetscScalar *weight; 582*fa0c1127SMark Adams PetscQuadrature xquad; 583*fa0c1127SMark Adams const PetscReal *xq, *xwq; 584*fa0c1127SMark Adams const PetscInt order = 5; 585*fa0c1127SMark Adams PetscReal xi0[3]; 586*fa0c1127SMark Adams PetscReal xwtot = 0., pwtot = 0.; 587*fa0c1127SMark Adams PetscInt xNq; 588*fa0c1127SMark Adams PetscInt dim, Ns, xcStart, xcEnd, vcStart, vcEnd, debug = ((DM_Swarm *)sw->data)->printWeights; 589*fa0c1127SMark Adams MPI_Comm comm; 590*fa0c1127SMark Adams PetscMPIInt rank; 591*fa0c1127SMark Adams 592*fa0c1127SMark Adams PetscFunctionBegin; 593*fa0c1127SMark Adams PetscCall(PetscObjectGetComm((PetscObject)sw, &comm)); 594*fa0c1127SMark Adams PetscCallMPI(MPI_Comm_rank(comm, &rank)); 595*fa0c1127SMark Adams PetscCall(DMGetDimension(sw, &dim)); 596*fa0c1127SMark Adams PetscCall(DMSwarmGetCellDM(sw, &xdm)); 597*fa0c1127SMark Adams PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 598*fa0c1127SMark Adams PetscCall(DMPlexGetHeightStratum(xdm, 0, &xcStart, &xcEnd)); 599*fa0c1127SMark Adams PetscCall(DMSwarmGetCellDMByName(sw, "velocity", &celldm)); 600*fa0c1127SMark Adams PetscCall(DMSwarmCellDMGetDM(celldm, &vdm)); 601*fa0c1127SMark Adams PetscCall(DMPlexGetHeightStratum(vdm, 0, &vcStart, &vcEnd)); 602*fa0c1127SMark Adams 603*fa0c1127SMark Adams // Setup Quadrature for spatial and velocity weight calculations 604*fa0c1127SMark Adams PetscCall(PetscDTGaussTensorQuadrature(dim, 1, order, -1.0, 1.0, &xquad)); 605*fa0c1127SMark Adams PetscCall(PetscQuadratureGetData(xquad, NULL, NULL, &xNq, &xq, &xwq)); 606*fa0c1127SMark Adams for (PetscInt d = 0; d < dim; ++d) xi0[d] = -1.0; 607*fa0c1127SMark Adams 608*fa0c1127SMark Adams // Integrate the density function to get the weights of particles in each cell 609*fa0c1127SMark Adams PetscCall(DMGetCoordinatesLocalSetUp(vdm)); 610*fa0c1127SMark Adams PetscCall(DMSwarmSortGetAccess(sw)); 611*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&weight)); 612*fa0c1127SMark Adams for (PetscInt c = xcStart; c < xcEnd; ++c) { 613*fa0c1127SMark Adams PetscReal xv0[3], xJ[9], xinvJ[9], xdetJ, xqr[3], xden, xw = 0.; 614*fa0c1127SMark Adams PetscInt *pidx, Npc; 615*fa0c1127SMark Adams PetscInt xNc; 616*fa0c1127SMark Adams const PetscScalar *xarray; 617*fa0c1127SMark Adams PetscScalar *xcoords = NULL; 618*fa0c1127SMark Adams PetscBool xisDG; 619*fa0c1127SMark Adams 620*fa0c1127SMark Adams PetscCall(DMPlexGetCellCoordinates(xdm, c, &xisDG, &xNc, &xarray, &xcoords)); 621*fa0c1127SMark Adams PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Npc, &pidx)); 622*fa0c1127SMark Adams PetscCheck(Npc == (vcEnd - vcStart) * Ns, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of particles %" PetscInt_FMT " in cell (rank %d) != %" PetscInt_FMT " number of velocity vertices", Npc, rank, (vcEnd - vcStart) * Ns); 623*fa0c1127SMark Adams PetscCall(DMPlexComputeCellGeometryFEM(xdm, c, NULL, xv0, xJ, xinvJ, &xdetJ)); 624*fa0c1127SMark Adams for (PetscInt q = 0; q < xNq; ++q) { 625*fa0c1127SMark Adams // Transform quadrature points from ref space to real space 626*fa0c1127SMark Adams CoordinatesRefToReal(dim, dim, xi0, xv0, xJ, &xq[q * dim], xqr); 627*fa0c1127SMark Adams // Get probability density at quad point 628*fa0c1127SMark Adams // No need to scale xqr since PDF will be periodic 629*fa0c1127SMark Adams PetscCall((*func)(xqr, param, &xden)); 630*fa0c1127SMark Adams xw += xden * (xwq[q] * xdetJ); 631*fa0c1127SMark Adams } 632*fa0c1127SMark Adams xwtot += xw; 633*fa0c1127SMark Adams if (debug) { 634*fa0c1127SMark Adams IS globalOrdering; 635*fa0c1127SMark Adams const PetscInt *ordering; 636*fa0c1127SMark Adams 637*fa0c1127SMark Adams PetscCall(DMPlexGetCellNumbering(xdm, &globalOrdering)); 638*fa0c1127SMark Adams PetscCall(ISGetIndices(globalOrdering, &ordering)); 639*fa0c1127SMark Adams PetscCall(PetscSynchronizedPrintf(comm, "c:%" PetscInt_FMT " [x_a,x_b] = %1.15f,%1.15f -> cell weight = %1.15f\n", ordering[c], (double)PetscRealPart(xcoords[0]), (double)PetscRealPart(xcoords[0 + dim]), (double)xw)); 640*fa0c1127SMark Adams PetscCall(ISRestoreIndices(globalOrdering, &ordering)); 641*fa0c1127SMark Adams } 642*fa0c1127SMark Adams // Set weights to be Gaussian in velocity cells 643*fa0c1127SMark Adams for (PetscInt vc = vcStart; vc < vcEnd; ++vc) { 644*fa0c1127SMark Adams const PetscInt p = pidx[vc * Ns + 0]; 645*fa0c1127SMark Adams PetscReal vw = 0.; 646*fa0c1127SMark Adams PetscInt vNc; 647*fa0c1127SMark Adams const PetscScalar *varray; 648*fa0c1127SMark Adams PetscScalar *vcoords = NULL; 649*fa0c1127SMark Adams PetscBool visDG; 650*fa0c1127SMark Adams 651*fa0c1127SMark Adams PetscCall(DMPlexGetCellCoordinates(vdm, vc, &visDG, &vNc, &varray, &vcoords)); 652*fa0c1127SMark Adams PetscCheck(vNc > 0 && vNc % dim == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Velocity cell %" PetscInt_FMT " has invalid coordinates (vNc=%" PetscInt_FMT ", dim=%" PetscInt_FMT ")", vc, vNc, dim); 653*fa0c1127SMark Adams { 654*fa0c1127SMark Adams PetscReal vmin[3] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL}; 655*fa0c1127SMark Adams PetscReal vmax[3] = {-PETSC_MAX_REAL, -PETSC_MAX_REAL, -PETSC_MAX_REAL}; 656*fa0c1127SMark Adams PetscInt numVert = vNc / dim; 657*fa0c1127SMark Adams for (PetscInt i = 0; i < numVert; ++i) { 658*fa0c1127SMark Adams for (PetscInt d = 0; d < dim; ++d) { 659*fa0c1127SMark Adams PetscReal coord = PetscRealPart(vcoords[i * dim + d]); 660*fa0c1127SMark Adams vmin[d] = PetscMin(vmin[d], coord); 661*fa0c1127SMark Adams vmax[d] = PetscMax(vmax[d], coord); 662*fa0c1127SMark Adams } 663*fa0c1127SMark Adams } 664*fa0c1127SMark Adams vw = 1.0; 665*fa0c1127SMark Adams for (PetscInt d = 0; d < dim; ++d) vw *= 0.5 * (PetscErfReal(vmax[d] / PetscSqrtReal(2.)) - PetscErfReal(vmin[d] / PetscSqrtReal(2.))); 666*fa0c1127SMark Adams PetscCheck(PetscIsNormalReal(vw), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Particle %" PetscInt_FMT " velocity weight is not normal: vw=%g, vmin=(%g,%g,%g), vmax=(%g,%g,%g)", p, vw, vmin[0], vmin[1], vmin[2], vmax[0], vmax[1], vmax[2]); 667*fa0c1127SMark Adams } 668*fa0c1127SMark Adams 669*fa0c1127SMark Adams weight[p] = totalWeight * vw * xw; 670*fa0c1127SMark Adams pwtot += weight[p]; 671*fa0c1127SMark Adams PetscCheck(weight[p] <= 10., PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Particle %" PetscInt_FMT " weight exceeded 10: weight=%g, xw=%g, vw=%g, totalWeight=%g", p, weight[p], xw, vw, totalWeight); 672*fa0c1127SMark Adams PetscCall(DMPlexRestoreCellCoordinates(vdm, vc, &visDG, &vNc, &varray, &vcoords)); 673*fa0c1127SMark Adams if (debug > 1) PetscCall(PetscPrintf(comm, "particle %" PetscInt_FMT ": %g, vw: %g xw: %g\n", p, weight[p], vw, xw)); 674*fa0c1127SMark Adams } 675*fa0c1127SMark Adams PetscCall(DMPlexRestoreCellCoordinates(xdm, c, &xisDG, &xNc, &xarray, &xcoords)); 676*fa0c1127SMark Adams PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Npc, &pidx)); 677*fa0c1127SMark Adams } 678*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight)); 679*fa0c1127SMark Adams PetscCall(DMSwarmSortRestoreAccess(sw)); 680*fa0c1127SMark Adams PetscCall(PetscQuadratureDestroy(&xquad)); 681*fa0c1127SMark Adams 682*fa0c1127SMark Adams if (debug) { 683*fa0c1127SMark Adams PetscReal wtot[2] = {pwtot, xwtot}, gwtot[2]; 684*fa0c1127SMark Adams 685*fa0c1127SMark Adams PetscCall(PetscSynchronizedFlush(comm, NULL)); 686*fa0c1127SMark Adams PetscCallMPI(MPIU_Allreduce(wtot, gwtot, 2, MPIU_REAL, MPIU_SUM, PETSC_COMM_WORLD)); 687*fa0c1127SMark Adams PetscCall(PetscPrintf(comm, "particle weight sum = %1.10f cell weight sum = %1.10f\n", (double)gwtot[0], (double)gwtot[1])); 688*fa0c1127SMark Adams } 689*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 690*fa0c1127SMark Adams } 691*fa0c1127SMark Adams 692*fa0c1127SMark Adams static PetscErrorCode InitializeParticles_PerturbedWeights(DM sw, AppCtx *user) 693*fa0c1127SMark Adams { 694*fa0c1127SMark Adams PetscReal scale[2] = {user->cosine_coefficients[0], user->cosine_coefficients[1]}; 695*fa0c1127SMark Adams PetscInt dim; 696*fa0c1127SMark Adams 697*fa0c1127SMark Adams PetscFunctionBegin; 698*fa0c1127SMark Adams PetscCall(DMGetDimension(sw, &dim)); 699*fa0c1127SMark Adams PetscCall(InitializeParticles_Centroid(sw)); 700*fa0c1127SMark Adams PetscCall(InitializeWeights(sw, user->totalWeight, dim == 1 ? PetscPDFCosine1D : (dim == 2 ? PetscPDFCosine2D : PetscPDFCosine3D), scale)); 701*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 702*fa0c1127SMark Adams } 703*fa0c1127SMark Adams 704*fa0c1127SMark Adams static PetscErrorCode InitializeConstants(DM sw, AppCtx *user) 705*fa0c1127SMark Adams { 706*fa0c1127SMark Adams DM dm; 707*fa0c1127SMark Adams PetscInt *species; 708*fa0c1127SMark Adams PetscReal *weight, totalCharge = 0., totalWeight = 0., gmin[3], gmax[3], global_charge, global_weight; 709*fa0c1127SMark Adams PetscInt Np, dim; 710*fa0c1127SMark Adams 711*fa0c1127SMark Adams PetscFunctionBegin; 712*fa0c1127SMark Adams PetscCall(DMSwarmGetCellDM(sw, &dm)); 713*fa0c1127SMark Adams PetscCall(DMGetDimension(sw, &dim)); 714*fa0c1127SMark Adams PetscCall(DMSwarmGetLocalSize(sw, &Np)); 715*fa0c1127SMark Adams PetscCall(DMGetBoundingBox(dm, gmin, gmax)); 716*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&weight)); 717*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 718*fa0c1127SMark Adams for (PetscInt p = 0; p < Np; ++p) { 719*fa0c1127SMark Adams totalWeight += weight[p]; 720*fa0c1127SMark Adams totalCharge += user->charges[species[p]] * weight[p]; 721*fa0c1127SMark Adams } 722*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight)); 723*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 724*fa0c1127SMark Adams { 725*fa0c1127SMark Adams Parameter *param; 726*fa0c1127SMark Adams PetscReal Area; 727*fa0c1127SMark Adams 728*fa0c1127SMark Adams PetscCall(PetscBagGetData(user->bag, (void **)¶m)); 729*fa0c1127SMark Adams switch (dim) { 730*fa0c1127SMark Adams case 1: 731*fa0c1127SMark Adams Area = (gmax[0] - gmin[0]); 732*fa0c1127SMark Adams break; 733*fa0c1127SMark Adams case 2: 734*fa0c1127SMark Adams Area = (gmax[0] - gmin[0]) * (gmax[1] - gmin[1]); 735*fa0c1127SMark Adams break; 736*fa0c1127SMark Adams case 3: 737*fa0c1127SMark Adams Area = (gmax[0] - gmin[0]) * (gmax[1] - gmin[1]) * (gmax[2] - gmin[2]); 738*fa0c1127SMark Adams break; 739*fa0c1127SMark Adams default: 740*fa0c1127SMark Adams SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Dimension %" PetscInt_FMT " not supported", dim); 741*fa0c1127SMark Adams } 742*fa0c1127SMark Adams PetscCallMPI(MPIU_Allreduce(&totalWeight, &global_weight, 1, MPIU_REAL, MPIU_SUM, PETSC_COMM_WORLD)); 743*fa0c1127SMark Adams PetscCallMPI(MPIU_Allreduce(&totalCharge, &global_charge, 1, MPIU_REAL, MPIU_SUM, PETSC_COMM_WORLD)); 744*fa0c1127SMark Adams PetscCall(PetscPrintf(PETSC_COMM_WORLD, "dim = %" PetscInt_FMT "\ttotalWeight = %f, user->charges[species[0]] = %f\ttotalCharge = %f, Total Area = %f\n", dim, (double)global_weight, (double)user->charges[0], (double)global_charge, (double)Area)); 745*fa0c1127SMark Adams param->sigma = PetscAbsReal(global_charge / (Area)); 746*fa0c1127SMark Adams 747*fa0c1127SMark Adams PetscCall(PetscPrintf(PETSC_COMM_WORLD, "sigma: %g\n", (double)param->sigma)); 748*fa0c1127SMark Adams PetscCall(PetscPrintf(PETSC_COMM_WORLD, "(x0,v0,t0,m0,q0,phi0): (%e, %e, %e, %e, %e, %e) - (P, V) = (%e, %e)\n", (double)param->x0, (double)param->v0, (double)param->t0, (double)param->m0, (double)param->q0, (double)param->phi0, (double)param->poissonNumber, 749*fa0c1127SMark Adams (double)param->vlasovNumber)); 750*fa0c1127SMark Adams } 751*fa0c1127SMark Adams /* Setup Constants */ 752*fa0c1127SMark Adams { 753*fa0c1127SMark Adams PetscDS ds; 754*fa0c1127SMark Adams Parameter *param; 755*fa0c1127SMark Adams PetscCall(PetscBagGetData(user->bag, (void **)¶m)); 756*fa0c1127SMark Adams PetscScalar constants[NUM_CONSTANTS]; 757*fa0c1127SMark Adams constants[SIGMA] = param->sigma; 758*fa0c1127SMark Adams constants[V0] = param->v0; 759*fa0c1127SMark Adams constants[T0] = param->t0; 760*fa0c1127SMark Adams constants[X0] = param->x0; 761*fa0c1127SMark Adams constants[M0] = param->m0; 762*fa0c1127SMark Adams constants[Q0] = param->q0; 763*fa0c1127SMark Adams constants[PHI0] = param->phi0; 764*fa0c1127SMark Adams constants[POISSON] = param->poissonNumber; 765*fa0c1127SMark Adams constants[VLASOV] = param->vlasovNumber; 766*fa0c1127SMark Adams PetscCall(DMGetDS(dm, &ds)); 767*fa0c1127SMark Adams PetscCall(PetscDSSetConstants(ds, NUM_CONSTANTS, constants)); 768*fa0c1127SMark Adams } 769*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 770*fa0c1127SMark Adams } 771*fa0c1127SMark Adams 772*fa0c1127SMark Adams static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *ctx) 773*fa0c1127SMark Adams { 774*fa0c1127SMark Adams PetscBag bag; 775*fa0c1127SMark Adams Parameter *p; 776*fa0c1127SMark Adams 777*fa0c1127SMark Adams PetscFunctionBeginUser; 778*fa0c1127SMark Adams /* setup PETSc parameter bag */ 779*fa0c1127SMark Adams PetscCall(PetscBagGetData(ctx->bag, (void **)&p)); 780*fa0c1127SMark Adams PetscCall(PetscBagSetName(ctx->bag, "par", "Vlasov-Poisson Parameters")); 781*fa0c1127SMark Adams bag = ctx->bag; 782*fa0c1127SMark Adams PetscCall(PetscBagRegisterScalar(bag, &p->v0, 1.0, "v0", "Velocity scale, m/s")); 783*fa0c1127SMark Adams PetscCall(PetscBagRegisterScalar(bag, &p->t0, 1.0, "t0", "Time scale, s")); 784*fa0c1127SMark Adams PetscCall(PetscBagRegisterScalar(bag, &p->x0, 1.0, "x0", "Space scale, m")); 785*fa0c1127SMark Adams PetscCall(PetscBagRegisterScalar(bag, &p->phi0, 1.0, "phi0", "Potential scale, kg*m^2/A*s^3")); 786*fa0c1127SMark Adams PetscCall(PetscBagRegisterScalar(bag, &p->q0, 1.0, "q0", "Charge Scale, A*s")); 787*fa0c1127SMark Adams PetscCall(PetscBagRegisterScalar(bag, &p->m0, 1.0, "m0", "Mass Scale, kg")); 788*fa0c1127SMark Adams PetscCall(PetscBagRegisterScalar(bag, &p->epsi0, 1.0, "epsi0", "Permittivity of Free Space, kg")); 789*fa0c1127SMark Adams PetscCall(PetscBagRegisterScalar(bag, &p->kb, 1.0, "kb", "Boltzmann Constant, m^2 kg/s^2 K^1")); 790*fa0c1127SMark Adams 791*fa0c1127SMark Adams PetscCall(PetscBagRegisterScalar(bag, &p->sigma, 1.0, "sigma", "Charge per unit area, C/m^3")); 792*fa0c1127SMark Adams PetscCall(PetscBagRegisterScalar(bag, &p->poissonNumber, 1.0, "poissonNumber", "Non-Dimensional Poisson Number")); 793*fa0c1127SMark Adams PetscCall(PetscBagRegisterScalar(bag, &p->vlasovNumber, 1.0, "vlasovNumber", "Non-Dimensional Vlasov Number")); 794*fa0c1127SMark Adams PetscCall(PetscBagSetFromOptions(bag)); 795*fa0c1127SMark Adams { 796*fa0c1127SMark Adams PetscViewer viewer; 797*fa0c1127SMark Adams PetscViewerFormat format; 798*fa0c1127SMark Adams PetscBool flg; 799*fa0c1127SMark Adams 800*fa0c1127SMark Adams PetscCall(PetscOptionsCreateViewer(comm, NULL, NULL, "-param_view", &viewer, &format, &flg)); 801*fa0c1127SMark Adams if (flg) { 802*fa0c1127SMark Adams PetscCall(PetscViewerPushFormat(viewer, format)); 803*fa0c1127SMark Adams PetscCall(PetscBagView(bag, viewer)); 804*fa0c1127SMark Adams PetscCall(PetscViewerFlush(viewer)); 805*fa0c1127SMark Adams PetscCall(PetscViewerPopFormat(viewer)); 806*fa0c1127SMark Adams PetscCall(PetscViewerDestroy(&viewer)); 807*fa0c1127SMark Adams } 808*fa0c1127SMark Adams } 809*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 810*fa0c1127SMark Adams } 811*fa0c1127SMark Adams 812*fa0c1127SMark Adams static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 813*fa0c1127SMark Adams { 814*fa0c1127SMark Adams const char *prefix = "x"; 815*fa0c1127SMark Adams 816*fa0c1127SMark Adams PetscFunctionBeginUser; 817*fa0c1127SMark Adams PetscCall(DMCreate(comm, dm)); 818*fa0c1127SMark Adams PetscCall(DMSetType(*dm, DMPLEX)); 819*fa0c1127SMark Adams PetscCall(DMPlexSetOptionsPrefix(*dm, prefix)); 820*fa0c1127SMark Adams PetscCall(DMSetFromOptions(*dm)); 821*fa0c1127SMark Adams PetscCall(PetscObjectSetName((PetscObject)*dm, "space")); 822*fa0c1127SMark Adams PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 823*fa0c1127SMark Adams 824*fa0c1127SMark Adams // Cache the mesh geometry 825*fa0c1127SMark Adams DMField coordField; 826*fa0c1127SMark Adams IS cellIS; 827*fa0c1127SMark Adams PetscQuadrature quad; 828*fa0c1127SMark Adams PetscReal *wt, *pt; 829*fa0c1127SMark Adams PetscInt cdim, cStart, cEnd; 830*fa0c1127SMark Adams 831*fa0c1127SMark Adams PetscCall(DMGetCoordinateField(*dm, &coordField)); 832*fa0c1127SMark Adams PetscCheck(coordField, comm, PETSC_ERR_USER, "DM must have a coordinate field"); 833*fa0c1127SMark Adams PetscCall(DMGetCoordinateDim(*dm, &cdim)); 834*fa0c1127SMark Adams PetscCall(DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd)); 835*fa0c1127SMark Adams PetscCall(ISCreateStride(PETSC_COMM_SELF, cEnd - cStart, cStart, 1, &cellIS)); 836*fa0c1127SMark Adams PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &quad)); 837*fa0c1127SMark Adams PetscCall(PetscMalloc1(1, &wt)); 838*fa0c1127SMark Adams PetscCall(PetscMalloc1(cdim, &pt)); 839*fa0c1127SMark Adams wt[0] = 1.; 840*fa0c1127SMark Adams for (PetscInt d = 0; d < cdim; ++d) pt[d] = -1.; 841*fa0c1127SMark Adams PetscCall(PetscQuadratureSetData(quad, cdim, 1, 1, pt, wt)); 842*fa0c1127SMark Adams PetscCall(DMFieldCreateFEGeom(coordField, cellIS, quad, PETSC_FEGEOM_BASIC, &user->fegeom)); 843*fa0c1127SMark Adams PetscCall(PetscQuadratureDestroy(&quad)); 844*fa0c1127SMark Adams PetscCall(ISDestroy(&cellIS)); 845*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 846*fa0c1127SMark Adams } 847*fa0c1127SMark Adams 848*fa0c1127SMark Adams static void ion_f0(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 849*fa0c1127SMark Adams { 850*fa0c1127SMark Adams f0[0] = -constants[SIGMA]; 851*fa0c1127SMark Adams } 852*fa0c1127SMark Adams 853*fa0c1127SMark Adams static void laplacian_f1(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 854*fa0c1127SMark Adams { 855*fa0c1127SMark Adams PetscInt d; 856*fa0c1127SMark Adams for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 857*fa0c1127SMark Adams } 858*fa0c1127SMark Adams 859*fa0c1127SMark Adams static void laplacian_g3(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 860*fa0c1127SMark Adams { 861*fa0c1127SMark Adams PetscInt d; 862*fa0c1127SMark Adams for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0; 863*fa0c1127SMark Adams } 864*fa0c1127SMark Adams 865*fa0c1127SMark Adams static PetscErrorCode CreateFEM(DM dm, AppCtx *user) 866*fa0c1127SMark Adams { 867*fa0c1127SMark Adams PetscFE fephi; 868*fa0c1127SMark Adams PetscDS ds; 869*fa0c1127SMark Adams PetscBool simplex; 870*fa0c1127SMark Adams PetscInt dim; 871*fa0c1127SMark Adams MatNullSpace nullsp; 872*fa0c1127SMark Adams 873*fa0c1127SMark Adams PetscFunctionBeginUser; 874*fa0c1127SMark Adams PetscCall(DMGetDimension(dm, &dim)); 875*fa0c1127SMark Adams PetscCall(DMPlexIsSimplex(dm, &simplex)); 876*fa0c1127SMark Adams 877*fa0c1127SMark Adams PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, PETSC_DETERMINE, &fephi)); 878*fa0c1127SMark Adams PetscCall(PetscObjectSetName((PetscObject)fephi, "potential")); 879*fa0c1127SMark Adams PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fephi)); 880*fa0c1127SMark Adams PetscCall(DMCreateDS(dm)); 881*fa0c1127SMark Adams PetscCall(DMGetDS(dm, &ds)); 882*fa0c1127SMark Adams PetscCall(PetscDSSetResidual(ds, 0, ion_f0, laplacian_f1)); 883*fa0c1127SMark Adams PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, laplacian_g3)); 884*fa0c1127SMark Adams PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_TRUE, 0, NULL, &nullsp)); 885*fa0c1127SMark Adams PetscCall(PetscObjectCompose((PetscObject)fephi, "nullspace", (PetscObject)nullsp)); 886*fa0c1127SMark Adams PetscCall(MatNullSpaceDestroy(&nullsp)); 887*fa0c1127SMark Adams PetscCall(PetscFEDestroy(&fephi)); 888*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 889*fa0c1127SMark Adams } 890*fa0c1127SMark Adams 891*fa0c1127SMark Adams static PetscErrorCode CreatePoisson(DM dm, AppCtx *user) 892*fa0c1127SMark Adams { 893*fa0c1127SMark Adams SNES snes; 894*fa0c1127SMark Adams Mat J; 895*fa0c1127SMark Adams MatNullSpace nullSpace; 896*fa0c1127SMark Adams 897*fa0c1127SMark Adams PetscFunctionBeginUser; 898*fa0c1127SMark Adams PetscCall(CreateFEM(dm, user)); 899*fa0c1127SMark Adams PetscCall(SNESCreate(PetscObjectComm((PetscObject)dm), &snes)); 900*fa0c1127SMark Adams PetscCall(SNESSetOptionsPrefix(snes, "em_")); 901*fa0c1127SMark Adams PetscCall(SNESSetDM(snes, dm)); 902*fa0c1127SMark Adams PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, user)); 903*fa0c1127SMark Adams PetscCall(SNESSetFromOptions(snes)); 904*fa0c1127SMark Adams 905*fa0c1127SMark Adams PetscCall(DMCreateMatrix(dm, &J)); 906*fa0c1127SMark Adams PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_TRUE, 0, NULL, &nullSpace)); 907*fa0c1127SMark Adams PetscCall(MatSetNullSpace(J, nullSpace)); 908*fa0c1127SMark Adams PetscCall(MatNullSpaceDestroy(&nullSpace)); 909*fa0c1127SMark Adams PetscCall(SNESSetJacobian(snes, J, J, NULL, NULL)); 910*fa0c1127SMark Adams PetscCall(MatDestroy(&J)); 911*fa0c1127SMark Adams 912*fa0c1127SMark Adams user->dmPot = dm; 913*fa0c1127SMark Adams PetscCall(PetscObjectReference((PetscObject)user->dmPot)); 914*fa0c1127SMark Adams 915*fa0c1127SMark Adams PetscCall(DMCreateMassMatrix(user->dmPot, user->dmPot, &user->M)); 916*fa0c1127SMark Adams PetscCall(DMPlexCreateClosureIndex(dm, NULL)); 917*fa0c1127SMark Adams user->snes = snes; 918*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 919*fa0c1127SMark Adams } 920*fa0c1127SMark Adams 921*fa0c1127SMark Adams static PetscErrorCode CreateSwarm(DM dm, AppCtx *user, DM *sw) 922*fa0c1127SMark Adams { 923*fa0c1127SMark Adams DMSwarmCellDM celldm; 924*fa0c1127SMark Adams PetscInt dim; 925*fa0c1127SMark Adams 926*fa0c1127SMark Adams PetscFunctionBeginUser; 927*fa0c1127SMark Adams PetscCall(DMGetDimension(dm, &dim)); 928*fa0c1127SMark Adams PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw)); 929*fa0c1127SMark Adams PetscCall(DMSetType(*sw, DMSWARM)); 930*fa0c1127SMark Adams PetscCall(DMSetDimension(*sw, dim)); 931*fa0c1127SMark Adams PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC)); 932*fa0c1127SMark Adams PetscCall(DMSetApplicationContext(*sw, user)); 933*fa0c1127SMark Adams 934*fa0c1127SMark Adams PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR)); 935*fa0c1127SMark Adams PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "velocity", dim, PETSC_REAL)); 936*fa0c1127SMark Adams PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "species", 1, PETSC_INT)); 937*fa0c1127SMark Adams PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "E_field", dim, PETSC_REAL)); 938*fa0c1127SMark Adams 939*fa0c1127SMark Adams const char *fieldnames[2] = {DMSwarmPICField_coor, "velocity"}; 940*fa0c1127SMark Adams 941*fa0c1127SMark Adams PetscCall(DMSwarmCellDMCreate(dm, 2, fieldnames, 1, fieldnames, &celldm)); 942*fa0c1127SMark Adams PetscCall(DMSwarmAddCellDM(*sw, celldm)); 943*fa0c1127SMark Adams PetscCall(DMSwarmCellDMDestroy(&celldm)); 944*fa0c1127SMark Adams 945*fa0c1127SMark Adams const char *vfieldnames[2] = {"w_q"}; 946*fa0c1127SMark Adams DM vdm; 947*fa0c1127SMark Adams 948*fa0c1127SMark Adams PetscCall(CreateVelocityDM(*sw, &vdm)); 949*fa0c1127SMark Adams PetscCall(DMSwarmCellDMCreate(vdm, 1, vfieldnames, 1, &fieldnames[1], &celldm)); 950*fa0c1127SMark Adams PetscCall(DMSwarmAddCellDM(*sw, celldm)); 951*fa0c1127SMark Adams PetscCall(DMSwarmCellDMDestroy(&celldm)); 952*fa0c1127SMark Adams PetscCall(DMDestroy(&vdm)); 953*fa0c1127SMark Adams 954*fa0c1127SMark Adams DM mdm; 955*fa0c1127SMark Adams 956*fa0c1127SMark Adams PetscCall(DMClone(dm, &mdm)); 957*fa0c1127SMark Adams PetscCall(PetscObjectSetName((PetscObject)mdm, "moments")); 958*fa0c1127SMark Adams PetscCall(DMCopyDisc(dm, mdm)); 959*fa0c1127SMark Adams PetscCall(DMSwarmCellDMCreate(mdm, 1, vfieldnames, 1, fieldnames, &celldm)); 960*fa0c1127SMark Adams PetscCall(DMDestroy(&mdm)); 961*fa0c1127SMark Adams PetscCall(DMSwarmAddCellDM(*sw, celldm)); 962*fa0c1127SMark Adams PetscCall(DMSwarmCellDMDestroy(&celldm)); 963*fa0c1127SMark Adams 964*fa0c1127SMark Adams PetscCall(DMSetFromOptions(*sw)); 965*fa0c1127SMark Adams PetscCall(DMSetUp(*sw)); 966*fa0c1127SMark Adams 967*fa0c1127SMark Adams PetscCall(DMSwarmSetCellDMActive(*sw, "space")); 968*fa0c1127SMark Adams user->swarm = *sw; 969*fa0c1127SMark Adams // TODO: This is redundant init since it is done in InitializeSolveAndSwarm, however DMSetUp() requires the local size be set 970*fa0c1127SMark Adams PetscCall(InitializeParticles_PerturbedWeights(*sw, user)); 971*fa0c1127SMark Adams PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles")); 972*fa0c1127SMark Adams PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view")); 973*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 974*fa0c1127SMark Adams } 975*fa0c1127SMark Adams 976*fa0c1127SMark Adams static PetscErrorCode ComputeFieldAtParticles_Primal(SNES snes, DM sw, Mat M_p, PetscReal E[]) 977*fa0c1127SMark Adams { 978*fa0c1127SMark Adams DM dm; 979*fa0c1127SMark Adams AppCtx *user; 980*fa0c1127SMark Adams PetscDS ds; 981*fa0c1127SMark Adams PetscFE fe; 982*fa0c1127SMark Adams KSP ksp; 983*fa0c1127SMark Adams Vec rhoRhs; // Weak charge density, \int phi_i rho 984*fa0c1127SMark Adams Vec rho; // Charge density, M^{-1} rhoRhs 985*fa0c1127SMark Adams Vec phi, locPhi; // Potential 986*fa0c1127SMark Adams Vec f; // Particle weights 987*fa0c1127SMark Adams PetscReal *coords; 988*fa0c1127SMark Adams PetscInt dim, cStart, cEnd, Np; 989*fa0c1127SMark Adams 990*fa0c1127SMark Adams PetscFunctionBegin; 991*fa0c1127SMark Adams PetscCall(DMGetApplicationContext(sw, (void *)&user)); 992*fa0c1127SMark Adams PetscCall(PetscLogEventBegin(user->ESolveEvent, snes, sw, 0, 0)); 993*fa0c1127SMark Adams PetscCall(DMGetDimension(sw, &dim)); 994*fa0c1127SMark Adams PetscCall(DMSwarmGetLocalSize(sw, &Np)); 995*fa0c1127SMark Adams 996*fa0c1127SMark Adams PetscCall(SNESGetDM(snes, &dm)); 997*fa0c1127SMark Adams PetscCall(DMGetGlobalVector(dm, &rhoRhs)); 998*fa0c1127SMark Adams PetscCall(PetscObjectSetName((PetscObject)rhoRhs, "Weak charge density")); 999*fa0c1127SMark Adams PetscCall(DMGetNamedGlobalVector(user->dmPot, "rho", &rho)); 1000*fa0c1127SMark Adams PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f)); 1001*fa0c1127SMark Adams PetscCall(PetscObjectSetName((PetscObject)f, "particle weight")); 1002*fa0c1127SMark Adams 1003*fa0c1127SMark Adams PetscCall(MatViewFromOptions(M_p, NULL, "-mp_view")); 1004*fa0c1127SMark Adams PetscCall(MatViewFromOptions(user->M, NULL, "-m_view")); 1005*fa0c1127SMark Adams PetscCall(VecViewFromOptions(f, NULL, "-weights_view")); 1006*fa0c1127SMark Adams 1007*fa0c1127SMark Adams PetscCall(MatMultTranspose(M_p, f, rhoRhs)); 1008*fa0c1127SMark Adams PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f)); 1009*fa0c1127SMark Adams 1010*fa0c1127SMark Adams // Low-pass filter rhoRhs 1011*fa0c1127SMark Adams PetscInt window = 0; 1012*fa0c1127SMark Adams PetscCall(PetscOptionsGetInt(NULL, NULL, "-rho_average", &window, NULL)); 1013*fa0c1127SMark Adams if (window) { 1014*fa0c1127SMark Adams PetscScalar *a; 1015*fa0c1127SMark Adams PetscInt n; 1016*fa0c1127SMark Adams PetscReal width = 2. * window + 1.; 1017*fa0c1127SMark Adams 1018*fa0c1127SMark Adams // This will only work for P_1 1019*fa0c1127SMark Adams // This works because integration against a test function is linear 1020*fa0c1127SMark Adams // Do a real integral against weight function for higher order 1021*fa0c1127SMark Adams PetscCall(VecGetLocalSize(rhoRhs, &n)); 1022*fa0c1127SMark Adams PetscCall(VecGetArrayWrite(rhoRhs, &a)); 1023*fa0c1127SMark Adams for (PetscInt i = 0; i < n; ++i) { 1024*fa0c1127SMark Adams PetscScalar avg = a[i]; 1025*fa0c1127SMark Adams for (PetscInt j = 1; j <= window; ++j) avg += a[(i - j + n) % n] + a[(i + j) % n]; 1026*fa0c1127SMark Adams a[i] = avg / width; 1027*fa0c1127SMark Adams //a[i] = (a[(i - 1 + n) % n] + a[i] + a[(i + 1) % n]) / 3.; 1028*fa0c1127SMark Adams } 1029*fa0c1127SMark Adams PetscCall(VecRestoreArrayWrite(rhoRhs, &a)); 1030*fa0c1127SMark Adams } 1031*fa0c1127SMark Adams 1032*fa0c1127SMark Adams PetscCall(KSPCreate(PetscObjectComm((PetscObject)dm), &ksp)); 1033*fa0c1127SMark Adams PetscCall(KSPSetOptionsPrefix(ksp, "em_proj_")); 1034*fa0c1127SMark Adams PetscCall(KSPSetOperators(ksp, user->M, user->M)); 1035*fa0c1127SMark Adams PetscCall(KSPSetFromOptions(ksp)); 1036*fa0c1127SMark Adams PetscCall(KSPSolve(ksp, rhoRhs, rho)); 1037*fa0c1127SMark Adams 1038*fa0c1127SMark Adams PetscCall(VecScale(rhoRhs, -1.0)); 1039*fa0c1127SMark Adams 1040*fa0c1127SMark Adams PetscCall(VecViewFromOptions(rhoRhs, NULL, "-rho_view")); 1041*fa0c1127SMark Adams PetscCall(DMRestoreNamedGlobalVector(user->dmPot, "rho", &rho)); 1042*fa0c1127SMark Adams PetscCall(KSPDestroy(&ksp)); 1043*fa0c1127SMark Adams 1044*fa0c1127SMark Adams PetscCall(DMGetNamedGlobalVector(user->dmPot, "phi", &phi)); 1045*fa0c1127SMark Adams PetscCall(VecSet(phi, 0.0)); 1046*fa0c1127SMark Adams PetscCall(SNESSolve(snes, rhoRhs, phi)); 1047*fa0c1127SMark Adams PetscCall(DMRestoreGlobalVector(dm, &rhoRhs)); 1048*fa0c1127SMark Adams PetscCall(VecViewFromOptions(phi, NULL, "-phi_view")); 1049*fa0c1127SMark Adams 1050*fa0c1127SMark Adams PetscCall(DMGetLocalVector(dm, &locPhi)); 1051*fa0c1127SMark Adams PetscCall(DMGlobalToLocalBegin(dm, phi, INSERT_VALUES, locPhi)); 1052*fa0c1127SMark Adams PetscCall(DMGlobalToLocalEnd(dm, phi, INSERT_VALUES, locPhi)); 1053*fa0c1127SMark Adams PetscCall(DMRestoreNamedGlobalVector(user->dmPot, "phi", &phi)); 1054*fa0c1127SMark Adams PetscCall(PetscLogEventEnd(user->ESolveEvent, snes, sw, 0, 0)); 1055*fa0c1127SMark Adams 1056*fa0c1127SMark Adams PetscCall(DMGetDS(dm, &ds)); 1057*fa0c1127SMark Adams PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&fe)); 1058*fa0c1127SMark Adams PetscCall(DMSwarmSortGetAccess(sw)); 1059*fa0c1127SMark Adams PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 1060*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 1061*fa0c1127SMark Adams 1062*fa0c1127SMark Adams PetscCall(PetscLogEventBegin(user->ETabEvent, snes, sw, 0, 0)); 1063*fa0c1127SMark Adams PetscTabulation tab; 1064*fa0c1127SMark Adams PetscReal *pcoord, *refcoord; 1065*fa0c1127SMark Adams PetscFEGeom *chunkgeom = NULL; 1066*fa0c1127SMark Adams PetscInt maxNcp = 0; 1067*fa0c1127SMark Adams 1068*fa0c1127SMark Adams for (PetscInt c = cStart; c < cEnd; ++c) { 1069*fa0c1127SMark Adams PetscInt Ncp; 1070*fa0c1127SMark Adams 1071*fa0c1127SMark Adams PetscCall(DMSwarmSortGetNumberOfPointsPerCell(sw, c, &Ncp)); 1072*fa0c1127SMark Adams maxNcp = PetscMax(maxNcp, Ncp); 1073*fa0c1127SMark Adams } 1074*fa0c1127SMark Adams PetscCall(DMGetWorkArray(dm, maxNcp * dim, MPIU_REAL, &refcoord)); 1075*fa0c1127SMark Adams PetscCall(PetscArrayzero(refcoord, maxNcp * dim)); 1076*fa0c1127SMark Adams PetscCall(DMGetWorkArray(dm, maxNcp * dim, MPIU_REAL, &pcoord)); 1077*fa0c1127SMark Adams PetscCall(PetscFECreateTabulation(fe, 1, maxNcp, refcoord, 1, &tab)); 1078*fa0c1127SMark Adams for (PetscInt c = cStart; c < cEnd; ++c) { 1079*fa0c1127SMark Adams PetscScalar *clPhi = NULL; 1080*fa0c1127SMark Adams PetscInt *points; 1081*fa0c1127SMark Adams PetscInt Ncp; 1082*fa0c1127SMark Adams 1083*fa0c1127SMark Adams PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Ncp, &points)); 1084*fa0c1127SMark Adams for (PetscInt cp = 0; cp < Ncp; ++cp) 1085*fa0c1127SMark Adams for (PetscInt d = 0; d < dim; ++d) pcoord[cp * dim + d] = coords[points[cp] * dim + d]; 1086*fa0c1127SMark Adams { 1087*fa0c1127SMark Adams PetscCall(PetscFEGeomGetChunk(user->fegeom, c - cStart, c - cStart + 1, &chunkgeom)); 1088*fa0c1127SMark Adams for (PetscInt i = 0; i < Ncp; ++i) { 1089*fa0c1127SMark Adams const PetscReal x0[3] = {-1., -1., -1.}; 1090*fa0c1127SMark Adams CoordinatesRealToRef(dim, dim, x0, chunkgeom->v, chunkgeom->invJ, &pcoord[dim * i], &refcoord[dim * i]); 1091*fa0c1127SMark Adams } 1092*fa0c1127SMark Adams } 1093*fa0c1127SMark Adams PetscCall(PetscFEComputeTabulation(fe, Ncp, refcoord, 1, tab)); 1094*fa0c1127SMark Adams PetscCall(DMPlexVecGetClosure(dm, NULL, locPhi, c, NULL, &clPhi)); 1095*fa0c1127SMark Adams for (PetscInt cp = 0; cp < Ncp; ++cp) { 1096*fa0c1127SMark Adams const PetscReal *basisDer = tab->T[1]; 1097*fa0c1127SMark Adams const PetscInt p = points[cp]; 1098*fa0c1127SMark Adams 1099*fa0c1127SMark Adams for (PetscInt d = 0; d < dim; ++d) E[p * dim + d] = 0.; 1100*fa0c1127SMark Adams PetscCall(PetscFEFreeInterpolateGradient_Static(fe, basisDer, clPhi, dim, chunkgeom->invJ, NULL, cp, &E[p * dim])); 1101*fa0c1127SMark Adams for (PetscInt d = 0; d < dim; ++d) E[p * dim + d] *= -1.0; 1102*fa0c1127SMark Adams } 1103*fa0c1127SMark Adams PetscCall(DMPlexVecRestoreClosure(dm, NULL, locPhi, c, NULL, &clPhi)); 1104*fa0c1127SMark Adams PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Ncp, &points)); 1105*fa0c1127SMark Adams } 1106*fa0c1127SMark Adams PetscCall(DMRestoreWorkArray(dm, maxNcp * dim, MPIU_REAL, &pcoord)); 1107*fa0c1127SMark Adams PetscCall(DMRestoreWorkArray(dm, maxNcp * dim, MPIU_REAL, &refcoord)); 1108*fa0c1127SMark Adams PetscCall(PetscTabulationDestroy(&tab)); 1109*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 1110*fa0c1127SMark Adams PetscCall(DMSwarmSortRestoreAccess(sw)); 1111*fa0c1127SMark Adams PetscCall(DMRestoreLocalVector(dm, &locPhi)); 1112*fa0c1127SMark Adams PetscCall(PetscFEGeomRestoreChunk(user->fegeom, 0, 1, &chunkgeom)); 1113*fa0c1127SMark Adams PetscCall(PetscLogEventEnd(user->ETabEvent, snes, sw, 0, 0)); 1114*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 1115*fa0c1127SMark Adams } 1116*fa0c1127SMark Adams 1117*fa0c1127SMark Adams static PetscErrorCode ComputeFieldAtParticles(SNES snes, DM sw) 1118*fa0c1127SMark Adams { 1119*fa0c1127SMark Adams AppCtx *user; 1120*fa0c1127SMark Adams Mat M_p; 1121*fa0c1127SMark Adams PetscReal *E; 1122*fa0c1127SMark Adams PetscInt dim, Np; 1123*fa0c1127SMark Adams 1124*fa0c1127SMark Adams PetscFunctionBegin; 1125*fa0c1127SMark Adams PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1126*fa0c1127SMark Adams PetscValidHeaderSpecific(sw, DM_CLASSID, 2); 1127*fa0c1127SMark Adams PetscCall(DMGetDimension(sw, &dim)); 1128*fa0c1127SMark Adams PetscCall(DMSwarmGetLocalSize(sw, &Np)); 1129*fa0c1127SMark Adams PetscCall(DMGetApplicationContext(sw, &user)); 1130*fa0c1127SMark Adams 1131*fa0c1127SMark Adams PetscCall(DMSwarmSetCellDMActive(sw, "moments")); 1132*fa0c1127SMark Adams // TODO: Could share sort context with space cellDM 1133*fa0c1127SMark Adams PetscCall(DMSwarmMigrate(sw, PETSC_FALSE)); 1134*fa0c1127SMark Adams PetscCall(DMCreateMassMatrix(sw, user->dmPot, &M_p)); 1135*fa0c1127SMark Adams PetscCall(DMSwarmSetCellDMActive(sw, "space")); 1136*fa0c1127SMark Adams 1137*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, "E_field", NULL, NULL, (void **)&E)); 1138*fa0c1127SMark Adams PetscCall(PetscArrayzero(E, Np * dim)); 1139*fa0c1127SMark Adams user->validE = PETSC_TRUE; 1140*fa0c1127SMark Adams 1141*fa0c1127SMark Adams PetscCall(ComputeFieldAtParticles_Primal(snes, sw, M_p, E)); 1142*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, "E_field", NULL, NULL, (void **)&E)); 1143*fa0c1127SMark Adams PetscCall(MatDestroy(&M_p)); 1144*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 1145*fa0c1127SMark Adams } 1146*fa0c1127SMark Adams 1147*fa0c1127SMark Adams static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec G, void *ctx) 1148*fa0c1127SMark Adams { 1149*fa0c1127SMark Adams DM sw; 1150*fa0c1127SMark Adams SNES snes = ((AppCtx *)ctx)->snes; 1151*fa0c1127SMark Adams const PetscScalar *u; 1152*fa0c1127SMark Adams PetscScalar *g; 1153*fa0c1127SMark Adams PetscReal *E, m_p = 1., q_p = -1.; 1154*fa0c1127SMark Adams PetscInt dim, d, Np, p; 1155*fa0c1127SMark Adams 1156*fa0c1127SMark Adams PetscFunctionBeginUser; 1157*fa0c1127SMark Adams PetscCall(TSGetDM(ts, &sw)); 1158*fa0c1127SMark Adams PetscCall(ComputeFieldAtParticles(snes, sw)); 1159*fa0c1127SMark Adams 1160*fa0c1127SMark Adams PetscCall(DMGetDimension(sw, &dim)); 1161*fa0c1127SMark Adams PetscCall(DMSwarmGetLocalSize(sw, &Np)); 1162*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, "E_field", NULL, NULL, (void **)&E)); 1163*fa0c1127SMark Adams PetscCall(VecGetArrayRead(U, &u)); 1164*fa0c1127SMark Adams PetscCall(VecGetArray(G, &g)); 1165*fa0c1127SMark Adams Np /= 2 * dim; 1166*fa0c1127SMark Adams for (p = 0; p < Np; ++p) { 1167*fa0c1127SMark Adams for (d = 0; d < dim; ++d) { 1168*fa0c1127SMark Adams g[(p * 2 + 0) * dim + d] = u[(p * 2 + 1) * dim + d]; 1169*fa0c1127SMark Adams g[(p * 2 + 1) * dim + d] = q_p * E[p * dim + d] / m_p; 1170*fa0c1127SMark Adams } 1171*fa0c1127SMark Adams } 1172*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, "E_field", NULL, NULL, (void **)&E)); 1173*fa0c1127SMark Adams PetscCall(VecRestoreArrayRead(U, &u)); 1174*fa0c1127SMark Adams PetscCall(VecRestoreArray(G, &g)); 1175*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 1176*fa0c1127SMark Adams } 1177*fa0c1127SMark Adams 1178*fa0c1127SMark Adams /* J_{ij} = dF_i/dx_j 1179*fa0c1127SMark Adams J_p = ( 0 1) 1180*fa0c1127SMark Adams (-w^2 0) 1181*fa0c1127SMark Adams TODO Now there is another term with w^2 from the electric field. I think we will need to invert the operator. 1182*fa0c1127SMark Adams Perhaps we can approximate the Jacobian using only the cellwise P-P gradient from Coulomb 1183*fa0c1127SMark Adams */ 1184*fa0c1127SMark Adams static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat J, Mat P, void *ctx) 1185*fa0c1127SMark Adams { 1186*fa0c1127SMark Adams DM sw; 1187*fa0c1127SMark Adams const PetscReal *coords, *vel; 1188*fa0c1127SMark Adams PetscInt dim, d, Np, p, rStart; 1189*fa0c1127SMark Adams 1190*fa0c1127SMark Adams PetscFunctionBeginUser; 1191*fa0c1127SMark Adams PetscCall(TSGetDM(ts, &sw)); 1192*fa0c1127SMark Adams PetscCall(DMGetDimension(sw, &dim)); 1193*fa0c1127SMark Adams PetscCall(DMSwarmGetLocalSize(sw, &Np)); 1194*fa0c1127SMark Adams PetscCall(MatGetOwnershipRange(J, &rStart, NULL)); 1195*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 1196*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&vel)); 1197*fa0c1127SMark Adams Np /= 2 * dim; 1198*fa0c1127SMark Adams for (p = 0; p < Np; ++p) { 1199*fa0c1127SMark Adams // TODO This is not right because dv/dx has the electric field in it 1200*fa0c1127SMark Adams PetscScalar vals[4] = {0., 1., -1, 0.}; 1201*fa0c1127SMark Adams 1202*fa0c1127SMark Adams for (d = 0; d < dim; ++d) { 1203*fa0c1127SMark Adams const PetscInt rows[2] = {(p * 2 + 0) * dim + d + rStart, (p * 2 + 1) * dim + d + rStart}; 1204*fa0c1127SMark Adams PetscCall(MatSetValues(J, 2, rows, 2, rows, vals, INSERT_VALUES)); 1205*fa0c1127SMark Adams } 1206*fa0c1127SMark Adams } 1207*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 1208*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&vel)); 1209*fa0c1127SMark Adams PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 1210*fa0c1127SMark Adams PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 1211*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 1212*fa0c1127SMark Adams } 1213*fa0c1127SMark Adams 1214*fa0c1127SMark Adams static PetscErrorCode RHSFunctionX(TS ts, PetscReal t, Vec V, Vec Xres, void *ctx) 1215*fa0c1127SMark Adams { 1216*fa0c1127SMark Adams AppCtx *user = (AppCtx *)ctx; 1217*fa0c1127SMark Adams DM sw; 1218*fa0c1127SMark Adams const PetscScalar *v; 1219*fa0c1127SMark Adams PetscScalar *xres; 1220*fa0c1127SMark Adams PetscInt Np, p, d, dim; 1221*fa0c1127SMark Adams 1222*fa0c1127SMark Adams PetscFunctionBeginUser; 1223*fa0c1127SMark Adams PetscCall(PetscLogEventBegin(user->RhsXEvent, ts, 0, 0, 0)); 1224*fa0c1127SMark Adams PetscCall(TSGetDM(ts, &sw)); 1225*fa0c1127SMark Adams PetscCall(DMGetDimension(sw, &dim)); 1226*fa0c1127SMark Adams PetscCall(VecGetLocalSize(Xres, &Np)); 1227*fa0c1127SMark Adams PetscCall(VecGetArrayRead(V, &v)); 1228*fa0c1127SMark Adams PetscCall(VecGetArray(Xres, &xres)); 1229*fa0c1127SMark Adams Np /= dim; 1230*fa0c1127SMark Adams for (p = 0; p < Np; ++p) { 1231*fa0c1127SMark Adams for (d = 0; d < dim; ++d) xres[p * dim + d] = v[p * dim + d]; 1232*fa0c1127SMark Adams } 1233*fa0c1127SMark Adams PetscCall(VecRestoreArrayRead(V, &v)); 1234*fa0c1127SMark Adams PetscCall(VecRestoreArray(Xres, &xres)); 1235*fa0c1127SMark Adams PetscCall(PetscLogEventEnd(user->RhsXEvent, ts, 0, 0, 0)); 1236*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 1237*fa0c1127SMark Adams } 1238*fa0c1127SMark Adams 1239*fa0c1127SMark Adams static PetscErrorCode RHSFunctionV(TS ts, PetscReal t, Vec X, Vec Vres, void *ctx) 1240*fa0c1127SMark Adams { 1241*fa0c1127SMark Adams DM sw; 1242*fa0c1127SMark Adams AppCtx *user = (AppCtx *)ctx; 1243*fa0c1127SMark Adams SNES snes = ((AppCtx *)ctx)->snes; 1244*fa0c1127SMark Adams const PetscScalar *x; 1245*fa0c1127SMark Adams PetscScalar *vres; 1246*fa0c1127SMark Adams PetscReal *E, m_p, q_p; 1247*fa0c1127SMark Adams PetscInt Np, p, dim, d; 1248*fa0c1127SMark Adams Parameter *param; 1249*fa0c1127SMark Adams 1250*fa0c1127SMark Adams PetscFunctionBeginUser; 1251*fa0c1127SMark Adams PetscCall(PetscLogEventBegin(user->RhsVEvent, ts, 0, 0, 0)); 1252*fa0c1127SMark Adams PetscCall(TSGetDM(ts, &sw)); 1253*fa0c1127SMark Adams PetscCall(ComputeFieldAtParticles(snes, sw)); 1254*fa0c1127SMark Adams 1255*fa0c1127SMark Adams PetscCall(DMGetDimension(sw, &dim)); 1256*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, "E_field", NULL, NULL, (void **)&E)); 1257*fa0c1127SMark Adams PetscCall(PetscBagGetData(user->bag, (void **)¶m)); 1258*fa0c1127SMark Adams m_p = user->masses[0] * param->m0; 1259*fa0c1127SMark Adams q_p = user->charges[0] * param->q0; 1260*fa0c1127SMark Adams PetscCall(VecGetLocalSize(Vres, &Np)); 1261*fa0c1127SMark Adams PetscCall(VecGetArrayRead(X, &x)); 1262*fa0c1127SMark Adams PetscCall(VecGetArray(Vres, &vres)); 1263*fa0c1127SMark Adams Np /= dim; 1264*fa0c1127SMark Adams for (p = 0; p < Np; ++p) { 1265*fa0c1127SMark Adams for (d = 0; d < dim; ++d) vres[p * dim + d] = q_p * E[p * dim + d] / m_p; 1266*fa0c1127SMark Adams } 1267*fa0c1127SMark Adams PetscCall(VecRestoreArrayRead(X, &x)); 1268*fa0c1127SMark Adams /* 1269*fa0c1127SMark Adams Synchronized, ordered output for parallel/sequential test cases. 1270*fa0c1127SMark Adams In the 1D (on the 2D mesh) case, every y component should be zero. 1271*fa0c1127SMark Adams */ 1272*fa0c1127SMark Adams if (user->checkVRes) { 1273*fa0c1127SMark Adams PetscBool pr = user->checkVRes > 1 ? PETSC_TRUE : PETSC_FALSE; 1274*fa0c1127SMark Adams PetscInt step; 1275*fa0c1127SMark Adams 1276*fa0c1127SMark Adams PetscCall(TSGetStepNumber(ts, &step)); 1277*fa0c1127SMark Adams if (pr) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "step: %" PetscInt_FMT "\n", step)); 1278*fa0c1127SMark Adams for (PetscInt p = 0; p < Np; ++p) { 1279*fa0c1127SMark Adams if (pr) PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "Residual: %.12g %.12g\n", (double)PetscRealPart(vres[p * dim + 0]), (double)PetscRealPart(vres[p * dim + 1]))); 1280*fa0c1127SMark Adams PetscCheck(PetscAbsScalar(vres[p * dim + 1]) < PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Y velocity should be 0., not %g", (double)PetscRealPart(vres[p * dim + 1])); 1281*fa0c1127SMark Adams } 1282*fa0c1127SMark Adams if (pr) PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 1283*fa0c1127SMark Adams } 1284*fa0c1127SMark Adams PetscCall(VecRestoreArray(Vres, &vres)); 1285*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, "E_field", NULL, NULL, (void **)&E)); 1286*fa0c1127SMark Adams PetscCall(PetscLogEventEnd(user->RhsVEvent, ts, 0, 0, 0)); 1287*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 1288*fa0c1127SMark Adams } 1289*fa0c1127SMark Adams 1290*fa0c1127SMark Adams /* Discrete Gradients Formulation: S, F, gradF (G) */ 1291*fa0c1127SMark Adams PetscErrorCode RHSJacobianS(TS ts, PetscReal t, Vec U, Mat S, void *ctx) 1292*fa0c1127SMark Adams { 1293*fa0c1127SMark Adams PetscScalar vals[4] = {0., 1., -1., 0.}; 1294*fa0c1127SMark Adams DM sw; 1295*fa0c1127SMark Adams PetscInt dim, d, Np, p, rStart; 1296*fa0c1127SMark Adams 1297*fa0c1127SMark Adams PetscFunctionBeginUser; 1298*fa0c1127SMark Adams PetscCall(TSGetDM(ts, &sw)); 1299*fa0c1127SMark Adams PetscCall(DMGetDimension(sw, &dim)); 1300*fa0c1127SMark Adams PetscCall(VecGetLocalSize(U, &Np)); 1301*fa0c1127SMark Adams PetscCall(MatGetOwnershipRange(S, &rStart, NULL)); 1302*fa0c1127SMark Adams Np /= 2 * dim; 1303*fa0c1127SMark Adams for (p = 0; p < Np; ++p) { 1304*fa0c1127SMark Adams for (d = 0; d < dim; ++d) { 1305*fa0c1127SMark Adams const PetscInt rows[2] = {(p * 2 + 0) * dim + d + rStart, (p * 2 + 1) * dim + d + rStart}; 1306*fa0c1127SMark Adams PetscCall(MatSetValues(S, 2, rows, 2, rows, vals, INSERT_VALUES)); 1307*fa0c1127SMark Adams } 1308*fa0c1127SMark Adams } 1309*fa0c1127SMark Adams PetscCall(MatAssemblyBegin(S, MAT_FINAL_ASSEMBLY)); 1310*fa0c1127SMark Adams PetscCall(MatAssemblyEnd(S, MAT_FINAL_ASSEMBLY)); 1311*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 1312*fa0c1127SMark Adams } 1313*fa0c1127SMark Adams 1314*fa0c1127SMark Adams PetscErrorCode RHSObjectiveF(TS ts, PetscReal t, Vec U, PetscScalar *F, void *ctx) 1315*fa0c1127SMark Adams { 1316*fa0c1127SMark Adams AppCtx *user = (AppCtx *)ctx; 1317*fa0c1127SMark Adams DM sw; 1318*fa0c1127SMark Adams Vec phi; 1319*fa0c1127SMark Adams const PetscScalar *u; 1320*fa0c1127SMark Adams PetscInt dim, Np, cStart, cEnd; 1321*fa0c1127SMark Adams PetscReal *vel, *coords, m_p = 1.; 1322*fa0c1127SMark Adams 1323*fa0c1127SMark Adams PetscFunctionBeginUser; 1324*fa0c1127SMark Adams PetscCall(TSGetDM(ts, &sw)); 1325*fa0c1127SMark Adams PetscCall(DMGetDimension(sw, &dim)); 1326*fa0c1127SMark Adams PetscCall(DMPlexGetHeightStratum(user->dmPot, 0, &cStart, &cEnd)); 1327*fa0c1127SMark Adams 1328*fa0c1127SMark Adams PetscCall(DMGetNamedGlobalVector(user->dmPot, "phi", &phi)); 1329*fa0c1127SMark Adams PetscCall(VecViewFromOptions(phi, NULL, "-phi_view_dg")); 1330*fa0c1127SMark Adams PetscCall(computeFieldEnergy(user->dmPot, phi, F)); 1331*fa0c1127SMark Adams PetscCall(DMRestoreNamedGlobalVector(user->dmPot, "phi", &phi)); 1332*fa0c1127SMark Adams 1333*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 1334*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&vel)); 1335*fa0c1127SMark Adams PetscCall(DMSwarmSortGetAccess(sw)); 1336*fa0c1127SMark Adams PetscCall(VecGetArrayRead(U, &u)); 1337*fa0c1127SMark Adams PetscCall(VecGetLocalSize(U, &Np)); 1338*fa0c1127SMark Adams Np /= 2 * dim; 1339*fa0c1127SMark Adams for (PetscInt c = cStart; c < cEnd; ++c) { 1340*fa0c1127SMark Adams PetscInt *points; 1341*fa0c1127SMark Adams PetscInt Ncp; 1342*fa0c1127SMark Adams 1343*fa0c1127SMark Adams PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Ncp, &points)); 1344*fa0c1127SMark Adams for (PetscInt cp = 0; cp < Ncp; ++cp) { 1345*fa0c1127SMark Adams const PetscInt p = points[cp]; 1346*fa0c1127SMark Adams const PetscReal v2 = DMPlex_DotRealD_Internal(dim, &u[(p * 2 + 1) * dim], &u[(p * 2 + 1) * dim]); 1347*fa0c1127SMark Adams 1348*fa0c1127SMark Adams *F += 0.5 * m_p * v2; 1349*fa0c1127SMark Adams } 1350*fa0c1127SMark Adams PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Ncp, &points)); 1351*fa0c1127SMark Adams } 1352*fa0c1127SMark Adams PetscCall(VecRestoreArrayRead(U, &u)); 1353*fa0c1127SMark Adams PetscCall(DMSwarmSortRestoreAccess(sw)); 1354*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 1355*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&vel)); 1356*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 1357*fa0c1127SMark Adams } 1358*fa0c1127SMark Adams 1359*fa0c1127SMark Adams /* dF/dx = q E dF/dv = v */ 1360*fa0c1127SMark Adams PetscErrorCode RHSFunctionG(TS ts, PetscReal t, Vec U, Vec G, void *ctx) 1361*fa0c1127SMark Adams { 1362*fa0c1127SMark Adams DM sw; 1363*fa0c1127SMark Adams SNES snes = ((AppCtx *)ctx)->snes; 1364*fa0c1127SMark Adams const PetscReal *coords, *vel, *E; 1365*fa0c1127SMark Adams const PetscScalar *u; 1366*fa0c1127SMark Adams PetscScalar *g; 1367*fa0c1127SMark Adams PetscReal m_p = 1., q_p = -1.; 1368*fa0c1127SMark Adams PetscInt dim, d, Np, p; 1369*fa0c1127SMark Adams 1370*fa0c1127SMark Adams PetscFunctionBeginUser; 1371*fa0c1127SMark Adams PetscCall(TSGetDM(ts, &sw)); 1372*fa0c1127SMark Adams PetscCall(DMGetDimension(sw, &dim)); 1373*fa0c1127SMark Adams PetscCall(DMSwarmGetLocalSize(sw, &Np)); 1374*fa0c1127SMark Adams PetscCall(VecGetArrayRead(U, &u)); 1375*fa0c1127SMark Adams PetscCall(VecGetArray(G, &g)); 1376*fa0c1127SMark Adams 1377*fa0c1127SMark Adams PetscLogEvent COMPUTEFIELD; 1378*fa0c1127SMark Adams PetscCall(PetscLogEventRegister("COMPFIELDATPART", TS_CLASSID, &COMPUTEFIELD)); 1379*fa0c1127SMark Adams PetscCall(PetscLogEventBegin(COMPUTEFIELD, 0, 0, 0, 0)); 1380*fa0c1127SMark Adams PetscCall(ComputeFieldAtParticles(snes, sw)); 1381*fa0c1127SMark Adams PetscCall(PetscLogEventEnd(COMPUTEFIELD, 0, 0, 0, 0)); 1382*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 1383*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&vel)); 1384*fa0c1127SMark Adams PetscCall(DMSwarmGetField(sw, "E_field", NULL, NULL, (void **)&E)); 1385*fa0c1127SMark Adams for (p = 0; p < Np; ++p) { 1386*fa0c1127SMark Adams for (d = 0; d < dim; ++d) { 1387*fa0c1127SMark Adams g[(p * 2 + 0) * dim + d] = -(q_p / m_p) * E[p * dim + d]; 1388*fa0c1127SMark Adams g[(p * 2 + 1) * dim + d] = u[(p * 2 + 1) * dim + d]; 1389*fa0c1127SMark Adams } 1390*fa0c1127SMark Adams } 1391*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, "E_field", NULL, NULL, (void **)&E)); 1392*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 1393*fa0c1127SMark Adams PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&vel)); 1394*fa0c1127SMark Adams PetscCall(VecRestoreArrayRead(U, &u)); 1395*fa0c1127SMark Adams PetscCall(VecRestoreArray(G, &g)); 1396*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 1397*fa0c1127SMark Adams } 1398*fa0c1127SMark Adams 1399*fa0c1127SMark Adams static PetscErrorCode CreateSolution(TS ts) 1400*fa0c1127SMark Adams { 1401*fa0c1127SMark Adams DM sw; 1402*fa0c1127SMark Adams Vec u; 1403*fa0c1127SMark Adams PetscInt dim, Np; 1404*fa0c1127SMark Adams 1405*fa0c1127SMark Adams PetscFunctionBegin; 1406*fa0c1127SMark Adams PetscCall(TSGetDM(ts, &sw)); 1407*fa0c1127SMark Adams PetscCall(DMGetDimension(sw, &dim)); 1408*fa0c1127SMark Adams PetscCall(DMSwarmGetLocalSize(sw, &Np)); 1409*fa0c1127SMark Adams PetscCall(VecCreate(PETSC_COMM_WORLD, &u)); 1410*fa0c1127SMark Adams PetscCall(VecSetBlockSize(u, dim)); 1411*fa0c1127SMark Adams PetscCall(VecSetSizes(u, 2 * Np * dim, PETSC_DECIDE)); 1412*fa0c1127SMark Adams PetscCall(VecSetUp(u)); 1413*fa0c1127SMark Adams PetscCall(TSSetSolution(ts, u)); 1414*fa0c1127SMark Adams PetscCall(VecDestroy(&u)); 1415*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 1416*fa0c1127SMark Adams } 1417*fa0c1127SMark Adams 1418*fa0c1127SMark Adams static PetscErrorCode SetProblem(TS ts) 1419*fa0c1127SMark Adams { 1420*fa0c1127SMark Adams AppCtx *user; 1421*fa0c1127SMark Adams DM sw; 1422*fa0c1127SMark Adams 1423*fa0c1127SMark Adams PetscFunctionBegin; 1424*fa0c1127SMark Adams PetscCall(TSGetDM(ts, &sw)); 1425*fa0c1127SMark Adams PetscCall(DMGetApplicationContext(sw, (void **)&user)); 1426*fa0c1127SMark Adams // Define unified system for (X, V) 1427*fa0c1127SMark Adams { 1428*fa0c1127SMark Adams Mat J; 1429*fa0c1127SMark Adams PetscInt dim, Np; 1430*fa0c1127SMark Adams 1431*fa0c1127SMark Adams PetscCall(DMGetDimension(sw, &dim)); 1432*fa0c1127SMark Adams PetscCall(DMSwarmGetLocalSize(sw, &Np)); 1433*fa0c1127SMark Adams PetscCall(MatCreate(PETSC_COMM_WORLD, &J)); 1434*fa0c1127SMark Adams PetscCall(MatSetSizes(J, 2 * Np * dim, 2 * Np * dim, PETSC_DECIDE, PETSC_DECIDE)); 1435*fa0c1127SMark Adams PetscCall(MatSetBlockSize(J, 2 * dim)); 1436*fa0c1127SMark Adams PetscCall(MatSetFromOptions(J)); 1437*fa0c1127SMark Adams PetscCall(MatSetUp(J)); 1438*fa0c1127SMark Adams PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, user)); 1439*fa0c1127SMark Adams PetscCall(TSSetRHSJacobian(ts, J, J, RHSJacobian, user)); 1440*fa0c1127SMark Adams PetscCall(MatDestroy(&J)); 1441*fa0c1127SMark Adams } 1442*fa0c1127SMark Adams /* Define split system for X and V */ 1443*fa0c1127SMark Adams { 1444*fa0c1127SMark Adams Vec u; 1445*fa0c1127SMark Adams IS isx, isv, istmp; 1446*fa0c1127SMark Adams const PetscInt *idx; 1447*fa0c1127SMark Adams PetscInt dim, Np, rstart; 1448*fa0c1127SMark Adams 1449*fa0c1127SMark Adams PetscCall(TSGetSolution(ts, &u)); 1450*fa0c1127SMark Adams PetscCall(DMGetDimension(sw, &dim)); 1451*fa0c1127SMark Adams PetscCall(DMSwarmGetLocalSize(sw, &Np)); 1452*fa0c1127SMark Adams PetscCall(VecGetOwnershipRange(u, &rstart, NULL)); 1453*fa0c1127SMark Adams PetscCall(ISCreateStride(PETSC_COMM_WORLD, Np, (rstart / dim) + 0, 2, &istmp)); 1454*fa0c1127SMark Adams PetscCall(ISGetIndices(istmp, &idx)); 1455*fa0c1127SMark Adams PetscCall(ISCreateBlock(PETSC_COMM_WORLD, dim, Np, idx, PETSC_COPY_VALUES, &isx)); 1456*fa0c1127SMark Adams PetscCall(ISRestoreIndices(istmp, &idx)); 1457*fa0c1127SMark Adams PetscCall(ISDestroy(&istmp)); 1458*fa0c1127SMark Adams PetscCall(ISCreateStride(PETSC_COMM_WORLD, Np, (rstart / dim) + 1, 2, &istmp)); 1459*fa0c1127SMark Adams PetscCall(ISGetIndices(istmp, &idx)); 1460*fa0c1127SMark Adams PetscCall(ISCreateBlock(PETSC_COMM_WORLD, dim, Np, idx, PETSC_COPY_VALUES, &isv)); 1461*fa0c1127SMark Adams PetscCall(ISRestoreIndices(istmp, &idx)); 1462*fa0c1127SMark Adams PetscCall(ISDestroy(&istmp)); 1463*fa0c1127SMark Adams PetscCall(TSRHSSplitSetIS(ts, "position", isx)); 1464*fa0c1127SMark Adams PetscCall(TSRHSSplitSetIS(ts, "momentum", isv)); 1465*fa0c1127SMark Adams PetscCall(ISDestroy(&isx)); 1466*fa0c1127SMark Adams PetscCall(ISDestroy(&isv)); 1467*fa0c1127SMark Adams PetscCall(TSRHSSplitSetRHSFunction(ts, "position", NULL, RHSFunctionX, user)); 1468*fa0c1127SMark Adams PetscCall(TSRHSSplitSetRHSFunction(ts, "momentum", NULL, RHSFunctionV, user)); 1469*fa0c1127SMark Adams } 1470*fa0c1127SMark Adams // Define symplectic formulation U_t = S . G, where G = grad F 1471*fa0c1127SMark Adams { 1472*fa0c1127SMark Adams PetscCall(TSDiscGradSetFormulation(ts, RHSJacobianS, RHSObjectiveF, RHSFunctionG, user)); 1473*fa0c1127SMark Adams } 1474*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 1475*fa0c1127SMark Adams } 1476*fa0c1127SMark Adams 1477*fa0c1127SMark Adams static PetscErrorCode DMSwarmTSRedistribute(TS ts) 1478*fa0c1127SMark Adams { 1479*fa0c1127SMark Adams DM sw; 1480*fa0c1127SMark Adams Vec u; 1481*fa0c1127SMark Adams PetscReal t, maxt, dt; 1482*fa0c1127SMark Adams PetscInt n, maxn; 1483*fa0c1127SMark Adams 1484*fa0c1127SMark Adams PetscFunctionBegin; 1485*fa0c1127SMark Adams PetscCall(TSGetDM(ts, &sw)); 1486*fa0c1127SMark Adams PetscCall(TSGetTime(ts, &t)); 1487*fa0c1127SMark Adams PetscCall(TSGetMaxTime(ts, &maxt)); 1488*fa0c1127SMark Adams PetscCall(TSGetTimeStep(ts, &dt)); 1489*fa0c1127SMark Adams PetscCall(TSGetStepNumber(ts, &n)); 1490*fa0c1127SMark Adams PetscCall(TSGetMaxSteps(ts, &maxn)); 1491*fa0c1127SMark Adams 1492*fa0c1127SMark Adams PetscCall(TSReset(ts)); 1493*fa0c1127SMark Adams PetscCall(TSSetDM(ts, sw)); 1494*fa0c1127SMark Adams PetscCall(TSSetFromOptions(ts)); 1495*fa0c1127SMark Adams PetscCall(TSSetTime(ts, t)); 1496*fa0c1127SMark Adams PetscCall(TSSetMaxTime(ts, maxt)); 1497*fa0c1127SMark Adams PetscCall(TSSetTimeStep(ts, dt)); 1498*fa0c1127SMark Adams PetscCall(TSSetStepNumber(ts, n)); 1499*fa0c1127SMark Adams PetscCall(TSSetMaxSteps(ts, maxn)); 1500*fa0c1127SMark Adams 1501*fa0c1127SMark Adams PetscCall(CreateSolution(ts)); 1502*fa0c1127SMark Adams PetscCall(SetProblem(ts)); 1503*fa0c1127SMark Adams PetscCall(TSGetSolution(ts, &u)); 1504*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 1505*fa0c1127SMark Adams } 1506*fa0c1127SMark Adams 1507*fa0c1127SMark Adams PetscErrorCode line(PetscInt dim, PetscReal time, const PetscReal dummy[], PetscInt p, PetscScalar x[], void *ctx) 1508*fa0c1127SMark Adams { 1509*fa0c1127SMark Adams DM sw, cdm; 1510*fa0c1127SMark Adams PetscInt Np; 1511*fa0c1127SMark Adams PetscReal low[2], high[2]; 1512*fa0c1127SMark Adams AppCtx *user = (AppCtx *)ctx; 1513*fa0c1127SMark Adams 1514*fa0c1127SMark Adams sw = user->swarm; 1515*fa0c1127SMark Adams PetscCall(DMSwarmGetCellDM(sw, &cdm)); 1516*fa0c1127SMark Adams // Get the bounding box so we can equally space the particles 1517*fa0c1127SMark Adams PetscCall(DMGetLocalBoundingBox(cdm, low, high)); 1518*fa0c1127SMark Adams PetscCall(DMSwarmGetLocalSize(sw, &Np)); 1519*fa0c1127SMark Adams // shift it by h/2 so nothing is initialized directly on a boundary 1520*fa0c1127SMark Adams x[0] = ((high[0] - low[0]) / Np) * (p + 0.5); 1521*fa0c1127SMark Adams x[1] = 0.; 1522*fa0c1127SMark Adams return PETSC_SUCCESS; 1523*fa0c1127SMark Adams } 1524*fa0c1127SMark Adams 1525*fa0c1127SMark Adams /* 1526*fa0c1127SMark Adams InitializeSolveAndSwarm - Set the solution values to the swarm coordinates and velocities, and also possibly set the initial values. 1527*fa0c1127SMark Adams 1528*fa0c1127SMark Adams Input Parameters: 1529*fa0c1127SMark Adams + ts - The TS 1530*fa0c1127SMark Adams - useInitial - Flag to also set the initial conditions to the current coordinates and velocities and setup the problem 1531*fa0c1127SMark Adams 1532*fa0c1127SMark Adams Output Parameters: 1533*fa0c1127SMark Adams . u - The initialized solution vector 1534*fa0c1127SMark Adams 1535*fa0c1127SMark Adams Level: advanced 1536*fa0c1127SMark Adams 1537*fa0c1127SMark Adams .seealso: InitializeSolve() 1538*fa0c1127SMark Adams */ 1539*fa0c1127SMark Adams static PetscErrorCode InitializeSolveAndSwarm(TS ts, PetscBool useInitial) 1540*fa0c1127SMark Adams { 1541*fa0c1127SMark Adams DM sw; 1542*fa0c1127SMark Adams Vec u, gc, gv; 1543*fa0c1127SMark Adams IS isx, isv; 1544*fa0c1127SMark Adams PetscInt dim; 1545*fa0c1127SMark Adams AppCtx *user; 1546*fa0c1127SMark Adams 1547*fa0c1127SMark Adams PetscFunctionBeginUser; 1548*fa0c1127SMark Adams PetscCall(TSGetDM(ts, &sw)); 1549*fa0c1127SMark Adams PetscCall(DMGetApplicationContext(sw, &user)); 1550*fa0c1127SMark Adams PetscCall(DMGetDimension(sw, &dim)); 1551*fa0c1127SMark Adams if (useInitial) { 1552*fa0c1127SMark Adams PetscCall(InitializeParticles_PerturbedWeights(sw, user)); 1553*fa0c1127SMark Adams PetscCall(DMSwarmMigrate(sw, PETSC_TRUE)); 1554*fa0c1127SMark Adams PetscCall(DMSwarmTSRedistribute(ts)); 1555*fa0c1127SMark Adams } 1556*fa0c1127SMark Adams PetscCall(DMSetUp(sw)); 1557*fa0c1127SMark Adams PetscCall(TSGetSolution(ts, &u)); 1558*fa0c1127SMark Adams PetscCall(TSRHSSplitGetIS(ts, "position", &isx)); 1559*fa0c1127SMark Adams PetscCall(TSRHSSplitGetIS(ts, "momentum", &isv)); 1560*fa0c1127SMark Adams PetscCall(DMSwarmCreateGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc)); 1561*fa0c1127SMark Adams PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "velocity", &gv)); 1562*fa0c1127SMark Adams PetscCall(VecISCopy(u, isx, SCATTER_FORWARD, gc)); 1563*fa0c1127SMark Adams PetscCall(VecISCopy(u, isv, SCATTER_FORWARD, gv)); 1564*fa0c1127SMark Adams PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc)); 1565*fa0c1127SMark Adams PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "velocity", &gv)); 1566*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 1567*fa0c1127SMark Adams } 1568*fa0c1127SMark Adams 1569*fa0c1127SMark Adams static PetscErrorCode InitializeSolve(TS ts, Vec u) 1570*fa0c1127SMark Adams { 1571*fa0c1127SMark Adams PetscFunctionBegin; 1572*fa0c1127SMark Adams PetscCall(TSSetSolution(ts, u)); 1573*fa0c1127SMark Adams PetscCall(InitializeSolveAndSwarm(ts, PETSC_TRUE)); 1574*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 1575*fa0c1127SMark Adams } 1576*fa0c1127SMark Adams 1577*fa0c1127SMark Adams static PetscErrorCode MigrateParticles(TS ts) 1578*fa0c1127SMark Adams { 1579*fa0c1127SMark Adams DM sw, cdm; 1580*fa0c1127SMark Adams const PetscReal *L; 1581*fa0c1127SMark Adams AppCtx *ctx; 1582*fa0c1127SMark Adams 1583*fa0c1127SMark Adams PetscFunctionBeginUser; 1584*fa0c1127SMark Adams PetscCall(TSGetDM(ts, &sw)); 1585*fa0c1127SMark Adams PetscCall(DMGetApplicationContext(sw, &ctx)); 1586*fa0c1127SMark Adams PetscCall(DMViewFromOptions(sw, NULL, "-migrate_view_pre")); 1587*fa0c1127SMark Adams { 1588*fa0c1127SMark Adams Vec u, gc, gv, position, momentum; 1589*fa0c1127SMark Adams IS isx, isv; 1590*fa0c1127SMark Adams PetscReal *pos, *mom; 1591*fa0c1127SMark Adams 1592*fa0c1127SMark Adams PetscCall(TSGetSolution(ts, &u)); 1593*fa0c1127SMark Adams PetscCall(TSRHSSplitGetIS(ts, "position", &isx)); 1594*fa0c1127SMark Adams PetscCall(TSRHSSplitGetIS(ts, "momentum", &isv)); 1595*fa0c1127SMark Adams PetscCall(VecGetSubVector(u, isx, &position)); 1596*fa0c1127SMark Adams PetscCall(VecGetSubVector(u, isv, &momentum)); 1597*fa0c1127SMark Adams PetscCall(VecGetArray(position, &pos)); 1598*fa0c1127SMark Adams PetscCall(VecGetArray(momentum, &mom)); 1599*fa0c1127SMark Adams PetscCall(DMSwarmCreateGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc)); 1600*fa0c1127SMark Adams PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "velocity", &gv)); 1601*fa0c1127SMark Adams PetscCall(VecISCopy(u, isx, SCATTER_REVERSE, gc)); 1602*fa0c1127SMark Adams PetscCall(VecISCopy(u, isv, SCATTER_REVERSE, gv)); 1603*fa0c1127SMark Adams 1604*fa0c1127SMark Adams PetscCall(DMSwarmGetCellDM(sw, &cdm)); 1605*fa0c1127SMark Adams PetscCall(DMGetPeriodicity(cdm, NULL, NULL, &L)); 1606*fa0c1127SMark Adams if ((L[0] || L[1]) >= 0.) { 1607*fa0c1127SMark Adams PetscReal *x, *v, upper[3], lower[3]; 1608*fa0c1127SMark Adams PetscInt Np, dim; 1609*fa0c1127SMark Adams 1610*fa0c1127SMark Adams PetscCall(DMSwarmGetLocalSize(sw, &Np)); 1611*fa0c1127SMark Adams PetscCall(DMGetDimension(cdm, &dim)); 1612*fa0c1127SMark Adams PetscCall(DMGetBoundingBox(cdm, lower, upper)); 1613*fa0c1127SMark Adams PetscCall(VecGetArray(gc, &x)); 1614*fa0c1127SMark Adams PetscCall(VecGetArray(gv, &v)); 1615*fa0c1127SMark Adams for (PetscInt p = 0; p < Np; ++p) { 1616*fa0c1127SMark Adams for (PetscInt d = 0; d < dim; ++d) { 1617*fa0c1127SMark Adams if (pos[p * dim + d] < lower[d]) { 1618*fa0c1127SMark Adams x[p * dim + d] = pos[p * dim + d] + (upper[d] - lower[d]); 1619*fa0c1127SMark Adams } else if (pos[p * dim + d] > upper[d]) { 1620*fa0c1127SMark Adams x[p * dim + d] = pos[p * dim + d] - (upper[d] - lower[d]); 1621*fa0c1127SMark Adams } else { 1622*fa0c1127SMark Adams x[p * dim + d] = pos[p * dim + d]; 1623*fa0c1127SMark Adams } 1624*fa0c1127SMark Adams PetscCheck(x[p * dim + d] >= lower[d] && x[p * dim + d] <= upper[d], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "p: %" PetscInt_FMT "x[%" PetscInt_FMT "] %g", p, d, (double)x[p * dim + d]); 1625*fa0c1127SMark Adams v[p * dim + d] = mom[p * dim + d]; 1626*fa0c1127SMark Adams } 1627*fa0c1127SMark Adams } 1628*fa0c1127SMark Adams PetscCall(VecRestoreArray(gc, &x)); 1629*fa0c1127SMark Adams PetscCall(VecRestoreArray(gv, &v)); 1630*fa0c1127SMark Adams } 1631*fa0c1127SMark Adams PetscCall(VecRestoreArray(position, &pos)); 1632*fa0c1127SMark Adams PetscCall(VecRestoreArray(momentum, &mom)); 1633*fa0c1127SMark Adams PetscCall(VecRestoreSubVector(u, isx, &position)); 1634*fa0c1127SMark Adams PetscCall(VecRestoreSubVector(u, isv, &momentum)); 1635*fa0c1127SMark Adams PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "velocity", &gv)); 1636*fa0c1127SMark Adams PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc)); 1637*fa0c1127SMark Adams } 1638*fa0c1127SMark Adams PetscCall(DMSwarmMigrate(sw, PETSC_TRUE)); 1639*fa0c1127SMark Adams // This MUST come last, since it recreates the subswarms and they must DMClone() the new swarm 1640*fa0c1127SMark Adams PetscCall(DMSwarmTSRedistribute(ts)); 1641*fa0c1127SMark Adams PetscCall(InitializeSolveAndSwarm(ts, PETSC_FALSE)); 1642*fa0c1127SMark Adams { 1643*fa0c1127SMark Adams const char *fieldnames[2] = {DMSwarmPICField_coor, "velocity"}; 1644*fa0c1127SMark Adams PetscCall(DMSwarmVectorDefineFields(sw, 2, fieldnames)); 1645*fa0c1127SMark Adams } 1646*fa0c1127SMark Adams PetscFunctionReturn(PETSC_SUCCESS); 1647*fa0c1127SMark Adams } 1648*fa0c1127SMark Adams 1649*fa0c1127SMark Adams int main(int argc, char **argv) 1650*fa0c1127SMark Adams { 1651*fa0c1127SMark Adams DM dm, sw; 1652*fa0c1127SMark Adams TS ts; 1653*fa0c1127SMark Adams Vec u; 1654*fa0c1127SMark Adams PetscReal dt; 1655*fa0c1127SMark Adams PetscInt maxn; 1656*fa0c1127SMark Adams AppCtx user; 1657*fa0c1127SMark Adams 1658*fa0c1127SMark Adams PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1659*fa0c1127SMark Adams PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 1660*fa0c1127SMark Adams PetscCall(PetscBagCreate(PETSC_COMM_SELF, sizeof(Parameter), &user.bag)); 1661*fa0c1127SMark Adams PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 1662*fa0c1127SMark Adams PetscCall(CreatePoisson(dm, &user)); 1663*fa0c1127SMark Adams PetscCall(CreateSwarm(dm, &user, &sw)); 1664*fa0c1127SMark Adams PetscCall(SetupParameters(PETSC_COMM_WORLD, &user)); 1665*fa0c1127SMark Adams PetscCall(InitializeConstants(sw, &user)); 1666*fa0c1127SMark Adams PetscCall(DMSetApplicationContext(sw, &user)); 1667*fa0c1127SMark Adams 1668*fa0c1127SMark Adams PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1669*fa0c1127SMark Adams PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 1670*fa0c1127SMark Adams PetscCall(TSSetDM(ts, sw)); 1671*fa0c1127SMark Adams PetscCall(TSSetMaxTime(ts, 0.1)); 1672*fa0c1127SMark Adams PetscCall(TSSetTimeStep(ts, 0.00001)); 1673*fa0c1127SMark Adams PetscCall(TSSetMaxSteps(ts, 100)); 1674*fa0c1127SMark Adams PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 1675*fa0c1127SMark Adams 1676*fa0c1127SMark Adams if (user.efield_monitor) PetscCall(TSMonitorSet(ts, MonitorEField, &user, NULL)); 1677*fa0c1127SMark Adams if (user.moment_monitor) PetscCall(TSMonitorSet(ts, MonitorMoments, &user, NULL)); 1678*fa0c1127SMark Adams if (user.particle_monitor) PetscCall(TSMonitorSet(ts, MonitorParticles, &user, NULL)); 1679*fa0c1127SMark Adams 1680*fa0c1127SMark Adams PetscCall(TSSetFromOptions(ts)); 1681*fa0c1127SMark Adams PetscCall(TSGetTimeStep(ts, &dt)); 1682*fa0c1127SMark Adams PetscCall(TSGetMaxSteps(ts, &maxn)); 1683*fa0c1127SMark Adams user.steps = maxn; 1684*fa0c1127SMark Adams user.stepSize = dt; 1685*fa0c1127SMark Adams PetscCall(SetupContext(dm, sw, &user)); 1686*fa0c1127SMark Adams PetscCall(TSSetComputeInitialCondition(ts, InitializeSolve)); 1687*fa0c1127SMark Adams PetscCall(TSSetPostStep(ts, MigrateParticles)); 1688*fa0c1127SMark Adams PetscCall(CreateSolution(ts)); 1689*fa0c1127SMark Adams PetscCall(TSGetSolution(ts, &u)); 1690*fa0c1127SMark Adams PetscCall(TSComputeInitialCondition(ts, u)); 1691*fa0c1127SMark Adams PetscCall(CheckNonNegativeWeights(sw, &user)); 1692*fa0c1127SMark Adams PetscCall(TSSolve(ts, NULL)); 1693*fa0c1127SMark Adams 1694*fa0c1127SMark Adams PetscCall(SNESDestroy(&user.snes)); 1695*fa0c1127SMark Adams PetscCall(DMDestroy(&user.dmPot)); 1696*fa0c1127SMark Adams PetscCall(MatDestroy(&user.M)); 1697*fa0c1127SMark Adams PetscCall(PetscFEGeomDestroy(&user.fegeom)); 1698*fa0c1127SMark Adams PetscCall(TSDestroy(&ts)); 1699*fa0c1127SMark Adams PetscCall(DMDestroy(&sw)); 1700*fa0c1127SMark Adams PetscCall(DMDestroy(&dm)); 1701*fa0c1127SMark Adams PetscCall(DestroyContext(&user)); 1702*fa0c1127SMark Adams PetscCall(PetscFinalize()); 1703*fa0c1127SMark Adams return 0; 1704*fa0c1127SMark Adams } 1705*fa0c1127SMark Adams 1706*fa0c1127SMark Adams /*TEST 1707*fa0c1127SMark Adams 1708*fa0c1127SMark Adams build: 1709*fa0c1127SMark Adams requires: !complex double 1710*fa0c1127SMark Adams 1711*fa0c1127SMark Adams testset: 1712*fa0c1127SMark Adams nsize: 2 1713*fa0c1127SMark Adams args: -cosine_coefficients 0.01 -charges -1. -total_weight 1. -xdm_plex_hash_location -vpetscspace_degree 2 -petscspace_degree 1 -em_snes_atol 1.e-12 \ 1714*fa0c1127SMark Adams -em_snes_error_if_not_converged -em_ksp_error_if_not_converged -em_ksp_type cg -em_pc_type gamg -em_mg_coarse_ksp_type preonly -em_mg_coarse_pc_type svd -em_proj_ksp_type cg \ 1715*fa0c1127SMark Adams -em_proj_pc_type gamg -em_proj_mg_coarse_ksp_type preonly -em_proj_mg_coarse_pc_type svd -ts_time_step 0.03 -xdm_plex_simplex 0 \ 1716*fa0c1127SMark Adams -ts_max_time 100 -output_step 1 -ts_type basicsymplectic -ts_basicsymplectic_type 1 -check_weights -ts_max_steps 60 1717*fa0c1127SMark Adams 1718*fa0c1127SMark Adams test: 1719*fa0c1127SMark Adams suffix: landau_damping_1d 1720*fa0c1127SMark Adams args: -xdm_plex_dim 1 -xdm_plex_box_faces 60 -xdm_plex_box_lower 0. -xdm_plex_box_upper 12.5664 -xdm_plex_box_bd periodic -vdm_plex_dim 1 -vdm_plex_box_faces 60 \ 1721*fa0c1127SMark Adams -vdm_plex_box_lower -6 -vdm_plex_box_upper 6 -vdm_plex_hash_location -vdm_plex_simplex 0 -vdm_plex_box_bd none 1722*fa0c1127SMark Adams 1723*fa0c1127SMark Adams test: 1724*fa0c1127SMark Adams suffix: landau_damping_2d 1725*fa0c1127SMark Adams args: -xdm_plex_dim 2 -xdm_plex_box_bd periodic,periodic -vdm_plex_dim 2 -xdm_plex_box_lower 0.,-.5 -vdm_plex_box_lower -6,-6 -vdm_plex_box_upper 6,6 -xdm_plex_box_faces 6,3 \ 1726*fa0c1127SMark Adams -xdm_plex_box_upper 12.5664,.5 -vdm_plex_box_faces 15,15 -vdm_plex_box_bd none,none -vdm_plex_hash_location -vdm_plex_simplex 0 1727*fa0c1127SMark Adams 1728*fa0c1127SMark Adams test: 1729*fa0c1127SMark Adams suffix: landau_damping_3d 1730*fa0c1127SMark Adams args: -xdm_plex_dim 3 -xdm_plex_box_faces 6,3,3 -xdm_plex_box_lower 0.,-1,-1 -xdm_plex_box_upper 12.5664,1,1 -xdm_plex_box_bd periodic,periodic,periodic -vdm_plex_dim 3 -vdm_plex_box_faces 4,4,4 \ 1731*fa0c1127SMark Adams -vdm_plex_box_lower -6,-6,-6 -vdm_plex_box_upper 6,6,6 -vdm_plex_hash_location -vdm_plex_simplex 0 -vdm_plex_box_bd none,none,none 1732*fa0c1127SMark Adams 1733*fa0c1127SMark Adams test: 1734*fa0c1127SMark Adams requires: !defined(PETSC_USE_DMLANDAU_2D) 1735*fa0c1127SMark Adams suffix: sphere_3d 1736*fa0c1127SMark Adams nsize: 1 1737*fa0c1127SMark Adams args: -use_landau_velocity_space -xdm_plex_dim 3 -vdm_landau_thermal_temps 1 -vdm_landau_device_type cpu -xdm_plex_box_faces 6,3,3 -xdm_plex_box_lower 0.,-1,-1 \ 1738*fa0c1127SMark Adams -xdm_plex_box_upper 12.5664,1,1 -xdm_plex_box_bd periodic,periodic,periodic -vdm_landau_verbose 2 -vdm_landau_sphere -vdm_landau_map_sphere \ 1739*fa0c1127SMark Adams -vdm_landau_domain_radius 6,6,6 -vdm_landau_sphere_inner_radius_90degree_scale .35 -vdm_refine 1 1740*fa0c1127SMark Adams 1741*fa0c1127SMark Adams TEST*/ 1742