1 static char help[] = "Landau Damping/Two Stream instability test using Vlasov-Poisson equations\n"; 2 3 /* 4 TODO: 5 - Cache mesh geometry 6 - Move electrostatic solver to MG+SVD 7 8 To run the code with particles sinusoidally perturbed in x space use the test "pp_poisson_bsi_1d_4" or "pp_poisson_bsi_2d_4" 9 According to Lukas, good damping results come at ~16k particles per cell 10 11 To visualize the maximum electric field use 12 13 -efield_monitor 14 15 To monitor velocity moments of the distribution use 16 17 -ptof_pc_type lu -moments_monitor 18 19 To monitor the particle positions in phase space use 20 21 -positions_monitor 22 23 To monitor the charge density, E field, and potential use 24 25 -poisson_monitor 26 27 To monitor the remapping field use 28 29 -remap_uf_view draw 30 31 To visualize the swarm distribution use 32 33 -ts_monitor_hg_swarm 34 35 To visualize the particles, we can use 36 37 -ts_monitor_sp_swarm -ts_monitor_sp_swarm_retain 0 -ts_monitor_sp_swarm_phase 1 -draw_size 500,500 38 39 For a Landau Damping verification run, we use 40 41 -dm_plex_dim 2 -fake_1D -dm_plex_simplex 0 -dm_plex_box_faces 10,1 \ 42 -dm_plex_box_lower 0.,-0.5 -dm_plex_box_upper 12.5664,0.5 -dm_plex_box_bd periodic,none \ 43 -vdm_plex_dim 1 -vdm_plex_simplex 0 -vdm_plex_box_faces 2000 -vdm_plex_box_lower -10 -vdm_plex_box_upper 10 \ 44 -dm_swarm_num_species 1 -charges -1.,1. \ 45 -cosine_coefficients 0.01,0.5 -perturbed_weights -total_weight 1. \ 46 -ts_type basicsymplectic -ts_basicsymplectic_type 1 -ts_dt 0.03 -ts_max_time 500 -ts_max_steps 500 \ 47 -em_type primal -petscspace_degree 1 -em_snes_atol 1.e-12 -em_snes_error_if_not_converged -em_ksp_error_if_not_converged -em_pc_type svd \ 48 -output_step 100 -check_vel_res -efield_monitor -ts_monitor -log_view 49 50 */ 51 #include <petscts.h> 52 #include <petscdmplex.h> 53 #include <petscdmswarm.h> 54 #include <petscfe.h> 55 #include <petscds.h> 56 #include <petscbag.h> 57 #include <petscdraw.h> 58 #include <petsc/private/dmpleximpl.h> /* For norm and dot */ 59 #include <petsc/private/petscfeimpl.h> /* For interpolation */ 60 #include <petsc/private/dmswarmimpl.h> /* For swarm debugging */ 61 #include "petscdm.h" 62 #include "petscdmlabel.h" 63 64 PETSC_EXTERN PetscErrorCode stream(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *); 65 PETSC_EXTERN PetscErrorCode line(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *); 66 67 const char *EMTypes[] = {"primal", "mixed", "coulomb", "none", "EMType", "EM_", NULL}; 68 typedef enum { 69 EM_PRIMAL, 70 EM_MIXED, 71 EM_COULOMB, 72 EM_NONE 73 } EMType; 74 75 const char *RemapTypes[] = {"colella", "pfak", "RemapType", "RM_", NULL}; 76 typedef enum { 77 RM_COLELLA, 78 RM_PFAK, 79 } RemapType; 80 81 typedef enum { 82 V0, 83 X0, 84 T0, 85 M0, 86 Q0, 87 PHI0, 88 POISSON, 89 VLASOV, 90 SIGMA, 91 NUM_CONSTANTS 92 } ConstantType; 93 typedef struct { 94 PetscScalar v0; /* Velocity scale, often the thermal velocity */ 95 PetscScalar t0; /* Time scale */ 96 PetscScalar x0; /* Space scale */ 97 PetscScalar m0; /* Mass scale */ 98 PetscScalar q0; /* Charge scale */ 99 PetscScalar kb; 100 PetscScalar epsi0; 101 PetscScalar phi0; /* Potential scale */ 102 PetscScalar poissonNumber; /* Non-Dimensional Poisson Number */ 103 PetscScalar vlasovNumber; /* Non-Dimensional Vlasov Number */ 104 PetscReal sigma; /* Nondimensional charge per length in x */ 105 } Parameter; 106 107 typedef struct { 108 PetscBag bag; // Problem parameters 109 PetscBool error; // Flag for printing the error 110 PetscBool remap; // Flag to remap particles 111 RemapType remapType; // Remap algorithm 112 PetscInt remapFreq; // Number of timesteps between remapping 113 PetscBool efield_monitor; // Flag to show electric field monitor 114 PetscBool moment_monitor; // Flag to show distribution moment monitor 115 PetscBool weight_monitor; // Flag to show weight monitor 116 PetscBool positions_monitor; // Flag to show particle positins at each time step 117 PetscBool poisson_monitor; // Flag to display charge, E field, and potential at each solve 118 PetscBool initial_monitor; // Flag to monitor the initial conditions 119 PetscBool fake_1D; // Run simulation in 2D but zeroing second dimension 120 PetscBool perturbed_weights; // Uniformly sample x,v space with gaussian weights 121 PetscInt ostep; // Print the energy at each ostep time steps 122 PetscInt numParticles; 123 PetscReal timeScale; /* Nondimensionalizing time scale */ 124 PetscReal charges[2]; /* The charges of each species */ 125 PetscReal masses[2]; /* The masses of each species */ 126 PetscReal thermal_energy[2]; /* Thermal Energy (used to get other constants)*/ 127 PetscReal cosine_coefficients[2]; /*(alpha, k)*/ 128 PetscReal totalWeight; 129 PetscReal stepSize; 130 PetscInt steps; 131 PetscReal initVel; 132 EMType em; // Type of electrostatic model 133 SNES snes; // EM solver 134 Mat M; // The finite element mass matrix 135 PetscFEGeom *fegeom; // Geometric information for the DM cells 136 PetscDraw drawic_x; 137 PetscDraw drawic_v; 138 PetscDraw drawic_w; 139 PetscDrawHG drawhgic_x; 140 PetscDrawHG drawhgic_v; 141 PetscDrawHG drawhgic_w; 142 PetscBool validE; // Flag to indicate E-field in swarm is valid 143 PetscReal drawlgEmin; // The minimum lg(E) to plot 144 PetscDrawLG drawlgE; // Logarithm of maximum electric field 145 PetscDrawSP drawspE; // Electric field at particle positions 146 PetscDrawSP drawspX; // Particle positions 147 PetscViewer viewerRho; // Charge density viewer 148 PetscViewer viewerPhi; // Potential viewer 149 DM swarm; 150 PetscRandom random; 151 PetscBool twostream; 152 PetscBool checkweights; 153 PetscInt checkVRes; /* Flag to check/output velocity residuals for nightly tests */ 154 155 PetscLogEvent RhsXEvent, RhsVEvent, ESolveEvent, ETabEvent; 156 } AppCtx; 157 158 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 159 { 160 PetscFunctionBeginUser; 161 PetscInt d = 2; 162 PetscInt maxSpecies = 2; 163 options->error = PETSC_FALSE; 164 options->remap = PETSC_FALSE; 165 options->remapType = RM_COLELLA; 166 options->remapFreq = 1; 167 options->efield_monitor = PETSC_FALSE; 168 options->moment_monitor = PETSC_FALSE; 169 options->weight_monitor = PETSC_FALSE; 170 options->initial_monitor = PETSC_FALSE; 171 options->fake_1D = PETSC_FALSE; 172 options->perturbed_weights = PETSC_FALSE; 173 options->poisson_monitor = PETSC_FALSE; 174 options->positions_monitor = PETSC_FALSE; 175 options->ostep = 100; 176 options->timeScale = 2.0e-14; 177 options->charges[0] = -1.0; 178 options->charges[1] = 1.0; 179 options->masses[0] = 1.0; 180 options->masses[1] = 1000.0; 181 options->thermal_energy[0] = 1.0; 182 options->thermal_energy[1] = 1.0; 183 options->cosine_coefficients[0] = 0.01; 184 options->cosine_coefficients[1] = 0.5; 185 options->initVel = 1; 186 options->totalWeight = 1.0; 187 options->drawic_x = NULL; 188 options->drawic_v = NULL; 189 options->drawic_w = NULL; 190 options->drawhgic_x = NULL; 191 options->drawhgic_v = NULL; 192 options->drawhgic_w = NULL; 193 options->drawlgEmin = -6; 194 options->drawlgE = NULL; 195 options->drawspE = NULL; 196 options->drawspX = NULL; 197 options->viewerRho = NULL; 198 options->viewerPhi = NULL; 199 options->em = EM_COULOMB; 200 options->numParticles = 32768; 201 options->twostream = PETSC_FALSE; 202 options->checkweights = PETSC_FALSE; 203 options->checkVRes = 0; 204 205 PetscOptionsBegin(comm, "", "Landau Damping and Two Stream options", "DMSWARM"); 206 PetscCall(PetscOptionsBool("-error", "Flag to print the error", "ex2.c", options->error, &options->error, NULL)); 207 PetscCall(PetscOptionsBool("-remap", "Flag to remap the particles", "ex2.c", options->remap, &options->remap, NULL)); 208 PetscCall(PetscOptionsInt("-remap_freq", "Number", "ex2.c", options->remapFreq, &options->remapFreq, NULL)); 209 PetscCall(PetscOptionsEnum("-remap_type", "Remap algorithm", "ex2.c", RemapTypes, (PetscEnum)options->remapType, (PetscEnum *)&options->remapType, NULL)); 210 PetscCall(PetscOptionsBool("-efield_monitor", "Flag to plot log(max E) over time", "ex2.c", options->efield_monitor, &options->efield_monitor, NULL)); 211 PetscCall(PetscOptionsReal("-efield_min_monitor", "Minimum E field to plot", "ex2.c", options->drawlgEmin, &options->drawlgEmin, NULL)); 212 PetscCall(PetscOptionsBool("-moments_monitor", "Flag to show moments table", "ex2.c", options->moment_monitor, &options->moment_monitor, NULL)); 213 PetscCall(PetscOptionsBool("-weight_monitor", "The flag to show particle weights", "ex2.c", options->weight_monitor, &options->weight_monitor, NULL)); 214 PetscCall(PetscOptionsBool("-ics_monitor", "Flag to show initial condition histograms", "ex2.c", options->initial_monitor, &options->initial_monitor, NULL)); 215 PetscCall(PetscOptionsBool("-positions_monitor", "The flag to show particle positions", "ex2.c", options->positions_monitor, &options->positions_monitor, NULL)); 216 PetscCall(PetscOptionsBool("-poisson_monitor", "The flag to show charges, Efield and potential solve", "ex2.c", options->poisson_monitor, &options->poisson_monitor, NULL)); 217 PetscCall(PetscOptionsBool("-fake_1D", "Flag to run a 1D simulation (but really in 2D)", "ex2.c", options->fake_1D, &options->fake_1D, NULL)); 218 PetscCall(PetscOptionsBool("-twostream", "Run two stream instability", "ex2.c", options->twostream, &options->twostream, NULL)); 219 PetscCall(PetscOptionsBool("-perturbed_weights", "Flag to run uniform sampling with perturbed weights", "ex2.c", options->perturbed_weights, &options->perturbed_weights, NULL)); 220 PetscCall(PetscOptionsBool("-check_weights", "Ensure all particle weights are positive", "ex2.c", options->checkweights, &options->checkweights, NULL)); 221 PetscCall(PetscOptionsInt("-check_vel_res", "Check particle velocity residuals for nightly tests", "ex2.c", options->checkVRes, &options->checkVRes, NULL)); 222 PetscCall(PetscOptionsInt("-output_step", "Number of time steps between output", "ex2.c", options->ostep, &options->ostep, NULL)); 223 PetscCall(PetscOptionsReal("-timeScale", "Nondimensionalizing time scale", "ex2.c", options->timeScale, &options->timeScale, NULL)); 224 PetscCall(PetscOptionsReal("-initial_velocity", "Initial velocity of perturbed particle", "ex2.c", options->initVel, &options->initVel, NULL)); 225 PetscCall(PetscOptionsReal("-total_weight", "Total weight of all particles", "ex2.c", options->totalWeight, &options->totalWeight, NULL)); 226 PetscCall(PetscOptionsRealArray("-cosine_coefficients", "Amplitude and frequency of cosine equation used in initialization", "ex2.c", options->cosine_coefficients, &d, NULL)); 227 PetscCall(PetscOptionsRealArray("-charges", "Species charges", "ex2.c", options->charges, &maxSpecies, NULL)); 228 PetscCall(PetscOptionsEnum("-em_type", "Type of electrostatic solver", "ex2.c", EMTypes, (PetscEnum)options->em, (PetscEnum *)&options->em, NULL)); 229 PetscOptionsEnd(); 230 231 PetscCall(PetscLogEventRegister("RhsX", TS_CLASSID, &options->RhsXEvent)); 232 PetscCall(PetscLogEventRegister("RhsV", TS_CLASSID, &options->RhsVEvent)); 233 PetscCall(PetscLogEventRegister("ESolve", TS_CLASSID, &options->ESolveEvent)); 234 PetscCall(PetscLogEventRegister("ETab", TS_CLASSID, &options->ETabEvent)); 235 PetscFunctionReturn(PETSC_SUCCESS); 236 } 237 238 static PetscErrorCode SetupContext(DM dm, DM sw, AppCtx *user) 239 { 240 PetscFunctionBeginUser; 241 if (user->efield_monitor) { 242 PetscDraw draw; 243 PetscDrawAxis axis; 244 245 PetscCall(PetscDrawCreate(PETSC_COMM_WORLD, NULL, "Max Electric Field", 0, 300, 400, 300, &draw)); 246 PetscCall(PetscDrawSetSave(draw, "ex9_Efield")); 247 PetscCall(PetscDrawSetFromOptions(draw)); 248 PetscCall(PetscDrawLGCreate(draw, 1, &user->drawlgE)); 249 PetscCall(PetscDrawDestroy(&draw)); 250 PetscCall(PetscDrawLGGetAxis(user->drawlgE, &axis)); 251 PetscCall(PetscDrawAxisSetLabels(axis, "Electron Electric Field", "time", "E_max")); 252 PetscCall(PetscDrawLGSetLimits(user->drawlgE, 0., user->steps * user->stepSize, user->drawlgEmin, 0.)); 253 } 254 if (user->initial_monitor) { 255 PetscDrawAxis axis1, axis2, axis3; 256 PetscReal dmboxlower[2], dmboxupper[2]; 257 PetscInt dim, cStart, cEnd; 258 PetscCall(DMGetDimension(sw, &dim)); 259 PetscCall(DMGetBoundingBox(dm, dmboxlower, dmboxupper)); 260 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 261 262 PetscCall(PetscDrawCreate(PETSC_COMM_WORLD, NULL, "monitor_initial_conditions_x", 0, 300, 400, 300, &user->drawic_x)); 263 PetscCall(PetscDrawSetSave(user->drawic_x, "ex9_ic_x")); 264 PetscCall(PetscDrawSetFromOptions(user->drawic_x)); 265 PetscCall(PetscDrawHGCreate(user->drawic_x, (int)dim, &user->drawhgic_x)); 266 PetscCall(PetscDrawHGGetAxis(user->drawhgic_x, &axis1)); 267 PetscCall(PetscDrawHGSetNumberBins(user->drawhgic_x, (int)(cEnd - cStart))); 268 PetscCall(PetscDrawAxisSetLabels(axis1, "Initial X Distribution", "X", "counts")); 269 PetscCall(PetscDrawAxisSetLimits(axis1, dmboxlower[0], dmboxupper[0], 0, 1500)); 270 271 PetscCall(PetscDrawCreate(PETSC_COMM_WORLD, NULL, "monitor_initial_conditions_v", 400, 300, 400, 300, &user->drawic_v)); 272 PetscCall(PetscDrawSetSave(user->drawic_v, "ex9_ic_v")); 273 PetscCall(PetscDrawSetFromOptions(user->drawic_v)); 274 PetscCall(PetscDrawHGCreate(user->drawic_v, (int)dim, &user->drawhgic_v)); 275 PetscCall(PetscDrawHGGetAxis(user->drawhgic_v, &axis2)); 276 PetscCall(PetscDrawHGSetNumberBins(user->drawhgic_v, 1000)); 277 PetscCall(PetscDrawAxisSetLabels(axis2, "Initial V_x Distribution", "V", "counts")); 278 PetscCall(PetscDrawAxisSetLimits(axis2, -1, 1, 0, 1500)); 279 280 PetscCall(PetscDrawCreate(PETSC_COMM_WORLD, NULL, "monitor_initial_conditions_w", 800, 300, 400, 300, &user->drawic_w)); 281 PetscCall(PetscDrawSetSave(user->drawic_w, "ex9_ic_w")); 282 PetscCall(PetscDrawSetFromOptions(user->drawic_w)); 283 PetscCall(PetscDrawHGCreate(user->drawic_w, (int)dim, &user->drawhgic_w)); 284 PetscCall(PetscDrawHGGetAxis(user->drawhgic_w, &axis3)); 285 PetscCall(PetscDrawHGSetNumberBins(user->drawhgic_w, 10)); 286 PetscCall(PetscDrawAxisSetLabels(axis3, "Initial W Distribution", "weight", "counts")); 287 PetscCall(PetscDrawAxisSetLimits(axis3, 0, 0.01, 0, 5000)); 288 } 289 if (user->positions_monitor) { 290 PetscDraw draw; 291 PetscDrawAxis axis; 292 293 PetscCall(PetscDrawCreate(PETSC_COMM_WORLD, NULL, "Particle Position", 0, 0, 400, 300, &draw)); 294 PetscCall(PetscDrawSetSave(draw, "ex9_pos")); 295 PetscCall(PetscDrawSetFromOptions(draw)); 296 PetscCall(PetscDrawSPCreate(draw, 10, &user->drawspX)); 297 PetscCall(PetscDrawDestroy(&draw)); 298 PetscCall(PetscDrawSPSetDimension(user->drawspX, 1)); 299 PetscCall(PetscDrawSPGetAxis(user->drawspX, &axis)); 300 PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "x", "v")); 301 PetscCall(PetscDrawSPReset(user->drawspX)); 302 } 303 if (user->poisson_monitor) { 304 Vec rho, phi; 305 PetscDraw draw; 306 PetscDrawAxis axis; 307 308 PetscCall(PetscDrawCreate(PETSC_COMM_WORLD, NULL, "Electric_Field", 0, 0, 400, 300, &draw)); 309 PetscCall(PetscDrawSetFromOptions(draw)); 310 PetscCall(PetscDrawSetSave(draw, "ex9_E_spatial")); 311 PetscCall(PetscDrawSPCreate(draw, 10, &user->drawspE)); 312 PetscCall(PetscDrawDestroy(&draw)); 313 PetscCall(PetscDrawSPSetDimension(user->drawspE, 1)); 314 PetscCall(PetscDrawSPGetAxis(user->drawspE, &axis)); 315 PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "x", "E")); 316 PetscCall(PetscDrawSPReset(user->drawspE)); 317 318 PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD, NULL, "Charge Density", 0, 0, 400, 300, &user->viewerRho)); 319 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)user->viewerRho, "rho_")); 320 PetscCall(PetscViewerDrawGetDraw(user->viewerRho, 0, &draw)); 321 PetscCall(PetscDrawSetSave(draw, "ex9_rho_spatial")); 322 PetscCall(PetscViewerSetFromOptions(user->viewerRho)); 323 PetscCall(DMGetNamedGlobalVector(dm, "rho", &rho)); 324 PetscCall(PetscObjectSetName((PetscObject)rho, "charge_density")); 325 PetscCall(DMRestoreNamedGlobalVector(dm, "rho", &rho)); 326 327 PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD, NULL, "Potential", 400, 0, 400, 300, &user->viewerPhi)); 328 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)user->viewerPhi, "phi_")); 329 PetscCall(PetscViewerDrawGetDraw(user->viewerPhi, 0, &draw)); 330 PetscCall(PetscDrawSetSave(draw, "ex9_phi_spatial")); 331 PetscCall(PetscViewerSetFromOptions(user->viewerPhi)); 332 PetscCall(DMGetNamedGlobalVector(dm, "phi", &phi)); 333 PetscCall(PetscObjectSetName((PetscObject)phi, "potential")); 334 PetscCall(DMRestoreNamedGlobalVector(dm, "phi", &phi)); 335 } 336 PetscFunctionReturn(PETSC_SUCCESS); 337 } 338 339 static PetscErrorCode DestroyContext(AppCtx *user) 340 { 341 PetscFunctionBeginUser; 342 PetscCall(PetscDrawHGDestroy(&user->drawhgic_x)); 343 PetscCall(PetscDrawDestroy(&user->drawic_x)); 344 PetscCall(PetscDrawHGDestroy(&user->drawhgic_v)); 345 PetscCall(PetscDrawDestroy(&user->drawic_v)); 346 PetscCall(PetscDrawHGDestroy(&user->drawhgic_w)); 347 PetscCall(PetscDrawDestroy(&user->drawic_w)); 348 349 PetscCall(PetscDrawLGDestroy(&user->drawlgE)); 350 PetscCall(PetscDrawSPDestroy(&user->drawspE)); 351 PetscCall(PetscDrawSPDestroy(&user->drawspX)); 352 PetscCall(PetscViewerDestroy(&user->viewerRho)); 353 PetscCall(PetscViewerDestroy(&user->viewerPhi)); 354 355 PetscCall(PetscBagDestroy(&user->bag)); 356 PetscFunctionReturn(PETSC_SUCCESS); 357 } 358 359 static PetscErrorCode CheckNonNegativeWeights(DM sw, AppCtx *user) 360 { 361 const PetscScalar *w; 362 PetscInt Np; 363 364 PetscFunctionBeginUser; 365 if (!user->checkweights) PetscFunctionReturn(PETSC_SUCCESS); 366 PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&w)); 367 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 368 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]); 369 PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&w)); 370 PetscFunctionReturn(PETSC_SUCCESS); 371 } 372 373 static PetscErrorCode computeVelocityFEMMoments(DM sw, PetscReal moments[], AppCtx *user) 374 { 375 DMSwarmCellDM celldm; 376 DM vdm; 377 Vec u[1]; 378 const char *fields[1] = {"w_q"}; 379 380 PetscFunctionBegin; 381 PetscCall(DMSwarmSetCellDMActive(sw, "velocity")); 382 PetscCall(DMSwarmGetCellDMActive(sw, &celldm)); 383 PetscCall(DMSwarmCellDMGetDM(celldm, &vdm)); 384 PetscCall(DMGetGlobalVector(vdm, &u[0])); 385 PetscCall(DMSwarmProjectFields(sw, vdm, 1, fields, u, SCATTER_FORWARD)); 386 PetscCall(DMPlexComputeMoments(vdm, u[0], moments)); 387 PetscCall(DMRestoreGlobalVector(vdm, &u[0])); 388 PetscCall(DMSwarmSetCellDMActive(sw, "space")); 389 PetscFunctionReturn(PETSC_SUCCESS); 390 } 391 392 static PetscErrorCode MonitorEField(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) 393 { 394 AppCtx *user = (AppCtx *)ctx; 395 DM sw; 396 PetscReal *E, *x, *weight; 397 PetscReal Enorm = 0., lgEnorm, lgEmax, sum = 0., Emax = 0., chargesum = 0.; 398 PetscReal pmoments[4]; /* \int f, \int v f, \int v^2 f */ 399 PetscInt *species, dim, Np; 400 401 PetscFunctionBeginUser; 402 if (step < 0 || !user->validE) PetscFunctionReturn(PETSC_SUCCESS); 403 PetscCall(TSGetDM(ts, &sw)); 404 PetscCall(DMGetDimension(sw, &dim)); 405 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 406 PetscCall(DMSwarmSortGetAccess(sw)); 407 PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x)); 408 PetscCall(DMSwarmGetField(sw, "E_field", NULL, NULL, (void **)&E)); 409 PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 410 PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&weight)); 411 412 for (PetscInt p = 0; p < Np; ++p) { 413 for (PetscInt d = 0; d < 1; ++d) { 414 PetscReal temp = PetscAbsReal(E[p * dim + d]); 415 if (temp > Emax) Emax = temp; 416 } 417 Enorm += PetscSqrtReal(E[p * dim] * E[p * dim]); 418 sum += E[p * dim]; 419 chargesum += user->charges[0] * weight[p]; 420 } 421 lgEnorm = Enorm != 0 ? PetscLog10Real(Enorm) : -16.; 422 lgEmax = Emax != 0 ? PetscLog10Real(Emax) : user->drawlgEmin; 423 424 PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x)); 425 PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight)); 426 PetscCall(DMSwarmRestoreField(sw, "E_field", NULL, NULL, (void **)&E)); 427 PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 428 PetscCall(PetscDrawLGAddPoint(user->drawlgE, &t, &lgEmax)); 429 PetscCall(PetscDrawLGDraw(user->drawlgE)); 430 PetscDraw draw; 431 PetscCall(PetscDrawLGGetDraw(user->drawlgE, &draw)); 432 PetscCall(PetscDrawSave(draw)); 433 434 PetscCall(DMSwarmComputeMoments(sw, "velocity", "w_q", pmoments)); 435 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%f\t%+e\t%e\t%f\t%f\t%f\t%f\t%f\t%f\t%f\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])); 436 PetscCall(DMViewFromOptions(sw, NULL, "-sw_efield_view")); 437 PetscFunctionReturn(PETSC_SUCCESS); 438 } 439 440 static PetscErrorCode MonitorMoments(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) 441 { 442 AppCtx *user = (AppCtx *)ctx; 443 DM sw; 444 PetscReal pmoments[4], fmoments[4]; /* \int f, \int v f, \int v^2 f */ 445 446 PetscFunctionBeginUser; 447 if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); 448 PetscCall(TSGetDM(ts, &sw)); 449 450 PetscCall(DMSwarmComputeMoments(sw, "velocity", "w_q", pmoments)); 451 PetscCall(computeVelocityFEMMoments(sw, fmoments, user)); 452 453 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%f\t%f\t%f\t%f\t%f\t%f\t%f\n", (double)t, (double)pmoments[0], (double)pmoments[1], (double)pmoments[3], (double)fmoments[0], (double)fmoments[1], (double)fmoments[2])); 454 PetscFunctionReturn(PETSC_SUCCESS); 455 } 456 457 PetscErrorCode MonitorInitialConditions(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) 458 { 459 AppCtx *user = (AppCtx *)ctx; 460 DM dm, sw; 461 const PetscScalar *u; 462 PetscReal *weight, *pos, *vel; 463 PetscInt dim, p, Np, cStart, cEnd; 464 465 PetscFunctionBegin; 466 if (step < 0) PetscFunctionReturn(PETSC_SUCCESS); /* -1 indicates interpolated solution */ 467 PetscCall(TSGetDM(ts, &sw)); 468 PetscCall(DMSwarmGetCellDM(sw, &dm)); 469 PetscCall(DMGetDimension(sw, &dim)); 470 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 471 PetscCall(DMSwarmSortGetAccess(sw)); 472 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 473 474 if (step == 0) { 475 PetscCall(PetscDrawHGReset(user->drawhgic_x)); 476 PetscCall(PetscDrawHGGetDraw(user->drawhgic_x, &user->drawic_x)); 477 PetscCall(PetscDrawClear(user->drawic_x)); 478 PetscCall(PetscDrawFlush(user->drawic_x)); 479 480 PetscCall(PetscDrawHGReset(user->drawhgic_v)); 481 PetscCall(PetscDrawHGGetDraw(user->drawhgic_v, &user->drawic_v)); 482 PetscCall(PetscDrawClear(user->drawic_v)); 483 PetscCall(PetscDrawFlush(user->drawic_v)); 484 485 PetscCall(PetscDrawHGReset(user->drawhgic_w)); 486 PetscCall(PetscDrawHGGetDraw(user->drawhgic_w, &user->drawic_w)); 487 PetscCall(PetscDrawClear(user->drawic_w)); 488 PetscCall(PetscDrawFlush(user->drawic_w)); 489 490 PetscCall(VecGetArrayRead(U, &u)); 491 PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&vel)); 492 PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&weight)); 493 PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&pos)); 494 495 PetscCall(VecGetLocalSize(U, &Np)); 496 Np /= dim * 2; 497 for (p = 0; p < Np; ++p) { 498 PetscCall(PetscDrawHGAddValue(user->drawhgic_x, pos[p * dim])); 499 PetscCall(PetscDrawHGAddValue(user->drawhgic_v, vel[p * dim])); 500 PetscCall(PetscDrawHGAddValue(user->drawhgic_w, weight[p])); 501 } 502 503 PetscCall(VecRestoreArrayRead(U, &u)); 504 PetscCall(PetscDrawHGDraw(user->drawhgic_x)); 505 PetscCall(PetscDrawHGSave(user->drawhgic_x)); 506 507 PetscCall(PetscDrawHGDraw(user->drawhgic_v)); 508 PetscCall(PetscDrawHGSave(user->drawhgic_v)); 509 510 PetscCall(PetscDrawHGDraw(user->drawhgic_w)); 511 PetscCall(PetscDrawHGSave(user->drawhgic_w)); 512 513 PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&pos)); 514 PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&vel)); 515 PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight)); 516 } 517 PetscFunctionReturn(PETSC_SUCCESS); 518 } 519 520 static PetscErrorCode MonitorPositions_2D(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) 521 { 522 AppCtx *user = (AppCtx *)ctx; 523 DM dm, sw; 524 PetscScalar *x, *v, *weight; 525 PetscReal lower[3], upper[3], speed; 526 const PetscInt *s; 527 PetscInt dim, cStart, cEnd, c; 528 529 PetscFunctionBeginUser; 530 if (step > 0 && step % user->ostep == 0) { 531 PetscCall(TSGetDM(ts, &sw)); 532 PetscCall(DMSwarmGetCellDM(sw, &dm)); 533 PetscCall(DMGetDimension(dm, &dim)); 534 PetscCall(DMGetBoundingBox(dm, lower, upper)); 535 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 536 PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x)); 537 PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&v)); 538 PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&weight)); 539 PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&s)); 540 PetscCall(DMSwarmSortGetAccess(sw)); 541 PetscCall(PetscDrawSPReset(user->drawspX)); 542 PetscCall(PetscDrawSPSetLimits(user->drawspX, lower[0], upper[0], lower[1], upper[1])); 543 PetscCall(PetscDrawSPSetLimits(user->drawspX, lower[0], upper[0], -12, 12)); 544 for (c = 0; c < cEnd - cStart; ++c) { 545 PetscInt *pidx, Npc, q; 546 PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Npc, &pidx)); 547 for (q = 0; q < Npc; ++q) { 548 const PetscInt p = pidx[q]; 549 if (s[p] == 0) { 550 speed = 0.; 551 for (PetscInt d = 0; d < dim; ++d) speed += PetscSqr(v[p * dim + d]); 552 speed = PetscSqrtReal(speed); 553 if (dim == 1 || user->fake_1D) { 554 PetscCall(PetscDrawSPAddPointColorized(user->drawspX, &x[p * dim], &v[p * dim], &speed)); 555 } else { 556 PetscCall(PetscDrawSPAddPointColorized(user->drawspX, &x[p * dim], &x[p * dim + 1], &speed)); 557 } 558 } else if (s[p] == 1) { 559 PetscCall(PetscDrawSPAddPoint(user->drawspX, &x[p * dim], &v[p * dim])); 560 } 561 } 562 PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Npc, &pidx)); 563 } 564 PetscCall(PetscDrawSPDraw(user->drawspX, PETSC_TRUE)); 565 PetscDraw draw; 566 PetscCall(PetscDrawSPGetDraw(user->drawspX, &draw)); 567 PetscCall(PetscDrawSave(draw)); 568 PetscCall(DMSwarmSortRestoreAccess(sw)); 569 PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x)); 570 PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight)); 571 PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&v)); 572 PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&s)); 573 } 574 PetscFunctionReturn(PETSC_SUCCESS); 575 } 576 577 static PetscErrorCode MonitorPoisson(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) 578 { 579 AppCtx *user = (AppCtx *)ctx; 580 DM dm, sw; 581 582 PetscFunctionBeginUser; 583 if (step > 0 && step % user->ostep == 0) { 584 PetscCall(TSGetDM(ts, &sw)); 585 PetscCall(DMSwarmGetCellDM(sw, &dm)); 586 587 if (user->validE) { 588 PetscScalar *x, *E, *weight; 589 PetscReal lower[3], upper[3], xval; 590 PetscDraw draw; 591 PetscInt dim, cStart, cEnd; 592 593 PetscCall(DMGetDimension(dm, &dim)); 594 PetscCall(DMGetBoundingBox(dm, lower, upper)); 595 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 596 597 PetscCall(PetscDrawSPReset(user->drawspE)); 598 PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x)); 599 PetscCall(DMSwarmGetField(sw, "E_field", NULL, NULL, (void **)&E)); 600 PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&weight)); 601 602 PetscCall(DMSwarmSortGetAccess(sw)); 603 for (PetscInt c = 0; c < cEnd - cStart; ++c) { 604 PetscReal Eavg = 0.0; 605 PetscInt *pidx, Npc; 606 607 PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Npc, &pidx)); 608 for (PetscInt q = 0; q < Npc; ++q) { 609 const PetscInt p = pidx[q]; 610 Eavg += E[p * dim]; 611 } 612 Eavg /= Npc; 613 xval = (c + 0.5) * ((upper[0] - lower[0]) / (cEnd - cStart)); 614 PetscCall(PetscDrawSPAddPoint(user->drawspE, &xval, &Eavg)); 615 PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Npc, &pidx)); 616 } 617 PetscCall(PetscDrawSPDraw(user->drawspE, PETSC_TRUE)); 618 PetscCall(PetscDrawSPGetDraw(user->drawspE, &draw)); 619 PetscCall(PetscDrawSave(draw)); 620 PetscCall(DMSwarmSortRestoreAccess(sw)); 621 PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x)); 622 PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight)); 623 PetscCall(DMSwarmRestoreField(sw, "E_field", NULL, NULL, (void **)&E)); 624 } 625 626 Vec rho, phi; 627 628 PetscCall(DMGetNamedGlobalVector(dm, "rho", &rho)); 629 PetscCall(VecView(rho, user->viewerRho)); 630 PetscCall(DMRestoreNamedGlobalVector(dm, "rho", &rho)); 631 632 PetscCall(DMGetNamedGlobalVector(dm, "phi", &phi)); 633 PetscCall(VecView(phi, user->viewerPhi)); 634 PetscCall(DMRestoreNamedGlobalVector(dm, "phi", &phi)); 635 } 636 PetscFunctionReturn(PETSC_SUCCESS); 637 } 638 639 static PetscErrorCode MonitorSwarmWeights(DM sw, const char field[]) 640 { 641 PetscReal *w; 642 PetscInt bs; 643 644 PetscFunctionBegin; 645 PetscCall(DMSwarmGetField(sw, field, &bs, NULL, (void **)&w)); 646 PetscCall(DMSwarmRestoreField(sw, field, &bs, NULL, (void **)&w)); 647 PetscFunctionReturn(PETSC_SUCCESS); 648 } 649 650 static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *ctx) 651 { 652 PetscBag bag; 653 Parameter *p; 654 655 PetscFunctionBeginUser; 656 /* setup PETSc parameter bag */ 657 PetscCall(PetscBagGetData(ctx->bag, (void **)&p)); 658 PetscCall(PetscBagSetName(ctx->bag, "par", "Vlasov-Poisson Parameters")); 659 bag = ctx->bag; 660 PetscCall(PetscBagRegisterScalar(bag, &p->v0, 1.0, "v0", "Velocity scale, m/s")); 661 PetscCall(PetscBagRegisterScalar(bag, &p->t0, 1.0, "t0", "Time scale, s")); 662 PetscCall(PetscBagRegisterScalar(bag, &p->x0, 1.0, "x0", "Space scale, m")); 663 PetscCall(PetscBagRegisterScalar(bag, &p->v0, 1.0, "phi0", "Potential scale, kg*m^2/A*s^3")); 664 PetscCall(PetscBagRegisterScalar(bag, &p->q0, 1.0, "q0", "Charge Scale, A*s")); 665 PetscCall(PetscBagRegisterScalar(bag, &p->m0, 1.0, "m0", "Mass Scale, kg")); 666 PetscCall(PetscBagRegisterScalar(bag, &p->epsi0, 1.0, "epsi0", "Permittivity of Free Space, kg")); 667 PetscCall(PetscBagRegisterScalar(bag, &p->kb, 1.0, "kb", "Boltzmann Constant, m^2 kg/s^2 K^1")); 668 669 PetscCall(PetscBagRegisterScalar(bag, &p->sigma, 1.0, "sigma", "Charge per unit area, C/m^3")); 670 PetscCall(PetscBagRegisterScalar(bag, &p->poissonNumber, 1.0, "poissonNumber", "Non-Dimensional Poisson Number")); 671 PetscCall(PetscBagRegisterScalar(bag, &p->vlasovNumber, 1.0, "vlasovNumber", "Non-Dimensional Vlasov Number")); 672 PetscCall(PetscBagSetFromOptions(bag)); 673 { 674 PetscViewer viewer; 675 PetscViewerFormat format; 676 PetscBool flg; 677 678 PetscCall(PetscOptionsCreateViewer(comm, NULL, NULL, "-param_view", &viewer, &format, &flg)); 679 if (flg) { 680 PetscCall(PetscViewerPushFormat(viewer, format)); 681 PetscCall(PetscBagView(bag, viewer)); 682 PetscCall(PetscViewerFlush(viewer)); 683 PetscCall(PetscViewerPopFormat(viewer)); 684 PetscCall(PetscViewerDestroy(&viewer)); 685 } 686 } 687 PetscFunctionReturn(PETSC_SUCCESS); 688 } 689 690 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 691 { 692 PetscFunctionBeginUser; 693 PetscCall(DMCreate(comm, dm)); 694 PetscCall(DMSetType(*dm, DMPLEX)); 695 PetscCall(DMSetFromOptions(*dm)); 696 PetscCall(PetscObjectSetName((PetscObject)*dm, "space")); 697 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 698 699 // Cache the mesh geometry 700 DMField coordField; 701 IS cellIS; 702 PetscQuadrature quad; 703 PetscReal *wt, *pt; 704 PetscInt cdim, cStart, cEnd; 705 706 PetscCall(DMGetCoordinateField(*dm, &coordField)); 707 PetscCheck(coordField, comm, PETSC_ERR_USER, "DM must have a coordinate field"); 708 PetscCall(DMGetCoordinateDim(*dm, &cdim)); 709 PetscCall(DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd)); 710 PetscCall(ISCreateStride(PETSC_COMM_SELF, cEnd - cStart, cStart, 1, &cellIS)); 711 PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &quad)); 712 PetscCall(PetscMalloc1(1, &wt)); 713 PetscCall(PetscMalloc1(2, &pt)); 714 wt[0] = 1.; 715 pt[0] = -1.; 716 pt[1] = -1.; 717 PetscCall(PetscQuadratureSetData(quad, cdim, 1, 1, pt, wt)); 718 PetscCall(DMFieldCreateFEGeom(coordField, cellIS, quad, PETSC_FALSE, &user->fegeom)); 719 PetscCall(PetscQuadratureDestroy(&quad)); 720 PetscCall(ISDestroy(&cellIS)); 721 PetscFunctionReturn(PETSC_SUCCESS); 722 } 723 724 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[]) 725 { 726 f0[0] = -constants[SIGMA]; 727 } 728 729 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[]) 730 { 731 PetscInt d; 732 for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 733 } 734 735 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[]) 736 { 737 PetscInt d; 738 for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0; 739 } 740 741 static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 742 { 743 *u = 0.0; 744 return PETSC_SUCCESS; 745 } 746 747 /* 748 / I -grad\ / q \ = /0\ 749 \-div 0 / \phi/ \f/ 750 */ 751 static void f0_q(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[]) 752 { 753 for (PetscInt d = 0; d < dim; ++d) f0[d] += u[uOff[0] + d]; 754 } 755 756 static void f1_q(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[]) 757 { 758 for (PetscInt d = 0; d < dim; ++d) f1[d * dim + d] = u[uOff[1]]; 759 } 760 761 static void f0_phi_backgroundCharge(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[]) 762 { 763 f0[0] += constants[SIGMA]; 764 for (PetscInt d = 0; d < dim; ++d) f0[0] += u_x[uOff_x[0] + d * dim + d]; 765 } 766 767 /* Boundary residual. Dirichlet boundary for u means u_bdy=p*n */ 768 static void g0_qq(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 g0[]) 769 { 770 for (PetscInt d = 0; d < dim; ++d) g0[d * dim + d] = 1.0; 771 } 772 773 static void g2_qphi(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 g2[]) 774 { 775 for (PetscInt d = 0; d < dim; ++d) g2[d * dim + d] = 1.0; 776 } 777 778 static void g1_phiq(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 g1[]) 779 { 780 for (PetscInt d = 0; d < dim; ++d) g1[d * dim + d] = 1.0; 781 } 782 783 static PetscErrorCode CreateFEM(DM dm, AppCtx *user) 784 { 785 PetscFE fephi, feq; 786 PetscDS ds; 787 PetscBool simplex; 788 PetscInt dim; 789 790 PetscFunctionBeginUser; 791 PetscCall(DMGetDimension(dm, &dim)); 792 PetscCall(DMPlexIsSimplex(dm, &simplex)); 793 if (user->em == EM_MIXED) { 794 DMLabel label; 795 const PetscInt id = 1; 796 797 PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, dim, simplex, "field_", PETSC_DETERMINE, &feq)); 798 PetscCall(PetscObjectSetName((PetscObject)feq, "field")); 799 PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "potential_", PETSC_DETERMINE, &fephi)); 800 PetscCall(PetscObjectSetName((PetscObject)fephi, "potential")); 801 PetscCall(PetscFECopyQuadrature(feq, fephi)); 802 PetscCall(DMSetField(dm, 0, NULL, (PetscObject)feq)); 803 PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fephi)); 804 PetscCall(DMCreateDS(dm)); 805 PetscCall(PetscFEDestroy(&fephi)); 806 PetscCall(PetscFEDestroy(&feq)); 807 808 PetscCall(DMGetLabel(dm, "marker", &label)); 809 PetscCall(DMGetDS(dm, &ds)); 810 811 PetscCall(PetscDSSetResidual(ds, 0, f0_q, f1_q)); 812 PetscCall(PetscDSSetResidual(ds, 1, f0_phi_backgroundCharge, NULL)); 813 PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_qq, NULL, NULL, NULL)); 814 PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_qphi, NULL)); 815 PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_phiq, NULL, NULL)); 816 817 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))zero, NULL, NULL, NULL)); 818 819 } else { 820 MatNullSpace nullsp; 821 PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, PETSC_DETERMINE, &fephi)); 822 PetscCall(PetscObjectSetName((PetscObject)fephi, "potential")); 823 PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fephi)); 824 PetscCall(DMCreateDS(dm)); 825 PetscCall(DMGetDS(dm, &ds)); 826 PetscCall(PetscDSSetResidual(ds, 0, ion_f0, laplacian_f1)); 827 PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, laplacian_g3)); 828 PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_TRUE, 0, NULL, &nullsp)); 829 PetscCall(PetscObjectCompose((PetscObject)fephi, "nullspace", (PetscObject)nullsp)); 830 PetscCall(MatNullSpaceDestroy(&nullsp)); 831 PetscCall(PetscFEDestroy(&fephi)); 832 } 833 PetscFunctionReturn(PETSC_SUCCESS); 834 } 835 836 static PetscErrorCode CreatePoisson(DM dm, AppCtx *user) 837 { 838 SNES snes; 839 Mat J; 840 MatNullSpace nullSpace; 841 842 PetscFunctionBeginUser; 843 PetscCall(CreateFEM(dm, user)); 844 PetscCall(SNESCreate(PetscObjectComm((PetscObject)dm), &snes)); 845 PetscCall(SNESSetOptionsPrefix(snes, "em_")); 846 PetscCall(SNESSetDM(snes, dm)); 847 PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, user)); 848 PetscCall(SNESSetFromOptions(snes)); 849 850 PetscCall(DMCreateMatrix(dm, &J)); 851 PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_TRUE, 0, NULL, &nullSpace)); 852 PetscCall(MatSetNullSpace(J, nullSpace)); 853 PetscCall(MatNullSpaceDestroy(&nullSpace)); 854 PetscCall(SNESSetJacobian(snes, J, J, NULL, NULL)); 855 PetscCall(MatDestroy(&J)); 856 PetscCall(DMCreateMassMatrix(dm, dm, &user->M)); 857 PetscCall(DMPlexCreateClosureIndex(dm, NULL)); 858 user->snes = snes; 859 PetscFunctionReturn(PETSC_SUCCESS); 860 } 861 862 PetscErrorCode PetscPDFPertubedConstant2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 863 { 864 p[0] = (1 + 0.01 * PetscCosReal(0.5 * x[0])) / (2 * PETSC_PI); 865 p[1] = (1 + 0.01 * PetscCosReal(0.5 * x[1])) / (2 * PETSC_PI); 866 return PETSC_SUCCESS; 867 } 868 PetscErrorCode PetscPDFPertubedConstant1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 869 { 870 p[0] = (1. + 0.01 * PetscCosReal(0.5 * x[0])) / (2 * PETSC_PI); 871 return PETSC_SUCCESS; 872 } 873 874 PetscErrorCode PetscPDFCosine1D(const PetscReal x[], const PetscReal scale[], PetscReal p[]) 875 { 876 const PetscReal alpha = scale ? scale[0] : 0.0; 877 const PetscReal k = scale ? scale[1] : 1.; 878 p[0] = (1 + alpha * PetscCosReal(k * x[0])); 879 return PETSC_SUCCESS; 880 } 881 882 PetscErrorCode PetscPDFCosine2D(const PetscReal x[], const PetscReal scale[], PetscReal p[]) 883 { 884 const PetscReal alpha = scale ? scale[0] : 0.; 885 const PetscReal k = scale ? scale[0] : 1.; 886 p[0] = (1 + alpha * PetscCosReal(k * (x[0] + x[1]))); 887 return PETSC_SUCCESS; 888 } 889 890 static PetscErrorCode CreateVelocityDM(DM sw, DM *vdm) 891 { 892 PetscFE fe; 893 DMPolytopeType ct; 894 PetscInt dim, cStart; 895 const char *prefix = "v"; 896 897 PetscFunctionBegin; 898 PetscCall(DMCreate(PETSC_COMM_SELF, vdm)); 899 PetscCall(DMSetType(*vdm, DMPLEX)); 900 PetscCall(DMPlexSetOptionsPrefix(*vdm, prefix)); 901 PetscCall(DMSetFromOptions(*vdm)); 902 PetscCall(PetscObjectSetName((PetscObject)*vdm, "velocity")); 903 PetscCall(DMViewFromOptions(*vdm, NULL, "-dm_view")); 904 905 PetscCall(DMGetDimension(*vdm, &dim)); 906 PetscCall(DMPlexGetHeightStratum(*vdm, 0, &cStart, NULL)); 907 PetscCall(DMPlexGetCellType(*vdm, cStart, &ct)); 908 PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, 1, ct, prefix, PETSC_DETERMINE, &fe)); 909 PetscCall(PetscObjectSetName((PetscObject)fe, "distribution")); 910 PetscCall(DMSetField(*vdm, 0, NULL, (PetscObject)fe)); 911 PetscCall(DMCreateDS(*vdm)); 912 PetscCall(PetscFEDestroy(&fe)); 913 PetscFunctionReturn(PETSC_SUCCESS); 914 } 915 916 static PetscErrorCode CreateRemapDM(DM sw, DM *rdm) 917 { 918 PetscFE fe; 919 DMPolytopeType ct; 920 PetscInt dim, cStart; 921 const char *prefix = "remap_"; 922 923 PetscFunctionBegin; 924 PetscCall(DMCreate(PetscObjectComm((PetscObject)sw), rdm)); 925 PetscCall(DMSetType(*rdm, DMPLEX)); 926 PetscCall(DMPlexSetOptionsPrefix(*rdm, prefix)); 927 PetscCall(DMSetFromOptions(*rdm)); 928 PetscCall(PetscObjectSetName((PetscObject)*rdm, "remap")); 929 PetscCall(DMViewFromOptions(*rdm, NULL, "-dm_view")); 930 931 PetscCall(DMGetDimension(*rdm, &dim)); 932 PetscCall(DMPlexGetHeightStratum(*rdm, 0, &cStart, NULL)); 933 PetscCall(DMPlexGetCellType(*rdm, cStart, &ct)); 934 PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, 1, ct, prefix, PETSC_DETERMINE, &fe)); 935 PetscCall(PetscObjectSetName((PetscObject)fe, "distribution")); 936 PetscCall(DMSetField(*rdm, 0, NULL, (PetscObject)fe)); 937 PetscCall(DMCreateDS(*rdm)); 938 PetscCall(PetscFEDestroy(&fe)); 939 PetscFunctionReturn(PETSC_SUCCESS); 940 } 941 942 /* 943 InitializeParticles_Regular - Initialize a regular grid of particles in each cell 944 945 Input Parameters: 946 + sw - The `DMSWARM` 947 - n - The number of particles per dimension per species 948 949 Notes: 950 This functions sets the species, cellid, and cell DM coordinates. 951 952 It places n^d particles per species in each cell of the cell DM. 953 */ 954 static PetscErrorCode InitializeParticles_Regular(DM sw, PetscInt n) 955 { 956 DM_Swarm *swarm = (DM_Swarm *)sw->data; 957 DM dm; 958 DMSwarmCellDM celldm; 959 PetscInt dim, Ns, Npc, Np, cStart, cEnd, debug; 960 PetscBool flg; 961 MPI_Comm comm; 962 963 PetscFunctionBegin; 964 PetscCall(PetscObjectGetComm((PetscObject)sw, &comm)); 965 966 PetscOptionsBegin(comm, "", "DMSwarm Options", "DMSWARM"); 967 PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 968 PetscCall(PetscOptionsInt("-dm_swarm_num_species", "The number of species", "DMSwarmSetNumSpecies", Ns, &Ns, &flg)); 969 if (flg) PetscCall(DMSwarmSetNumSpecies(sw, Ns)); 970 PetscCall(PetscOptionsBoundedInt("-dm_swarm_print_coords", "Debug output level for particle coordinate computations", "InitializeParticles", 0, &swarm->printCoords, NULL, 0)); 971 PetscCall(PetscOptionsBoundedInt("-dm_swarm_print_weights", "Debug output level for particle weight computations", "InitializeWeights", 0, &swarm->printWeights, NULL, 0)); 972 PetscOptionsEnd(); 973 debug = swarm->printCoords; 974 975 // n^d particle per cell on the grid 976 PetscCall(DMSwarmGetCellDM(sw, &dm)); 977 PetscCall(DMGetDimension(dm, &dim)); 978 PetscCheck(!(dim % 2), comm, PETSC_ERR_SUP, "We only support even dimension, not %" PetscInt_FMT, dim); 979 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 980 Npc = Ns * PetscPowInt(n, dim); 981 Np = (cEnd - cStart) * Npc; 982 PetscCall(DMSwarmSetLocalSizes(sw, Np, 0)); 983 if (debug) { 984 PetscInt gNp; 985 PetscCallMPI(MPIU_Allreduce(&Np, &gNp, 1, MPIU_INT, MPIU_SUM, comm)); 986 PetscCall(PetscPrintf(comm, "Global Np = %" PetscInt_FMT "\n", gNp)); 987 } 988 PetscCall(PetscPrintf(comm, "Regular layout using %" PetscInt_FMT " particles per cell\n", Npc)); 989 990 // Set species and cellid 991 { 992 const char *cellidName; 993 PetscInt *species, *cellid; 994 995 PetscCall(DMSwarmGetCellDMActive(sw, &celldm)); 996 PetscCall(DMSwarmCellDMGetCellID(celldm, &cellidName)); 997 PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 998 PetscCall(DMSwarmGetField(sw, cellidName, NULL, NULL, (void **)&cellid)); 999 for (PetscInt c = 0, p = 0; c < cEnd - cStart; ++c) { 1000 for (PetscInt s = 0; s < Ns; ++s) { 1001 for (PetscInt q = 0; q < Npc / Ns; ++q, ++p) { 1002 species[p] = s; 1003 cellid[p] = c; 1004 } 1005 } 1006 } 1007 PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 1008 PetscCall(DMSwarmRestoreField(sw, cellidName, NULL, NULL, (void **)&cellid)); 1009 } 1010 1011 // Set particle coordinates 1012 { 1013 PetscReal *x, *v; 1014 const char **coordNames; 1015 PetscInt Ncoord; 1016 const PetscInt xdim = dim / 2, vdim = dim / 2; 1017 1018 PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Ncoord, &coordNames)); 1019 PetscCheck(Ncoord == 2, comm, PETSC_ERR_SUP, "We only support regular layout for 2 coordinate fields, not %" PetscInt_FMT, Ncoord); 1020 PetscCall(DMSwarmGetField(sw, coordNames[0], NULL, NULL, (void **)&x)); 1021 PetscCall(DMSwarmGetField(sw, coordNames[1], NULL, NULL, (void **)&v)); 1022 PetscCall(DMSwarmSortGetAccess(sw)); 1023 PetscCall(DMGetCoordinatesLocalSetUp(dm)); 1024 for (PetscInt c = 0; c < cEnd - cStart; ++c) { 1025 const PetscInt cell = c + cStart; 1026 const PetscScalar *a; 1027 PetscScalar *coords; 1028 PetscReal lower[6], upper[6]; 1029 PetscBool isDG; 1030 PetscInt *pidx, npc, Nc; 1031 1032 PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &npc, &pidx)); 1033 PetscCheck(Npc == npc, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of points per cell %" PetscInt_FMT " != %" PetscInt_FMT, npc, Npc); 1034 PetscCall(DMPlexGetCellCoordinates(dm, cell, &isDG, &Nc, &a, &coords)); 1035 for (PetscInt d = 0; d < dim; ++d) { 1036 lower[d] = coords[0 * dim + d]; 1037 upper[d] = coords[0 * dim + d]; 1038 } 1039 for (PetscInt i = 1; i < Nc / dim; ++i) { 1040 for (PetscInt d = 0; d < dim; ++d) { 1041 lower[d] = PetscMin(lower[d], coords[i * dim + d]); 1042 upper[d] = PetscMax(upper[d], coords[i * dim + d]); 1043 } 1044 } 1045 for (PetscInt s = 0; s < Ns; ++s) { 1046 for (PetscInt q = 0; q < Npc / Ns; ++q) { 1047 const PetscInt p = pidx[q * Ns + s]; 1048 PetscInt xi[3], vi[3]; 1049 1050 xi[0] = q % n; 1051 xi[1] = (q / n) % n; 1052 xi[2] = (q / PetscSqr(n)) % n; 1053 for (PetscInt d = 0; d < xdim; ++d) x[p * xdim + d] = lower[d] + (xi[d] + 0.5) * (upper[d] - lower[d]) / n; 1054 vi[0] = (q / PetscPowInt(n, xdim)) % n; 1055 vi[1] = (q / PetscPowInt(n, xdim + 1)) % n; 1056 vi[2] = (q / PetscPowInt(n, xdim + 2)); 1057 for (PetscInt d = 0; d < vdim; ++d) v[p * vdim + d] = lower[xdim + d] + (vi[d] + 0.5) * (upper[xdim + d] - lower[xdim + d]) / n; 1058 if (debug > 1) { 1059 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Particle %4" PetscInt_FMT " ", p)); 1060 PetscCall(PetscPrintf(PETSC_COMM_SELF, " x: (")); 1061 for (PetscInt d = 0; d < xdim; ++d) { 1062 if (d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 1063 PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", x[p * xdim + d])); 1064 } 1065 PetscCall(PetscPrintf(PETSC_COMM_SELF, ") v:(")); 1066 for (PetscInt d = 0; d < vdim; ++d) { 1067 if (d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 1068 PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", v[p * vdim + d])); 1069 } 1070 PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n")); 1071 } 1072 } 1073 } 1074 PetscCall(DMPlexRestoreCellCoordinates(dm, cell, &isDG, &Nc, &a, &coords)); 1075 PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Npc, &pidx)); 1076 } 1077 PetscCall(DMSwarmSortRestoreAccess(sw)); 1078 PetscCall(DMSwarmRestoreField(sw, coordNames[0], NULL, NULL, (void **)&x)); 1079 PetscCall(DMSwarmRestoreField(sw, coordNames[1], NULL, NULL, (void **)&v)); 1080 } 1081 PetscFunctionReturn(PETSC_SUCCESS); 1082 } 1083 1084 /* 1085 InitializeParticles_Centroid - Initialize a regular grid of particles. 1086 1087 Input Parameters: 1088 + sw - The `DMSWARM` 1089 - force1D - Treat the spatial domain as 1D 1090 1091 Notes: 1092 This functions sets the species, cellid, spatial coordinate, and velocity fields for all particles. 1093 1094 It places one particle in the centroid of each cell in the implicit tensor product of the spatial 1095 and velocity meshes. 1096 */ 1097 static PetscErrorCode InitializeParticles_Centroid(DM sw, PetscBool force1D) 1098 { 1099 DM_Swarm *swarm = (DM_Swarm *)sw->data; 1100 DMSwarmCellDM celldm; 1101 DM xdm, vdm; 1102 PetscReal vmin[3], vmax[3]; 1103 PetscReal *x, *v; 1104 PetscInt *species, *cellid; 1105 PetscInt dim, xcStart, xcEnd, vcStart, vcEnd, Ns, Np, Npc, debug; 1106 PetscBool flg; 1107 MPI_Comm comm; 1108 const char *cellidname; 1109 1110 PetscFunctionBegin; 1111 PetscCall(PetscObjectGetComm((PetscObject)sw, &comm)); 1112 1113 PetscOptionsBegin(comm, "", "DMSwarm Options", "DMSWARM"); 1114 PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 1115 PetscCall(PetscOptionsInt("-dm_swarm_num_species", "The number of species", "DMSwarmSetNumSpecies", Ns, &Ns, &flg)); 1116 if (flg) PetscCall(DMSwarmSetNumSpecies(sw, Ns)); 1117 PetscCall(PetscOptionsBoundedInt("-dm_swarm_print_coords", "Debug output level for particle coordinate computations", "InitializeParticles", 0, &swarm->printCoords, NULL, 0)); 1118 PetscCall(PetscOptionsBoundedInt("-dm_swarm_print_weights", "Debug output level for particle weight computations", "InitializeWeights", 0, &swarm->printWeights, NULL, 0)); 1119 PetscOptionsEnd(); 1120 debug = swarm->printCoords; 1121 1122 PetscCall(DMGetDimension(sw, &dim)); 1123 PetscCall(DMSwarmGetCellDM(sw, &xdm)); 1124 PetscCall(DMPlexGetHeightStratum(xdm, 0, &xcStart, &xcEnd)); 1125 1126 PetscCall(PetscObjectQuery((PetscObject)sw, "__vdm__", (PetscObject *)&vdm)); 1127 PetscCall(DMPlexGetHeightStratum(vdm, 0, &vcStart, &vcEnd)); 1128 1129 // One particle per centroid on the tensor product grid 1130 Npc = (vcEnd - vcStart) * Ns; 1131 Np = (xcEnd - xcStart) * Npc; 1132 PetscCall(DMSwarmSetLocalSizes(sw, Np, 0)); 1133 if (debug) { 1134 PetscInt gNp; 1135 PetscCallMPI(MPIU_Allreduce(&Np, &gNp, 1, MPIU_INT, MPIU_SUM, comm)); 1136 PetscCall(PetscPrintf(comm, "Global Np = %" PetscInt_FMT "\n", gNp)); 1137 } 1138 1139 // Set species and cellid 1140 PetscCall(DMSwarmGetCellDMActive(sw, &celldm)); 1141 PetscCall(DMSwarmCellDMGetCellID(celldm, &cellidname)); 1142 PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 1143 PetscCall(DMSwarmGetField(sw, cellidname, NULL, NULL, (void **)&cellid)); 1144 for (PetscInt c = 0, p = 0; c < xcEnd - xcStart; ++c) { 1145 for (PetscInt s = 0; s < Ns; ++s) { 1146 for (PetscInt q = 0; q < Npc / Ns; ++q, ++p) { 1147 species[p] = s; 1148 cellid[p] = c; 1149 } 1150 } 1151 } 1152 PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 1153 PetscCall(DMSwarmRestoreField(sw, cellidname, NULL, NULL, (void **)&cellid)); 1154 1155 // Set particle coordinates 1156 PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x)); 1157 PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&v)); 1158 PetscCall(DMSwarmSortGetAccess(sw)); 1159 PetscCall(DMGetBoundingBox(vdm, vmin, vmax)); 1160 PetscCall(DMGetCoordinatesLocalSetUp(xdm)); 1161 for (PetscInt c = 0; c < xcEnd - xcStart; ++c) { 1162 const PetscInt cell = c + xcStart; 1163 PetscInt *pidx, Npc; 1164 PetscReal centroid[3], volume; 1165 1166 PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Npc, &pidx)); 1167 PetscCall(DMPlexComputeCellGeometryFVM(xdm, cell, &volume, centroid, NULL)); 1168 for (PetscInt s = 0; s < Ns; ++s) { 1169 for (PetscInt q = 0; q < Npc / Ns; ++q) { 1170 const PetscInt p = pidx[q * Ns + s]; 1171 1172 for (PetscInt d = 0; d < dim; ++d) { 1173 x[p * dim + d] = centroid[d]; 1174 v[p * dim + d] = vmin[0] + (q + 0.5) * ((vmax[0] - vmin[0]) / (Npc / Ns)); 1175 if (force1D && d > 0) v[p * dim + d] = 0.; 1176 } 1177 if (debug > 1) { 1178 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Particle %4" PetscInt_FMT " ", p)); 1179 PetscCall(PetscPrintf(PETSC_COMM_SELF, " x: (")); 1180 for (PetscInt d = 0; d < dim; ++d) { 1181 if (d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 1182 PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", x[p * dim + d])); 1183 } 1184 PetscCall(PetscPrintf(PETSC_COMM_SELF, ") v:(")); 1185 for (PetscInt d = 0; d < dim; ++d) { 1186 if (d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 1187 PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", v[p * dim + d])); 1188 } 1189 PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n")); 1190 } 1191 } 1192 } 1193 PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Npc, &pidx)); 1194 } 1195 PetscCall(DMSwarmSortRestoreAccess(sw)); 1196 PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x)); 1197 PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&v)); 1198 PetscFunctionReturn(PETSC_SUCCESS); 1199 } 1200 1201 /* 1202 InitializeWeights - Compute weight for each local particle 1203 1204 Input Parameters: 1205 + sw - The `DMSwarm` 1206 . totalWeight - The sum of all particle weights 1207 . force1D - Flag to treat the problem as 1D 1208 . func - The PDF for the particle spatial distribution 1209 - param - The PDF parameters 1210 1211 Notes: 1212 The PDF for velocity is assumed to be a Gaussian 1213 1214 The particle weights are returned in the `w_q` field of `sw`. 1215 */ 1216 static PetscErrorCode InitializeWeights(DM sw, PetscReal totalWeight, PetscBool force1D, PetscProbFunc func, const PetscReal param[]) 1217 { 1218 DM xdm, vdm; 1219 PetscScalar *weight; 1220 PetscQuadrature xquad; 1221 const PetscReal *xq, *xwq; 1222 const PetscInt order = 5; 1223 PetscReal *xqd = NULL, xi0[3]; 1224 PetscReal xwtot = 0., pwtot = 0.; 1225 PetscInt xNq; 1226 PetscInt dim, Ns, xcStart, xcEnd, vcStart, vcEnd, debug = ((DM_Swarm *)sw->data)->printWeights; 1227 MPI_Comm comm; 1228 PetscMPIInt rank; 1229 1230 PetscFunctionBegin; 1231 PetscCall(PetscObjectGetComm((PetscObject)sw, &comm)); 1232 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1233 PetscCall(DMGetDimension(sw, &dim)); 1234 PetscCall(DMSwarmGetCellDM(sw, &xdm)); 1235 PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 1236 PetscCall(DMPlexGetHeightStratum(xdm, 0, &xcStart, &xcEnd)); 1237 PetscCall(PetscObjectQuery((PetscObject)sw, "__vdm__", (PetscObject *)&vdm)); 1238 PetscCall(DMPlexGetHeightStratum(vdm, 0, &vcStart, &vcEnd)); 1239 1240 // Setup Quadrature for spatial and velocity weight calculations 1241 if (force1D) PetscCall(PetscDTGaussTensorQuadrature(1, 1, order, -1.0, 1.0, &xquad)); 1242 else PetscCall(PetscDTGaussTensorQuadrature(dim, 1, order, -1.0, 1.0, &xquad)); 1243 PetscCall(PetscQuadratureGetData(xquad, NULL, NULL, &xNq, &xq, &xwq)); 1244 if (force1D) { 1245 PetscCall(PetscCalloc1(xNq * dim, &xqd)); 1246 for (PetscInt q = 0; q < xNq; ++q) xqd[q * dim] = xq[q]; 1247 } 1248 for (PetscInt d = 0; d < dim; ++d) xi0[d] = -1.0; 1249 1250 // Integrate the density function to get the weights of particles in each cell 1251 PetscCall(DMGetCoordinatesLocalSetUp(vdm)); 1252 PetscCall(DMSwarmSortGetAccess(sw)); 1253 PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&weight)); 1254 for (PetscInt c = xcStart; c < xcEnd; ++c) { 1255 PetscReal xv0[3], xJ[9], xinvJ[9], xdetJ, xqr[3], xden, xw = 0.; 1256 PetscInt *pidx, Npc; 1257 PetscInt xNc; 1258 const PetscScalar *xarray; 1259 PetscScalar *xcoords = NULL; 1260 PetscBool xisDG; 1261 1262 PetscCall(DMPlexGetCellCoordinates(xdm, c, &xisDG, &xNc, &xarray, &xcoords)); 1263 PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Npc, &pidx)); 1264 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); 1265 PetscCall(DMPlexComputeCellGeometryFEM(xdm, c, NULL, xv0, xJ, xinvJ, &xdetJ)); 1266 for (PetscInt q = 0; q < xNq; ++q) { 1267 // Transform quadrature points from ref space to real space 1268 if (force1D) CoordinatesRefToReal(dim, dim, xi0, xv0, xJ, &xqd[q * dim], xqr); 1269 else CoordinatesRefToReal(dim, dim, xi0, xv0, xJ, &xq[q * dim], xqr); 1270 // Get probability density at quad point 1271 // No need to scale xqr since PDF will be periodic 1272 PetscCall((*func)(xqr, param, &xden)); 1273 if (force1D) xdetJ = xJ[0]; // Only want x contribution 1274 xw += xden * (xwq[q] * xdetJ); 1275 } 1276 xwtot += xw; 1277 if (debug) { 1278 IS globalOrdering; 1279 const PetscInt *ordering; 1280 1281 PetscCall(DMPlexGetCellNumbering(xdm, &globalOrdering)); 1282 PetscCall(ISGetIndices(globalOrdering, &ordering)); 1283 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)); 1284 PetscCall(ISRestoreIndices(globalOrdering, &ordering)); 1285 } 1286 // Set weights to be Gaussian in velocity cells 1287 for (PetscInt vc = vcStart; vc < vcEnd; ++vc) { 1288 const PetscInt p = pidx[vc * Ns + 0]; 1289 PetscReal vw = 0.; 1290 PetscInt vNc; 1291 const PetscScalar *varray; 1292 PetscScalar *vcoords = NULL; 1293 PetscBool visDG; 1294 1295 PetscCall(DMPlexGetCellCoordinates(vdm, vc, &visDG, &vNc, &varray, &vcoords)); 1296 // TODO: Fix 2 stream Ask Joe 1297 // Two stream function from 1/2pi v^2 e^(-v^2/2) 1298 // vw = 1. / (PetscSqrtReal(2 * PETSC_PI)) * (((coords_v[0] * PetscExpReal(-PetscSqr(coords_v[0]) / 2.)) - (coords_v[1] * PetscExpReal(-PetscSqr(coords_v[1]) / 2.)))) - 0.5 * PetscErfReal(coords_v[0] / PetscSqrtReal(2.)) + 0.5 * (PetscErfReal(coords_v[1] / PetscSqrtReal(2.))); 1299 vw = 0.5 * (PetscErfReal(vcoords[1] / PetscSqrtReal(2.)) - PetscErfReal(vcoords[0] / PetscSqrtReal(2.))); 1300 1301 weight[p] = totalWeight * vw * xw; 1302 pwtot += weight[p]; 1303 PetscCheck(weight[p] <= 10., PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Particle %" PetscInt_FMT " weight exceeeded 1: %g, %g, %g", p, xw, vw, totalWeight); 1304 PetscCall(DMPlexRestoreCellCoordinates(vdm, vc, &visDG, &vNc, &varray, &vcoords)); 1305 if (debug > 1) PetscCall(PetscPrintf(comm, "particle %" PetscInt_FMT ": %g, vw: %g xw: %g\n", p, weight[p], vw, xw)); 1306 } 1307 PetscCall(DMPlexRestoreCellCoordinates(xdm, c, &xisDG, &xNc, &xarray, &xcoords)); 1308 PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Npc, &pidx)); 1309 } 1310 PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight)); 1311 PetscCall(DMSwarmSortRestoreAccess(sw)); 1312 PetscCall(PetscQuadratureDestroy(&xquad)); 1313 if (force1D) PetscCall(PetscFree(xqd)); 1314 1315 if (debug) { 1316 PetscReal wtot[2] = {pwtot, xwtot}, gwtot[2]; 1317 1318 PetscCall(PetscSynchronizedFlush(comm, NULL)); 1319 PetscCallMPI(MPIU_Allreduce(wtot, gwtot, 2, MPIU_REAL, MPIU_SUM, PETSC_COMM_WORLD)); 1320 PetscCall(PetscPrintf(comm, "particle weight sum = %1.10f cell weight sum = %1.10f\n", (double)gwtot[0], (double)gwtot[1])); 1321 } 1322 PetscFunctionReturn(PETSC_SUCCESS); 1323 } 1324 1325 static PetscErrorCode InitializeParticles_PerturbedWeights(DM sw, AppCtx *user) 1326 { 1327 PetscReal scale[2] = {user->cosine_coefficients[0], user->cosine_coefficients[1]}; 1328 PetscInt dim; 1329 1330 PetscFunctionBegin; 1331 PetscCall(DMGetDimension(sw, &dim)); 1332 PetscCall(InitializeParticles_Centroid(sw, user->fake_1D)); 1333 PetscCall(InitializeWeights(sw, user->totalWeight, user->fake_1D, dim == 1 || user->fake_1D ? PetscPDFCosine1D : PetscPDFCosine2D, scale)); 1334 PetscFunctionReturn(PETSC_SUCCESS); 1335 } 1336 1337 static PetscErrorCode InitializeConstants(DM sw, AppCtx *user) 1338 { 1339 DM dm; 1340 PetscInt *species; 1341 PetscReal *weight, totalCharge = 0., totalWeight = 0., gmin[3], gmax[3], global_charge, global_weight; 1342 PetscInt Np, dim; 1343 1344 PetscFunctionBegin; 1345 PetscCall(DMSwarmGetCellDM(sw, &dm)); 1346 PetscCall(DMGetDimension(sw, &dim)); 1347 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 1348 PetscCall(DMGetBoundingBox(dm, gmin, gmax)); 1349 PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&weight)); 1350 PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 1351 for (PetscInt p = 0; p < Np; ++p) { 1352 totalWeight += weight[p]; 1353 totalCharge += user->charges[species[p]] * weight[p]; 1354 } 1355 PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight)); 1356 PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 1357 { 1358 Parameter *param; 1359 PetscReal Area; 1360 1361 PetscCall(PetscBagGetData(user->bag, (void **)¶m)); 1362 switch (dim) { 1363 case 1: 1364 Area = (gmax[0] - gmin[0]); 1365 break; 1366 case 2: 1367 if (user->fake_1D) { 1368 Area = (gmax[0] - gmin[0]); 1369 } else { 1370 Area = (gmax[0] - gmin[0]) * (gmax[1] - gmin[1]); 1371 } 1372 break; 1373 case 3: 1374 Area = (gmax[0] - gmin[0]) * (gmax[1] - gmin[1]) * (gmax[2] - gmin[2]); 1375 break; 1376 default: 1377 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Dimension %" PetscInt_FMT " not supported", dim); 1378 } 1379 PetscCallMPI(MPIU_Allreduce(&totalWeight, &global_weight, 1, MPIU_REAL, MPIU_SUM, PETSC_COMM_WORLD)); 1380 PetscCallMPI(MPIU_Allreduce(&totalCharge, &global_charge, 1, MPIU_REAL, MPIU_SUM, PETSC_COMM_WORLD)); 1381 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)); 1382 param->sigma = PetscAbsReal(global_charge / (Area)); 1383 1384 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "sigma: %g\n", (double)param->sigma)); 1385 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, 1386 (double)param->vlasovNumber)); 1387 } 1388 /* Setup Constants */ 1389 { 1390 PetscDS ds; 1391 Parameter *param; 1392 PetscCall(PetscBagGetData(user->bag, (void **)¶m)); 1393 PetscScalar constants[NUM_CONSTANTS]; 1394 constants[SIGMA] = param->sigma; 1395 constants[V0] = param->v0; 1396 constants[T0] = param->t0; 1397 constants[X0] = param->x0; 1398 constants[M0] = param->m0; 1399 constants[Q0] = param->q0; 1400 constants[PHI0] = param->phi0; 1401 constants[POISSON] = param->poissonNumber; 1402 constants[VLASOV] = param->vlasovNumber; 1403 PetscCall(DMGetDS(dm, &ds)); 1404 PetscCall(PetscDSSetConstants(ds, NUM_CONSTANTS, constants)); 1405 } 1406 PetscFunctionReturn(PETSC_SUCCESS); 1407 } 1408 1409 static PetscErrorCode InitializeVelocities_Fake1D(DM sw, AppCtx *user) 1410 { 1411 DM dm; 1412 PetscReal *v; 1413 PetscInt *species, cStart, cEnd; 1414 PetscInt dim, Np; 1415 1416 PetscFunctionBegin; 1417 PetscCall(DMGetDimension(sw, &dim)); 1418 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 1419 PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&v)); 1420 PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 1421 PetscCall(DMSwarmGetCellDM(sw, &dm)); 1422 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 1423 PetscRandom rnd; 1424 PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)sw), &rnd)); 1425 PetscCall(PetscRandomSetInterval(rnd, 0, 1.)); 1426 PetscCall(PetscRandomSetFromOptions(rnd)); 1427 1428 for (PetscInt p = 0; p < Np; ++p) { 1429 PetscReal a[3] = {0., 0., 0.}, vel[3] = {0., 0., 0.}; 1430 1431 PetscCall(PetscRandomGetValueReal(rnd, &a[0])); 1432 if (user->perturbed_weights) { 1433 PetscCall(PetscPDFSampleConstant1D(a, NULL, vel)); 1434 } else { 1435 PetscCall(PetscPDFSampleGaussian1D(a, NULL, vel)); 1436 } 1437 v[p * dim] = vel[0]; 1438 } 1439 PetscCall(PetscRandomDestroy(&rnd)); 1440 PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&v)); 1441 PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 1442 PetscFunctionReturn(PETSC_SUCCESS); 1443 } 1444 1445 static PetscErrorCode CreateSwarm(DM dm, AppCtx *user, DM *sw) 1446 { 1447 DMSwarmCellDM celldm; 1448 DM vdm; 1449 PetscReal v0[2] = {1., 0.}; 1450 PetscInt dim; 1451 1452 PetscFunctionBeginUser; 1453 PetscCall(DMGetDimension(dm, &dim)); 1454 PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw)); 1455 PetscCall(DMSetType(*sw, DMSWARM)); 1456 PetscCall(DMSetDimension(*sw, dim)); 1457 PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC)); 1458 1459 PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR)); 1460 PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "velocity", dim, PETSC_REAL)); 1461 PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "species", 1, PETSC_INT)); 1462 PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "initCoordinates", dim, PETSC_REAL)); 1463 PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "initVelocity", dim, PETSC_REAL)); 1464 PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "E_field", dim, PETSC_REAL)); 1465 1466 const char *fieldnames[2] = {DMSwarmPICField_coor, "velocity"}; 1467 1468 PetscCall(DMSwarmCellDMCreate(dm, 2, fieldnames, 1, fieldnames, &celldm)); 1469 PetscCall(DMSwarmAddCellDM(*sw, celldm)); 1470 PetscCall(DMSwarmCellDMDestroy(&celldm)); 1471 1472 const char *vfieldnames[1] = {"w_q"}; 1473 1474 PetscCall(CreateVelocityDM(*sw, &vdm)); 1475 PetscCall(PetscObjectCompose((PetscObject)*sw, "__vdm__", (PetscObject)vdm)); 1476 PetscCall(DMSwarmCellDMCreate(vdm, 1, vfieldnames, 1, &fieldnames[1], &celldm)); 1477 PetscCall(DMSwarmAddCellDM(*sw, celldm)); 1478 PetscCall(DMSwarmCellDMDestroy(&celldm)); 1479 PetscCall(DMDestroy(&vdm)); 1480 1481 if (user->remap && user->remapType == RM_PFAK) { 1482 DM rdm; 1483 1484 PetscCall(CreateRemapDM(*sw, &rdm)); 1485 PetscCall(DMSwarmCellDMCreate(rdm, 1, vfieldnames, 2, fieldnames, &celldm)); 1486 PetscCall(DMSwarmAddCellDM(*sw, celldm)); 1487 PetscCall(DMSwarmCellDMDestroy(&celldm)); 1488 PetscCall(DMDestroy(&rdm)); 1489 } 1490 1491 PetscCall(DMSwarmFinalizeFieldRegister(*sw)); 1492 PetscCall(DMSetApplicationContext(*sw, user)); 1493 PetscCall(DMSetFromOptions(*sw)); 1494 PetscCall(DMSwarmSetCellDMActive(*sw, "space")); 1495 user->swarm = *sw; 1496 // TODO: This is redundant init since it is done in InitializeSolveAndSwarm. To remove it, we need to move the 1497 // creation of initCoordinates and initVelocity 1498 if (user->perturbed_weights) { 1499 PetscCall(InitializeParticles_PerturbedWeights(*sw, user)); 1500 } else { 1501 PetscCall(DMSwarmComputeLocalSizeFromOptions(*sw)); 1502 PetscCall(DMSwarmInitializeCoordinates(*sw)); 1503 if (user->fake_1D) { 1504 PetscCall(InitializeVelocities_Fake1D(*sw, user)); 1505 } else { 1506 PetscCall(DMSwarmInitializeVelocitiesFromOptions(*sw, v0)); 1507 } 1508 } 1509 PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles")); 1510 PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view")); 1511 { 1512 Vec gc, gc0, gv, gv0; 1513 1514 PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, DMSwarmPICField_coor, &gc)); 1515 PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, "initCoordinates", &gc0)); 1516 PetscCall(VecCopy(gc, gc0)); 1517 PetscCall(VecViewFromOptions(gc, NULL, "-ic_x_view")); 1518 PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, DMSwarmPICField_coor, &gc)); 1519 PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, "initCoordinates", &gc0)); 1520 PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, "velocity", &gv)); 1521 PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, "initVelocity", &gv0)); 1522 PetscCall(VecCopy(gv, gv0)); 1523 PetscCall(VecViewFromOptions(gv, NULL, "-ic_v_view")); 1524 PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, "velocity", &gv)); 1525 PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, "initVelocity", &gv0)); 1526 } 1527 PetscFunctionReturn(PETSC_SUCCESS); 1528 } 1529 1530 /* 1531 @article{MyersColellaVanStraalen2017, 1532 title = {A 4th-order particle-in-cell method with phase-space remapping for the {Vlasov-Poisson} equation}, 1533 author = {Andrew Myers and Phillip Colella and Brian Van Straalen}, 1534 journal = {SIAM Journal on Scientific Computing}, 1535 volume = {39}, 1536 issue = {3}, 1537 pages = {B467-B485}, 1538 doi = {10.1137/16M105962X}, 1539 issn = {10957197}, 1540 year = {2017}, 1541 } 1542 */ 1543 static PetscErrorCode W_3_Interpolation_Private(PetscReal x, PetscReal *w) 1544 { 1545 const PetscReal ax = PetscAbsReal(x); 1546 1547 PetscFunctionBegin; 1548 *w = 0.; 1549 // W_3(x) = 1 - 5/2 |x|^2 + 3/2 |x|^3 0 \le |x| \e 1 1550 if (ax <= 1.) *w = 1. - 2.5 * PetscSqr(ax) + 1.5 * PetscSqr(ax) * ax; 1551 // 1/2 (2 - |x|)^2 (1 - |x|) 1 \le |x| \le 2 1552 else if (ax <= 2.) *w = 0.5 * PetscSqr(2. - ax) * (1. - ax); 1553 //PetscCall(PetscPrintf(PETSC_COMM_SELF, " W_3 %g --> %g\n", x, *w)); 1554 PetscFunctionReturn(PETSC_SUCCESS); 1555 } 1556 1557 // Right now, we will assume that the spatial and velocity grids are regular, which will speed up point location immensely 1558 static PetscErrorCode DMSwarmRemap_Colella(DM sw, DM *rsw) 1559 { 1560 DM xdm, vdm; 1561 DMSwarmCellDM celldm; 1562 AppCtx *user; 1563 PetscReal xmin[3], xmax[3], vmin[3], vmax[3]; 1564 PetscInt xend[3], vend[3]; 1565 PetscReal *x, *v, *w, *rw; 1566 PetscReal hx[3], hv[3]; 1567 PetscInt dim, xcdim, vcdim, xcStart, xcEnd, vcStart, vcEnd, Np, Nfc; 1568 PetscInt debug = ((DM_Swarm *)sw->data)->printWeights; 1569 const char **coordFields; 1570 1571 PetscFunctionBegin; 1572 PetscCall(DMGetApplicationContext(sw, (void **)&user)); 1573 PetscCall(DMGetDimension(sw, &dim)); 1574 PetscCall(DMSwarmGetCellDM(sw, &xdm)); 1575 PetscCall(DMGetCoordinateDim(xdm, &xcdim)); 1576 // Create a new centroid swarm without weights 1577 PetscCall(CreateSwarm(xdm, user, rsw)); 1578 PetscCall(InitializeParticles_Centroid(*rsw, user->fake_1D)); 1579 PetscCall(DMSwarmGetLocalSize(*rsw, &Np)); 1580 // Assume quad mesh and calculate cell diameters (TODO this could be more robust) 1581 { 1582 const PetscScalar *array; 1583 PetscScalar *coords; 1584 PetscBool isDG; 1585 PetscInt Nc; 1586 1587 PetscCall(DMGetBoundingBox(xdm, xmin, xmax)); 1588 PetscCall(DMPlexGetHeightStratum(xdm, 0, &xcStart, &xcEnd)); 1589 PetscCall(DMPlexGetCellCoordinates(xdm, xcStart, &isDG, &Nc, &array, &coords)); 1590 hx[0] = coords[1 * xcdim + 0] - coords[0 * xcdim + 0]; 1591 hx[1] = xcdim > 1 ? coords[2 * xcdim + 1] - coords[1 * xcdim + 1] : 1.; 1592 PetscCall(DMPlexRestoreCellCoordinates(xdm, xcStart, &isDG, &Nc, &array, &coords)); 1593 PetscCall(PetscObjectQuery((PetscObject)sw, "__vdm__", (PetscObject *)&vdm)); 1594 PetscCall(DMGetCoordinateDim(vdm, &vcdim)); 1595 PetscCall(DMGetBoundingBox(vdm, vmin, vmax)); 1596 PetscCall(DMPlexGetHeightStratum(vdm, 0, &vcStart, &vcEnd)); 1597 PetscCall(DMPlexGetCellCoordinates(vdm, vcStart, &isDG, &Nc, &array, &coords)); 1598 hv[0] = coords[1 * vcdim + 0] - coords[0 * vcdim + 0]; 1599 hv[1] = vcdim > 1 ? coords[2 * vcdim + 1] - coords[1 * vcdim + 1] : 1.; 1600 PetscCall(DMPlexRestoreCellCoordinates(vdm, vcStart, &isDG, &Nc, &array, &coords)); 1601 1602 PetscCheck(dim == 1 || user->fake_1D, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONG, "Only support 1D distributions at this time"); 1603 xend[0] = xcEnd - xcStart; 1604 xend[1] = 1; 1605 vend[0] = vcEnd - vcStart; 1606 vend[1] = 1; 1607 if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Phase Grid (%g, %g, %g, %g) (%" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ")\n", hx[0], hx[1], hv[0], hv[1], xend[0], xend[1], vend[0], vend[1])); 1608 } 1609 // Iterate over particles in the original swarm 1610 PetscCall(DMSwarmGetCellDMActive(sw, &celldm)); 1611 PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields)); 1612 PetscCheck(Nfc == 1, PetscObjectComm((PetscObject)sw), PETSC_ERR_SUP, "We only support a single coordinate field right now, not %" PetscInt_FMT, Nfc); 1613 PetscCall(DMSwarmGetField(sw, coordFields[0], NULL, NULL, (void **)&x)); 1614 PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&v)); 1615 PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&w)); 1616 PetscCall(DMSwarmGetField(*rsw, "w_q", NULL, NULL, (void **)&rw)); 1617 PetscCall(DMSwarmSortGetAccess(sw)); 1618 PetscCall(DMSwarmSortGetAccess(*rsw)); 1619 PetscCall(DMGetBoundingBox(vdm, vmin, vmax)); 1620 PetscCall(DMGetCoordinatesLocalSetUp(xdm)); 1621 for (PetscInt i = 0; i < Np; ++i) rw[i] = 0.; 1622 for (PetscInt c = 0; c < xcEnd - xcStart; ++c) { 1623 PetscInt *pidx, Npc; 1624 PetscInt *rpidx, rNpc; 1625 1626 PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Npc, &pidx)); 1627 for (PetscInt q = 0; q < Npc; ++q) { 1628 const PetscInt p = pidx[q]; 1629 const PetscReal wp = w[p]; 1630 PetscReal Wx[3], Wv[3]; 1631 PetscInt xs[3], vs[3]; 1632 1633 // Determine the containing cell 1634 for (PetscInt d = 0; d < dim; ++d) { 1635 const PetscReal xp = x[p * dim + d]; 1636 const PetscReal vp = v[p * dim + d]; 1637 1638 xs[d] = PetscFloorReal((xp - xmin[d]) / hx[d]); 1639 vs[d] = PetscFloorReal((vp - vmin[d]) / hv[d]); 1640 } 1641 // Loop over all grid points within 2 spacings of the particle 1642 if (debug > 2) { 1643 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Interpolating particle %" PetscInt_FMT " wt %g (%g, %g, %g, %g) (%" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ")\n", p, wp, x[p * dim + 0], xcdim > 1 ? x[p * xcdim + 1] : 0., v[p * vcdim + 0], vcdim > 1 ? v[p * vcdim + 1] : 0., xs[0], xs[1], vs[0], vs[1])); 1644 } 1645 for (PetscInt xi = xs[0] - 1; xi < xs[0] + 3; ++xi) { 1646 // Treat xi as periodic 1647 const PetscInt xip = xi < 0 ? xi + xend[0] : (xi >= xend[0] ? xi - xend[0] : xi); 1648 PetscCall(W_3_Interpolation_Private((xmin[0] + (xi + 0.5) * hx[0] - x[p * dim + 0]) / hx[0], &Wx[0])); 1649 for (PetscInt xj = PetscMax(xs[1] - 1, 0); xj < PetscMin(xs[1] + 3, xend[1]); ++xj) { 1650 if (xcdim > 1) PetscCall(W_3_Interpolation_Private((xmin[1] + (xj + 0.5) * hx[1] - x[p * dim + 1]) / hx[1], &Wx[1])); 1651 else Wx[1] = 1.; 1652 for (PetscInt vi = PetscMax(vs[0] - 1, 0); vi < PetscMin(vs[0] + 3, vend[0]); ++vi) { 1653 PetscCall(W_3_Interpolation_Private((vmin[0] + (vi + 0.5) * hv[0] - v[p * dim + 0]) / hv[0], &Wv[0])); 1654 for (PetscInt vj = PetscMax(vs[1] - 1, 0); vj < PetscMin(vs[1] + 3, vend[1]); ++vj) { 1655 const PetscInt rc = xip * xend[1] + xj; 1656 const PetscInt rv = vi * vend[1] + vj; 1657 1658 PetscCall(DMSwarmSortGetPointsPerCell(*rsw, rc, &rNpc, &rpidx)); 1659 if (vcdim > 1) PetscCall(W_3_Interpolation_Private((vmin[1] + (vj + 0.5) * hv[1] - v[p * dim + 1]) / hv[1], &Wv[1])); 1660 else Wv[1] = 1.; 1661 if (debug > 2) 1662 PetscCall(PetscPrintf(PETSC_COMM_SELF, " Depositing on particle (%" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ") w = %g (%g, %g, %g, %g)\n", xi, xj, vi, vj, wp * Wx[0] * Wx[1] * Wv[0] * Wv[1], Wx[0], Wx[1], Wv[0], Wv[1])); 1663 // Add weight to new particles from original particle using interpolation function 1664 PetscCheck(rNpc == vend[0] * vend[1], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid particle velocity binning"); 1665 const PetscInt rp = rpidx[rv]; 1666 PetscCheck(rp >= 0 && rp < Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Particle index %" PetscInt_FMT " not in [0, %" PetscInt_FMT ")", rp, Np); 1667 rw[rp] += wp * Wx[0] * Wx[1] * Wv[0] * Wv[1]; 1668 if (debug > 2) PetscCall(PetscPrintf(PETSC_COMM_SELF, " Adding weight %g (%g) to particle %" PetscInt_FMT "\n", wp * Wx[0] * Wx[1] * Wv[0] * Wv[1], rw[rp], rp)); 1669 PetscCall(DMSwarmSortRestorePointsPerCell(*rsw, rc, &rNpc, &rpidx)); 1670 } 1671 } 1672 } 1673 } 1674 } 1675 PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Npc, &pidx)); 1676 } 1677 PetscCall(DMSwarmSortRestoreAccess(sw)); 1678 PetscCall(DMSwarmSortRestoreAccess(*rsw)); 1679 PetscCall(DMSwarmRestoreField(sw, coordFields[0], NULL, NULL, (void **)&x)); 1680 PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&v)); 1681 PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&w)); 1682 PetscCall(DMSwarmRestoreField(*rsw, "w_q", NULL, NULL, (void **)&rw)); 1683 1684 if (debug) { 1685 Vec w; 1686 1687 PetscCall(DMSwarmCreateGlobalVectorFromField(sw, coordFields[0], &w)); 1688 PetscCall(VecViewFromOptions(w, NULL, "-remap_view")); 1689 PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, coordFields[0], &w)); 1690 PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "velocity", &w)); 1691 PetscCall(VecViewFromOptions(w, NULL, "-remap_view")); 1692 PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "velocity", &w)); 1693 PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &w)); 1694 PetscCall(VecViewFromOptions(w, NULL, "-remap_view")); 1695 PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &w)); 1696 PetscCall(DMSwarmCreateGlobalVectorFromField(*rsw, coordFields[0], &w)); 1697 PetscCall(VecViewFromOptions(w, NULL, "-remap_view")); 1698 PetscCall(DMSwarmDestroyGlobalVectorFromField(*rsw, coordFields[0], &w)); 1699 PetscCall(DMSwarmCreateGlobalVectorFromField(*rsw, "velocity", &w)); 1700 PetscCall(VecViewFromOptions(w, NULL, "-remap_view")); 1701 PetscCall(DMSwarmDestroyGlobalVectorFromField(*rsw, "velocity", &w)); 1702 PetscCall(DMSwarmCreateGlobalVectorFromField(*rsw, "w_q", &w)); 1703 PetscCall(VecViewFromOptions(w, NULL, "-remap_view")); 1704 PetscCall(DMSwarmDestroyGlobalVectorFromField(*rsw, "w_q", &w)); 1705 } 1706 PetscFunctionReturn(PETSC_SUCCESS); 1707 } 1708 1709 static void f0_v2(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[]) 1710 { 1711 PetscInt d; 1712 1713 f0[0] = 0.0; 1714 for (d = dim / 2; d < dim; ++d) f0[0] += PetscSqr(x[d]) * u[0]; 1715 } 1716 1717 static PetscErrorCode DMSwarmRemap_PFAK(DM sw, DM *rsw) 1718 { 1719 DM xdm, vdm, rdm; 1720 DMSwarmCellDM rcelldm; 1721 Mat M_p, rM_p, rPM_p; 1722 Vec w, rw, rhs; 1723 PetscInt Nf; 1724 const char **fields; 1725 AppCtx *user; 1726 1727 PetscFunctionBegin; 1728 // Create a new centroid swarm without weights 1729 PetscCall(DMGetApplicationContext(sw, (void **)&user)); 1730 PetscCall(DMSwarmGetCellDM(sw, &xdm)); 1731 PetscCall(PetscObjectQuery((PetscObject)sw, "__vdm__", (PetscObject *)&vdm)); 1732 PetscCall(CreateSwarm(xdm, user, rsw)); 1733 // Set remap cell DM 1734 PetscCall(DMSwarmSetCellDMActive(sw, "remap")); 1735 PetscCall(DMSwarmGetCellDMActive(sw, &rcelldm)); 1736 PetscCall(DMSwarmCellDMGetFields(rcelldm, &Nf, &fields)); 1737 PetscCheck(Nf == 1, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONG, "We only allow a single weight field, not %" PetscInt_FMT, Nf); 1738 PetscCall(DMSwarmGetCellDM(sw, &rdm)); 1739 PetscCall(DMGetGlobalVector(rdm, &rhs)); 1740 PetscCall(DMSwarmMigrate(sw, PETSC_FALSE)); // Bin particles in remap mesh 1741 // Compute rhs = M_p w_p 1742 PetscCall(DMCreateMassMatrix(sw, rdm, &M_p)); 1743 PetscCall(DMSwarmCreateGlobalVectorFromField(sw, fields[0], &w)); 1744 PetscCall(VecViewFromOptions(w, NULL, "-remap_w_view")); 1745 PetscCall(MatMultTranspose(M_p, w, rhs)); 1746 PetscCall(VecViewFromOptions(rhs, NULL, "-remap_rhs_view")); 1747 PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, fields[0], &w)); 1748 PetscCall(MatDestroy(&M_p)); 1749 { 1750 KSP ksp; 1751 Mat M_f; 1752 Vec u_f; 1753 PetscReal mom[4]; 1754 PetscInt cdim; 1755 const char *prefix; 1756 1757 PetscCall(DMGetCoordinateDim(rdm, &cdim)); 1758 PetscCall(DMCreateMassMatrix(rdm, rdm, &M_f)); 1759 PetscCall(DMGetGlobalVector(rdm, &u_f)); 1760 1761 PetscCall(KSPCreate(PetscObjectComm((PetscObject)sw), &ksp)); 1762 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sw, &prefix)); 1763 PetscCall(KSPSetOptionsPrefix(ksp, prefix)); 1764 PetscCall(KSPAppendOptionsPrefix(ksp, "ptof_")); 1765 PetscCall(KSPSetFromOptions(ksp)); 1766 1767 PetscCall(KSPSetOperators(ksp, M_f, M_f)); 1768 PetscCall(KSPSolve(ksp, rhs, u_f)); 1769 PetscCall(KSPDestroy(&ksp)); 1770 PetscCall(VecViewFromOptions(u_f, NULL, "-remap_uf_view")); 1771 1772 PetscCall(DMPlexComputeMoments(rdm, u_f, mom)); 1773 // Energy is not correct since it uses (x^2 + v^2) 1774 PetscDS rds; 1775 PetscScalar rmom; 1776 void *ctx; 1777 1778 PetscCall(DMGetDS(rdm, &rds)); 1779 PetscCall(DMGetApplicationContext(rdm, &ctx)); 1780 PetscCall(PetscDSSetObjective(rds, 0, &f0_v2)); 1781 PetscCall(DMPlexComputeIntegralFEM(rdm, u_f, &rmom, ctx)); 1782 mom[1 + cdim] = PetscRealPart(rmom); 1783 1784 PetscCall(DMRestoreGlobalVector(rdm, &u_f)); 1785 PetscCall(PetscPrintf(PETSC_COMM_SELF, "========== PFAK u_f ==========\n")); 1786 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Mom 0: %g\n", mom[0])); 1787 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Mom x: %g\n", mom[1 + 0])); 1788 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Mom v: %g\n", mom[1 + 1])); 1789 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Mom 2: %g\n", mom[1 + cdim])); 1790 PetscCall(MatDestroy(&M_f)); 1791 } 1792 // Create Remap particle mass matrix M_p 1793 PetscInt xcStart, xcEnd, vcStart, vcEnd, cStart, cEnd, r; 1794 1795 PetscCall(DMSwarmSetCellDMActive(*rsw, "remap")); 1796 PetscCall(DMPlexGetHeightStratum(xdm, 0, &xcStart, &xcEnd)); 1797 PetscCall(DMPlexGetHeightStratum(vdm, 0, &vcStart, &vcEnd)); 1798 PetscCall(DMPlexGetHeightStratum(rdm, 0, &cStart, &cEnd)); 1799 r = (PetscInt)PetscSqrtReal(((xcEnd - xcStart) * (vcEnd - vcStart)) / (cEnd - cStart)); 1800 PetscCall(InitializeParticles_Regular(*rsw, r)); 1801 PetscCall(DMSwarmMigrate(*rsw, PETSC_FALSE)); // Bin particles in remap mesh 1802 PetscCall(DMCreateMassMatrix(*rsw, rdm, &rM_p)); 1803 PetscCall(MatViewFromOptions(rM_p, NULL, "-rM_p_view")); 1804 // Solve M_p 1805 { 1806 KSP ksp; 1807 PC pc; 1808 const char *prefix; 1809 PetscBool isBjacobi; 1810 1811 PetscCall(KSPCreate(PetscObjectComm((PetscObject)sw), &ksp)); 1812 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sw, &prefix)); 1813 PetscCall(KSPSetOptionsPrefix(ksp, prefix)); 1814 PetscCall(KSPAppendOptionsPrefix(ksp, "ftop_")); 1815 PetscCall(KSPSetFromOptions(ksp)); 1816 1817 PetscCall(KSPGetPC(ksp, &pc)); 1818 PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCBJACOBI, &isBjacobi)); 1819 if (isBjacobi) { 1820 PetscCall(DMSwarmCreateMassMatrixSquare(sw, rdm, &rPM_p)); 1821 } else { 1822 rPM_p = rM_p; 1823 PetscCall(PetscObjectReference((PetscObject)rPM_p)); 1824 } 1825 PetscCall(KSPSetOperators(ksp, rM_p, rPM_p)); 1826 PetscCall(DMSwarmCreateGlobalVectorFromField(*rsw, fields[0], &rw)); 1827 PetscCall(KSPSolveTranspose(ksp, rhs, rw)); 1828 PetscCall(VecViewFromOptions(rw, NULL, "-remap_rw_view")); 1829 PetscCall(DMSwarmDestroyGlobalVectorFromField(*rsw, fields[0], &rw)); 1830 PetscCall(KSPDestroy(&ksp)); 1831 PetscCall(MatDestroy(&rPM_p)); 1832 PetscCall(MatDestroy(&rM_p)); 1833 } 1834 PetscCall(DMRestoreGlobalVector(rdm, &rhs)); 1835 1836 // Restore original cell DM 1837 PetscCall(DMSwarmSetCellDMActive(sw, "space")); 1838 PetscCall(DMSwarmSetCellDMActive(*rsw, "space")); 1839 PetscCall(DMSwarmMigrate(*rsw, PETSC_FALSE)); // Bin particles in spatial mesh 1840 PetscFunctionReturn(PETSC_SUCCESS); 1841 } 1842 1843 static PetscErrorCode DMSwarmRemap(DM sw) 1844 { 1845 DM rsw; 1846 AppCtx *user; 1847 1848 PetscFunctionBegin; 1849 PetscCall(DMGetApplicationContext(sw, (void **)&user)); 1850 switch (user->remapType) { 1851 case RM_COLELLA: 1852 PetscCall(DMSwarmRemap_Colella(sw, &rsw)); 1853 break; 1854 case RM_PFAK: 1855 PetscCall(DMSwarmRemap_PFAK(sw, &rsw)); 1856 break; 1857 default: 1858 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No remap algorithmic %s", RemapTypes[user->remapType]); 1859 } 1860 1861 PetscReal mom[4], rmom[4]; 1862 PetscInt cdim; 1863 1864 PetscCall(DMGetCoordinateDim(sw, &cdim)); 1865 PetscCall(DMSwarmComputeMoments(sw, "velocity", "w_q", mom)); 1866 PetscCall(DMSwarmComputeMoments(rsw, "velocity", "w_q", rmom)); 1867 PetscCall(PetscPrintf(PETSC_COMM_SELF, "========== Remapped ==========\n")); 1868 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Mom 0: %g --> %g\n", mom[0], rmom[0])); 1869 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Mom 1: %g --> %g\n", mom[1], rmom[1])); 1870 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Mom 2: %g --> %g\n", mom[1 + cdim], rmom[1 + cdim])); 1871 if (user->weight_monitor) { 1872 PetscCall(MonitorSwarmWeights(sw, "w_q")); 1873 PetscCall(MonitorSwarmWeights(rsw, "w_q")); 1874 } 1875 PetscCall(DMSwarmReplace(sw, &rsw)); 1876 PetscFunctionReturn(PETSC_SUCCESS); 1877 } 1878 1879 static PetscErrorCode ComputeFieldAtParticles_Coulomb(SNES snes, DM sw, PetscReal E[]) 1880 { 1881 AppCtx *user; 1882 PetscReal *coords; 1883 PetscInt *species, dim, Np, Ns; 1884 PetscMPIInt size; 1885 1886 PetscFunctionBegin; 1887 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)snes), &size)); 1888 PetscCheck(size == 1, PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "Coulomb code only works in serial"); 1889 PetscCall(DMGetDimension(sw, &dim)); 1890 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 1891 PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 1892 PetscCall(DMGetApplicationContext(sw, (void *)&user)); 1893 1894 PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 1895 PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 1896 for (PetscInt p = 0; p < Np; ++p) { 1897 PetscReal *pcoord = &coords[p * dim]; 1898 PetscReal pE[3] = {0., 0., 0.}; 1899 1900 /* Calculate field at particle p due to particle q */ 1901 for (PetscInt q = 0; q < Np; ++q) { 1902 PetscReal *qcoord = &coords[q * dim]; 1903 PetscReal rpq[3], r, r3, q_q; 1904 1905 if (p == q) continue; 1906 q_q = user->charges[species[q]] * 1.; 1907 for (PetscInt d = 0; d < dim; ++d) rpq[d] = pcoord[d] - qcoord[d]; 1908 r = DMPlex_NormD_Internal(dim, rpq); 1909 if (r < PETSC_SQRT_MACHINE_EPSILON) continue; 1910 r3 = PetscPowRealInt(r, 3); 1911 for (PetscInt d = 0; d < dim; ++d) pE[d] += q_q * rpq[d] / r3; 1912 } 1913 for (PetscInt d = 0; d < dim; ++d) E[p * dim + d] = pE[d]; 1914 } 1915 PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 1916 PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 1917 PetscFunctionReturn(PETSC_SUCCESS); 1918 } 1919 1920 static PetscErrorCode ComputeFieldAtParticles_Primal(SNES snes, DM sw, PetscReal E[]) 1921 { 1922 DM dm; 1923 AppCtx *user; 1924 PetscDS ds; 1925 PetscFE fe; 1926 Mat M_p; 1927 Vec rhoRhs; // Weak charge density, \int phi_i rho 1928 Vec rho; // Charge density, M^{-1} rhoRhs 1929 Vec phi, locPhi; // Potential 1930 Vec f; // Particle weights 1931 PetscReal *coords; 1932 PetscInt dim, cStart, cEnd, Np; 1933 1934 PetscFunctionBegin; 1935 PetscCall(DMGetApplicationContext(sw, (void *)&user)); 1936 PetscCall(PetscLogEventBegin(user->ESolveEvent, snes, sw, 0, 0)); 1937 PetscCall(DMGetDimension(sw, &dim)); 1938 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 1939 1940 KSP ksp; 1941 const char **oldFields; 1942 PetscInt Nf; 1943 const char **tmp; 1944 1945 PetscCall(SNESGetDM(snes, &dm)); 1946 PetscCall(DMSwarmVectorGetField(sw, &Nf, &tmp)); 1947 PetscCall(PetscMalloc1(Nf, &oldFields)); 1948 for (PetscInt f = 0; f < Nf; ++f) PetscCall(PetscStrallocpy(tmp[f], (char **)&oldFields[f])); 1949 PetscCall(DMSwarmVectorDefineField(sw, "w_q")); 1950 PetscCall(DMCreateMassMatrix(sw, dm, &M_p)); 1951 PetscCall(DMSwarmVectorDefineFields(sw, Nf, oldFields)); 1952 for (PetscInt f = 0; f < Nf; ++f) PetscCall(PetscFree(oldFields[f])); 1953 PetscCall(PetscFree(oldFields)); 1954 1955 PetscCall(DMGetGlobalVector(dm, &rhoRhs)); 1956 PetscCall(PetscObjectSetName((PetscObject)rhoRhs, "Weak charge density")); 1957 PetscCall(DMGetNamedGlobalVector(dm, "rho", &rho)); 1958 PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f)); 1959 PetscCall(PetscObjectSetName((PetscObject)f, "particle weight")); 1960 1961 PetscCall(MatViewFromOptions(M_p, NULL, "-mp_view")); 1962 PetscCall(MatViewFromOptions(user->M, NULL, "-m_view")); 1963 PetscCall(VecViewFromOptions(f, NULL, "-weights_view")); 1964 1965 PetscCall(MatMultTranspose(M_p, f, rhoRhs)); 1966 PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f)); 1967 1968 PetscCall(KSPCreate(PetscObjectComm((PetscObject)dm), &ksp)); 1969 PetscCall(KSPSetOptionsPrefix(ksp, "em_proj_")); 1970 PetscCall(KSPSetOperators(ksp, user->M, user->M)); 1971 PetscCall(KSPSetFromOptions(ksp)); 1972 PetscCall(KSPSolve(ksp, rhoRhs, rho)); 1973 1974 PetscCall(VecScale(rhoRhs, -1.0)); 1975 1976 PetscCall(VecViewFromOptions(rhoRhs, NULL, "-rho_view")); 1977 PetscCall(DMRestoreNamedGlobalVector(dm, "rho", &rho)); 1978 PetscCall(KSPDestroy(&ksp)); 1979 PetscCall(MatDestroy(&M_p)); 1980 1981 PetscCall(DMGetNamedGlobalVector(dm, "phi", &phi)); 1982 PetscCall(VecSet(phi, 0.0)); 1983 PetscCall(SNESSolve(snes, rhoRhs, phi)); 1984 PetscCall(DMRestoreGlobalVector(dm, &rhoRhs)); 1985 PetscCall(VecViewFromOptions(phi, NULL, "-phi_view")); 1986 1987 PetscCall(DMGetLocalVector(dm, &locPhi)); 1988 PetscCall(DMGlobalToLocalBegin(dm, phi, INSERT_VALUES, locPhi)); 1989 PetscCall(DMGlobalToLocalEnd(dm, phi, INSERT_VALUES, locPhi)); 1990 PetscCall(DMRestoreNamedGlobalVector(dm, "phi", &phi)); 1991 PetscCall(PetscLogEventEnd(user->ESolveEvent, snes, sw, 0, 0)); 1992 1993 PetscCall(DMGetDS(dm, &ds)); 1994 PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&fe)); 1995 PetscCall(DMSwarmSortGetAccess(sw)); 1996 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 1997 PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 1998 1999 PetscCall(PetscLogEventBegin(user->ETabEvent, snes, sw, 0, 0)); 2000 PetscFEGeom *chunkgeom = NULL; 2001 for (PetscInt c = cStart; c < cEnd; ++c) { 2002 PetscTabulation tab; 2003 PetscScalar *clPhi = NULL; 2004 PetscReal *pcoord, *refcoord; 2005 PetscInt *points; 2006 PetscInt Ncp; 2007 2008 PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Ncp, &points)); 2009 PetscCall(DMGetWorkArray(dm, Ncp * dim, MPIU_REAL, &pcoord)); 2010 PetscCall(DMGetWorkArray(dm, Ncp * dim, MPIU_REAL, &refcoord)); 2011 for (PetscInt cp = 0; cp < Ncp; ++cp) 2012 for (PetscInt d = 0; d < dim; ++d) pcoord[cp * dim + d] = coords[points[cp] * dim + d]; 2013 { 2014 PetscCall(PetscFEGeomGetChunk(user->fegeom, c - cStart, c - cStart + 1, &chunkgeom)); 2015 for (PetscInt i = 0; i < Ncp; ++i) { 2016 const PetscReal x0[3] = {-1., -1., -1.}; 2017 CoordinatesRealToRef(dim, dim, x0, chunkgeom->v, chunkgeom->invJ, &pcoord[dim * i], &refcoord[dim * i]); 2018 } 2019 } 2020 PetscCall(PetscFECreateTabulation(fe, 1, Ncp, refcoord, 1, &tab)); 2021 PetscCall(DMPlexVecGetClosure(dm, NULL, locPhi, c, NULL, &clPhi)); 2022 for (PetscInt cp = 0; cp < Ncp; ++cp) { 2023 const PetscReal *basisDer = tab->T[1]; 2024 const PetscInt p = points[cp]; 2025 2026 for (PetscInt d = 0; d < dim; ++d) E[p * dim + d] = 0.; 2027 PetscCall(PetscFEFreeInterpolateGradient_Static(fe, basisDer, clPhi, dim, chunkgeom->invJ, NULL, cp, &E[p * dim])); 2028 for (PetscInt d = 0; d < dim; ++d) { 2029 E[p * dim + d] *= -1.0; 2030 if (user->fake_1D && d > 0) E[p * dim + d] = 0; 2031 } 2032 } 2033 PetscCall(DMPlexVecRestoreClosure(dm, NULL, locPhi, c, NULL, &clPhi)); 2034 PetscCall(PetscTabulationDestroy(&tab)); 2035 PetscCall(DMRestoreWorkArray(dm, Ncp * dim, MPIU_REAL, &pcoord)); 2036 PetscCall(DMRestoreWorkArray(dm, Ncp * dim, MPIU_REAL, &refcoord)); 2037 PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Ncp, &points)); 2038 } 2039 PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 2040 PetscCall(DMSwarmSortRestoreAccess(sw)); 2041 PetscCall(DMRestoreLocalVector(dm, &locPhi)); 2042 PetscCall(PetscFEGeomRestoreChunk(user->fegeom, 0, 1, &chunkgeom)); 2043 PetscCall(PetscLogEventEnd(user->ETabEvent, snes, sw, 0, 0)); 2044 PetscFunctionReturn(PETSC_SUCCESS); 2045 } 2046 2047 static PetscErrorCode ComputeFieldAtParticles_Mixed(SNES snes, DM sw, PetscReal E[]) 2048 { 2049 AppCtx *user; 2050 DM dm, potential_dm; 2051 KSP ksp; 2052 IS potential_IS; 2053 PetscDS ds; 2054 PetscFE fe; 2055 Mat M_p, M; 2056 Vec phi, locPhi, rho, f, temp_rho, rho0; 2057 PetscQuadrature q; 2058 PetscReal *coords; 2059 PetscInt dim, cStart, cEnd, Np, pot_field = 1; 2060 const char **oldFields; 2061 PetscInt Nf; 2062 const char **tmp; 2063 2064 PetscFunctionBegin; 2065 PetscCall(DMGetApplicationContext(sw, &user)); 2066 PetscCall(PetscLogEventBegin(user->ESolveEvent, snes, sw, 0, 0)); 2067 PetscCall(DMGetDimension(sw, &dim)); 2068 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 2069 2070 /* Create the charges rho */ 2071 PetscCall(SNESGetDM(snes, &dm)); 2072 PetscCall(DMGetGlobalVector(dm, &rho)); 2073 PetscCall(PetscObjectSetName((PetscObject)rho, "rho")); 2074 2075 PetscCall(DMCreateSubDM(dm, 1, &pot_field, &potential_IS, &potential_dm)); 2076 2077 PetscCall(DMSwarmVectorGetField(sw, &Nf, &tmp)); 2078 PetscCall(PetscMalloc1(Nf, &oldFields)); 2079 for (PetscInt f = 0; f < Nf; ++f) PetscCall(PetscStrallocpy(tmp[f], (char **)&oldFields[f])); 2080 PetscCall(DMSwarmVectorDefineField(sw, "w_q")); 2081 PetscCall(DMCreateMassMatrix(sw, potential_dm, &M_p)); 2082 PetscCall(DMSwarmVectorDefineFields(sw, Nf, oldFields)); 2083 for (PetscInt f = 0; f < Nf; ++f) PetscCall(PetscFree(oldFields[f])); 2084 PetscCall(PetscFree(oldFields)); 2085 2086 PetscCall(DMCreateMassMatrix(potential_dm, potential_dm, &M)); 2087 PetscCall(MatViewFromOptions(M_p, NULL, "-mp_view")); 2088 PetscCall(MatViewFromOptions(M, NULL, "-m_view")); 2089 PetscCall(DMGetGlobalVector(potential_dm, &temp_rho)); 2090 PetscCall(PetscObjectSetName((PetscObject)temp_rho, "Mf")); 2091 PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f)); 2092 PetscCall(PetscObjectSetName((PetscObject)f, "particle weight")); 2093 PetscCall(VecViewFromOptions(f, NULL, "-weights_view")); 2094 PetscCall(MatMultTranspose(M_p, f, temp_rho)); 2095 PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f)); 2096 PetscCall(DMGetGlobalVector(potential_dm, &rho0)); 2097 PetscCall(PetscObjectSetName((PetscObject)rho0, "Charge density (rho0) from Mixed Compute")); 2098 2099 PetscCall(KSPCreate(PetscObjectComm((PetscObject)dm), &ksp)); 2100 PetscCall(KSPSetOptionsPrefix(ksp, "em_proj")); 2101 PetscCall(KSPSetOperators(ksp, M, M)); 2102 PetscCall(KSPSetFromOptions(ksp)); 2103 PetscCall(KSPSolve(ksp, temp_rho, rho0)); 2104 PetscCall(VecViewFromOptions(rho0, NULL, "-rho0_view")); 2105 2106 PetscCall(VecISCopy(rho, potential_IS, SCATTER_FORWARD, temp_rho)); 2107 PetscCall(VecScale(rho, 0.25)); 2108 PetscCall(VecViewFromOptions(rho0, NULL, "-rho0_view")); 2109 PetscCall(VecViewFromOptions(temp_rho, NULL, "-temprho_view")); 2110 PetscCall(VecViewFromOptions(rho, NULL, "-rho_view")); 2111 PetscCall(DMRestoreGlobalVector(potential_dm, &temp_rho)); 2112 PetscCall(DMRestoreGlobalVector(potential_dm, &rho0)); 2113 2114 PetscCall(MatDestroy(&M_p)); 2115 PetscCall(MatDestroy(&M)); 2116 PetscCall(KSPDestroy(&ksp)); 2117 PetscCall(DMDestroy(&potential_dm)); 2118 PetscCall(ISDestroy(&potential_IS)); 2119 2120 PetscCall(DMGetGlobalVector(dm, &phi)); 2121 PetscCall(PetscObjectSetName((PetscObject)phi, "potential")); 2122 PetscCall(VecSet(phi, 0.0)); 2123 PetscCall(SNESSolve(snes, rho, phi)); 2124 PetscCall(DMRestoreGlobalVector(dm, &rho)); 2125 2126 PetscCall(VecViewFromOptions(phi, NULL, "-phi_view")); 2127 2128 PetscCall(DMGetLocalVector(dm, &locPhi)); 2129 PetscCall(DMGlobalToLocalBegin(dm, phi, INSERT_VALUES, locPhi)); 2130 PetscCall(DMGlobalToLocalEnd(dm, phi, INSERT_VALUES, locPhi)); 2131 PetscCall(DMRestoreGlobalVector(dm, &phi)); 2132 PetscCall(PetscLogEventEnd(user->ESolveEvent, snes, sw, 0, 0)); 2133 2134 PetscCall(PetscLogEventBegin(user->ETabEvent, snes, sw, 0, 0)); 2135 PetscCall(DMGetDS(dm, &ds)); 2136 PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&fe)); 2137 PetscCall(DMSwarmSortGetAccess(sw)); 2138 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 2139 PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 2140 PetscCall(PetscFEGetQuadrature(fe, &q)); 2141 PetscFEGeom *chunkgeom = NULL; 2142 for (PetscInt c = cStart; c < cEnd; ++c) { 2143 PetscTabulation tab; 2144 PetscScalar *clPhi = NULL; 2145 PetscReal *pcoord, *refcoord; 2146 PetscInt *points; 2147 PetscInt Ncp; 2148 2149 PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Ncp, &points)); 2150 PetscCall(DMGetWorkArray(dm, Ncp * dim, MPIU_REAL, &pcoord)); 2151 PetscCall(DMGetWorkArray(dm, Ncp * dim, MPIU_REAL, &refcoord)); 2152 for (PetscInt cp = 0; cp < Ncp; ++cp) 2153 for (PetscInt d = 0; d < dim; ++d) pcoord[cp * dim + d] = coords[points[cp] * dim + d]; 2154 { 2155 PetscCall(PetscFEGeomGetChunk(user->fegeom, c - cStart, c - cStart + 1, &chunkgeom)); 2156 for (PetscInt i = 0; i < Ncp; ++i) { 2157 // Apply the inverse affine transformation for each point 2158 const PetscReal x0[3] = {-1., -1., -1.}; 2159 CoordinatesRealToRef(dim, dim, x0, chunkgeom->v, chunkgeom->invJ, &pcoord[dim * i], &refcoord[dim * i]); 2160 } 2161 } 2162 PetscCall(PetscFECreateTabulation(fe, 1, Ncp, refcoord, 1, &tab)); 2163 PetscCall(DMPlexVecGetClosure(dm, NULL, locPhi, c, NULL, &clPhi)); 2164 2165 for (PetscInt cp = 0; cp < Ncp; ++cp) { 2166 const PetscInt p = points[cp]; 2167 2168 for (PetscInt d = 0; d < dim; ++d) E[p * dim + d] = 0.; 2169 PetscCall(PetscFEInterpolateAtPoints_Static(fe, tab, clPhi, chunkgeom, cp, &E[p * dim])); 2170 PetscCall(PetscFEPushforward(fe, chunkgeom, 1, &E[p * dim])); 2171 for (PetscInt d = 0; d < dim; ++d) { 2172 E[p * dim + d] *= -2.0; 2173 if (user->fake_1D && d > 0) E[p * dim + d] = 0; 2174 } 2175 } 2176 PetscCall(DMPlexVecRestoreClosure(dm, NULL, locPhi, c, NULL, &clPhi)); 2177 PetscCall(PetscTabulationDestroy(&tab)); 2178 PetscCall(DMRestoreWorkArray(dm, Ncp * dim, MPIU_REAL, &pcoord)); 2179 PetscCall(DMRestoreWorkArray(dm, Ncp * dim, MPIU_REAL, &refcoord)); 2180 PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Ncp, &points)); 2181 } 2182 PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 2183 PetscCall(DMSwarmSortRestoreAccess(sw)); 2184 PetscCall(DMRestoreLocalVector(dm, &locPhi)); 2185 PetscCall(PetscFEGeomRestoreChunk(user->fegeom, 0, 1, &chunkgeom)); 2186 PetscCall(PetscLogEventEnd(user->ETabEvent, snes, sw, 0, 0)); 2187 PetscFunctionReturn(PETSC_SUCCESS); 2188 } 2189 2190 static PetscErrorCode ComputeFieldAtParticles(SNES snes, DM sw, PetscReal E[]) 2191 { 2192 AppCtx *ctx; 2193 PetscInt dim, Np; 2194 2195 PetscFunctionBegin; 2196 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 2197 PetscValidHeaderSpecific(sw, DM_CLASSID, 2); 2198 PetscAssertPointer(E, 3); 2199 PetscCall(DMGetDimension(sw, &dim)); 2200 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 2201 PetscCall(DMGetApplicationContext(sw, &ctx)); 2202 PetscCall(PetscArrayzero(E, Np * dim)); 2203 ctx->validE = PETSC_TRUE; 2204 2205 switch (ctx->em) { 2206 case EM_PRIMAL: 2207 PetscCall(ComputeFieldAtParticles_Primal(snes, sw, E)); 2208 break; 2209 case EM_COULOMB: 2210 PetscCall(ComputeFieldAtParticles_Coulomb(snes, sw, E)); 2211 break; 2212 case EM_MIXED: 2213 PetscCall(ComputeFieldAtParticles_Mixed(snes, sw, E)); 2214 break; 2215 case EM_NONE: 2216 break; 2217 default: 2218 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No solver for electrostatic model %s", EMTypes[ctx->em]); 2219 } 2220 PetscFunctionReturn(PETSC_SUCCESS); 2221 } 2222 2223 static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec G, void *ctx) 2224 { 2225 DM sw; 2226 SNES snes = ((AppCtx *)ctx)->snes; 2227 const PetscReal *coords, *vel; 2228 const PetscScalar *u; 2229 PetscScalar *g; 2230 PetscReal *E, m_p = 1., q_p = -1.; 2231 PetscInt dim, d, Np, p; 2232 2233 PetscFunctionBeginUser; 2234 PetscCall(TSGetDM(ts, &sw)); 2235 PetscCall(DMGetDimension(sw, &dim)); 2236 PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **)&coords)); 2237 PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **)&vel)); 2238 PetscCall(DMSwarmGetField(sw, "E_field", NULL, NULL, (void **)&E)); 2239 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 2240 PetscCall(VecGetArrayRead(U, &u)); 2241 PetscCall(VecGetArray(G, &g)); 2242 2243 PetscCall(ComputeFieldAtParticles(snes, sw, E)); 2244 2245 Np /= 2 * dim; 2246 for (p = 0; p < Np; ++p) { 2247 for (d = 0; d < dim; ++d) { 2248 g[(p * 2 + 0) * dim + d] = u[(p * 2 + 1) * dim + d]; 2249 g[(p * 2 + 1) * dim + d] = q_p * E[p * dim + d] / m_p; 2250 } 2251 } 2252 PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **)&coords)); 2253 PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **)&vel)); 2254 PetscCall(DMSwarmRestoreField(sw, "E_field", NULL, NULL, (void **)&E)); 2255 PetscCall(VecRestoreArrayRead(U, &u)); 2256 PetscCall(VecRestoreArray(G, &g)); 2257 PetscFunctionReturn(PETSC_SUCCESS); 2258 } 2259 2260 /* J_{ij} = dF_i/dx_j 2261 J_p = ( 0 1) 2262 (-w^2 0) 2263 TODO Now there is another term with w^2 from the electric field. I think we will need to invert the operator. 2264 Perhaps we can approximate the Jacobian using only the cellwise P-P gradient from Coulomb 2265 */ 2266 static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat J, Mat P, void *ctx) 2267 { 2268 DM sw; 2269 const PetscReal *coords, *vel; 2270 PetscInt dim, d, Np, p, rStart; 2271 2272 PetscFunctionBeginUser; 2273 PetscCall(TSGetDM(ts, &sw)); 2274 PetscCall(DMGetDimension(sw, &dim)); 2275 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 2276 PetscCall(MatGetOwnershipRange(J, &rStart, NULL)); 2277 PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **)&coords)); 2278 PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **)&vel)); 2279 Np /= 2 * dim; 2280 for (p = 0; p < Np; ++p) { 2281 const PetscReal x0 = coords[p * dim + 0]; 2282 const PetscReal vy0 = vel[p * dim + 1]; 2283 const PetscReal omega = vy0 / x0; 2284 PetscScalar vals[4] = {0., 1., -PetscSqr(omega), 0.}; 2285 2286 for (d = 0; d < dim; ++d) { 2287 const PetscInt rows[2] = {(p * 2 + 0) * dim + d + rStart, (p * 2 + 1) * dim + d + rStart}; 2288 PetscCall(MatSetValues(J, 2, rows, 2, rows, vals, INSERT_VALUES)); 2289 } 2290 } 2291 PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **)&coords)); 2292 PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **)&vel)); 2293 PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 2294 PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 2295 PetscFunctionReturn(PETSC_SUCCESS); 2296 } 2297 2298 static PetscErrorCode RHSFunctionX(TS ts, PetscReal t, Vec V, Vec Xres, void *ctx) 2299 { 2300 AppCtx *user = (AppCtx *)ctx; 2301 DM sw; 2302 const PetscScalar *v; 2303 PetscScalar *xres; 2304 PetscInt Np, p, d, dim; 2305 2306 PetscFunctionBeginUser; 2307 PetscCall(PetscLogEventBegin(user->RhsXEvent, ts, 0, 0, 0)); 2308 PetscCall(TSGetDM(ts, &sw)); 2309 PetscCall(DMGetDimension(sw, &dim)); 2310 PetscCall(VecGetLocalSize(Xres, &Np)); 2311 PetscCall(VecGetArrayRead(V, &v)); 2312 PetscCall(VecGetArray(Xres, &xres)); 2313 Np /= dim; 2314 for (p = 0; p < Np; ++p) { 2315 for (d = 0; d < dim; ++d) { 2316 xres[p * dim + d] = v[p * dim + d]; 2317 if (user->fake_1D && d > 0) xres[p * dim + d] = 0; 2318 } 2319 } 2320 PetscCall(VecRestoreArrayRead(V, &v)); 2321 PetscCall(VecRestoreArray(Xres, &xres)); 2322 PetscCall(PetscLogEventEnd(user->RhsXEvent, ts, 0, 0, 0)); 2323 PetscFunctionReturn(PETSC_SUCCESS); 2324 } 2325 2326 static PetscErrorCode RHSFunctionV(TS ts, PetscReal t, Vec X, Vec Vres, void *ctx) 2327 { 2328 DM sw; 2329 AppCtx *user = (AppCtx *)ctx; 2330 SNES snes = ((AppCtx *)ctx)->snes; 2331 const PetscScalar *x; 2332 const PetscReal *coords, *vel; 2333 PetscReal *E, m_p, q_p; 2334 PetscScalar *vres; 2335 PetscInt Np, p, dim, d; 2336 Parameter *param; 2337 2338 PetscFunctionBeginUser; 2339 PetscCall(PetscLogEventBegin(user->RhsVEvent, ts, 0, 0, 0)); 2340 PetscCall(TSGetDM(ts, &sw)); 2341 PetscCall(DMGetDimension(sw, &dim)); 2342 PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **)&coords)); 2343 PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **)&vel)); 2344 PetscCall(DMSwarmGetField(sw, "E_field", NULL, NULL, (void **)&E)); 2345 PetscCall(PetscBagGetData(user->bag, (void **)¶m)); 2346 m_p = user->masses[0] * param->m0; 2347 q_p = user->charges[0] * param->q0; 2348 PetscCall(VecGetLocalSize(Vres, &Np)); 2349 PetscCall(VecGetArrayRead(X, &x)); 2350 PetscCall(VecGetArray(Vres, &vres)); 2351 PetscCall(ComputeFieldAtParticles(snes, sw, E)); 2352 2353 Np /= dim; 2354 for (p = 0; p < Np; ++p) { 2355 for (d = 0; d < dim; ++d) { 2356 vres[p * dim + d] = q_p * E[p * dim + d] / m_p; 2357 if (user->fake_1D && d > 0) vres[p * dim + d] = 0.; 2358 } 2359 } 2360 PetscCall(VecRestoreArrayRead(X, &x)); 2361 /* 2362 Synchronized, ordered output for parallel/sequential test cases. 2363 In the 1D (on the 2D mesh) case, every y component should be zero. 2364 */ 2365 if (user->checkVRes) { 2366 PetscBool pr = user->checkVRes > 1 ? PETSC_TRUE : PETSC_FALSE; 2367 PetscInt step; 2368 2369 PetscCall(TSGetStepNumber(ts, &step)); 2370 if (pr) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "step: %" PetscInt_FMT "\n", step)); 2371 for (PetscInt p = 0; p < Np; ++p) { 2372 if (pr) PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "Residual: %.12g %.12g\n", (double)PetscRealPart(vres[p * dim + 0]), (double)PetscRealPart(vres[p * dim + 1]))); 2373 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])); 2374 } 2375 if (pr) PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 2376 } 2377 PetscCall(VecRestoreArray(Vres, &vres)); 2378 PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **)&coords)); 2379 PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **)&vel)); 2380 PetscCall(DMSwarmRestoreField(sw, "E_field", NULL, NULL, (void **)&E)); 2381 PetscCall(PetscLogEventEnd(user->RhsVEvent, ts, 0, 0, 0)); 2382 PetscFunctionReturn(PETSC_SUCCESS); 2383 } 2384 2385 static PetscErrorCode CreateSolution(TS ts) 2386 { 2387 DM sw; 2388 Vec u; 2389 PetscInt dim, Np; 2390 2391 PetscFunctionBegin; 2392 PetscCall(TSGetDM(ts, &sw)); 2393 PetscCall(DMGetDimension(sw, &dim)); 2394 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 2395 PetscCall(VecCreate(PETSC_COMM_WORLD, &u)); 2396 PetscCall(VecSetBlockSize(u, dim)); 2397 PetscCall(VecSetSizes(u, 2 * Np * dim, PETSC_DECIDE)); 2398 PetscCall(VecSetUp(u)); 2399 PetscCall(TSSetSolution(ts, u)); 2400 PetscCall(VecDestroy(&u)); 2401 PetscFunctionReturn(PETSC_SUCCESS); 2402 } 2403 2404 static PetscErrorCode SetProblem(TS ts) 2405 { 2406 AppCtx *user; 2407 DM sw; 2408 2409 PetscFunctionBegin; 2410 PetscCall(TSGetDM(ts, &sw)); 2411 PetscCall(DMGetApplicationContext(sw, (void **)&user)); 2412 // Define unified system for (X, V) 2413 { 2414 Mat J; 2415 PetscInt dim, Np; 2416 2417 PetscCall(DMGetDimension(sw, &dim)); 2418 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 2419 PetscCall(MatCreate(PETSC_COMM_WORLD, &J)); 2420 PetscCall(MatSetSizes(J, 2 * Np * dim, 2 * Np * dim, PETSC_DECIDE, PETSC_DECIDE)); 2421 PetscCall(MatSetBlockSize(J, 2 * dim)); 2422 PetscCall(MatSetFromOptions(J)); 2423 PetscCall(MatSetUp(J)); 2424 PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, user)); 2425 PetscCall(TSSetRHSJacobian(ts, J, J, RHSJacobian, user)); 2426 PetscCall(MatDestroy(&J)); 2427 } 2428 /* Define split system for X and V */ 2429 { 2430 Vec u; 2431 IS isx, isv, istmp; 2432 const PetscInt *idx; 2433 PetscInt dim, Np, rstart; 2434 2435 PetscCall(TSGetSolution(ts, &u)); 2436 PetscCall(DMGetDimension(sw, &dim)); 2437 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 2438 PetscCall(VecGetOwnershipRange(u, &rstart, NULL)); 2439 PetscCall(ISCreateStride(PETSC_COMM_WORLD, Np, (rstart / dim) + 0, 2, &istmp)); 2440 PetscCall(ISGetIndices(istmp, &idx)); 2441 PetscCall(ISCreateBlock(PETSC_COMM_WORLD, dim, Np, idx, PETSC_COPY_VALUES, &isx)); 2442 PetscCall(ISRestoreIndices(istmp, &idx)); 2443 PetscCall(ISDestroy(&istmp)); 2444 PetscCall(ISCreateStride(PETSC_COMM_WORLD, Np, (rstart / dim) + 1, 2, &istmp)); 2445 PetscCall(ISGetIndices(istmp, &idx)); 2446 PetscCall(ISCreateBlock(PETSC_COMM_WORLD, dim, Np, idx, PETSC_COPY_VALUES, &isv)); 2447 PetscCall(ISRestoreIndices(istmp, &idx)); 2448 PetscCall(ISDestroy(&istmp)); 2449 PetscCall(TSRHSSplitSetIS(ts, "position", isx)); 2450 PetscCall(TSRHSSplitSetIS(ts, "momentum", isv)); 2451 PetscCall(ISDestroy(&isx)); 2452 PetscCall(ISDestroy(&isv)); 2453 PetscCall(TSRHSSplitSetRHSFunction(ts, "position", NULL, RHSFunctionX, user)); 2454 PetscCall(TSRHSSplitSetRHSFunction(ts, "momentum", NULL, RHSFunctionV, user)); 2455 } 2456 PetscFunctionReturn(PETSC_SUCCESS); 2457 } 2458 2459 static PetscErrorCode DMSwarmTSRedistribute(TS ts) 2460 { 2461 DM sw; 2462 Vec u; 2463 PetscReal t, maxt, dt; 2464 PetscInt n, maxn; 2465 2466 PetscFunctionBegin; 2467 PetscCall(TSGetDM(ts, &sw)); 2468 PetscCall(TSGetTime(ts, &t)); 2469 PetscCall(TSGetMaxTime(ts, &maxt)); 2470 PetscCall(TSGetTimeStep(ts, &dt)); 2471 PetscCall(TSGetStepNumber(ts, &n)); 2472 PetscCall(TSGetMaxSteps(ts, &maxn)); 2473 2474 PetscCall(TSReset(ts)); 2475 PetscCall(TSSetDM(ts, sw)); 2476 PetscCall(TSSetFromOptions(ts)); 2477 PetscCall(TSSetTime(ts, t)); 2478 PetscCall(TSSetMaxTime(ts, maxt)); 2479 PetscCall(TSSetTimeStep(ts, dt)); 2480 PetscCall(TSSetStepNumber(ts, n)); 2481 PetscCall(TSSetMaxSteps(ts, maxn)); 2482 2483 PetscCall(CreateSolution(ts)); 2484 PetscCall(SetProblem(ts)); 2485 PetscCall(TSGetSolution(ts, &u)); 2486 PetscFunctionReturn(PETSC_SUCCESS); 2487 } 2488 2489 PetscErrorCode line(PetscInt dim, PetscReal time, const PetscReal dummy[], PetscInt p, PetscScalar x[], void *ctx) 2490 { 2491 DM sw, cdm; 2492 PetscInt Np; 2493 PetscReal low[2], high[2]; 2494 AppCtx *user = (AppCtx *)ctx; 2495 2496 sw = user->swarm; 2497 PetscCall(DMSwarmGetCellDM(sw, &cdm)); 2498 // Get the bounding box so we can equally space the particles 2499 PetscCall(DMGetLocalBoundingBox(cdm, low, high)); 2500 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 2501 // shift it by h/2 so nothing is initialized directly on a boundary 2502 x[0] = ((high[0] - low[0]) / Np) * (p + 0.5); 2503 x[1] = 0.; 2504 return PETSC_SUCCESS; 2505 } 2506 2507 /* 2508 InitializeSolveAndSwarm - Set the solution values to the swarm coordinates and velocities, and also possibly set the initial values. 2509 2510 Input Parameters: 2511 + ts - The TS 2512 - useInitial - Flag to also set the initial conditions to the current coordinates and velocities and setup the problem 2513 2514 Output Parameters: 2515 . u - The initialized solution vector 2516 2517 Level: advanced 2518 2519 .seealso: InitializeSolve() 2520 */ 2521 static PetscErrorCode InitializeSolveAndSwarm(TS ts, PetscBool useInitial) 2522 { 2523 DM sw; 2524 Vec u, gc, gv, gc0, gv0; 2525 IS isx, isv; 2526 PetscInt dim; 2527 AppCtx *user; 2528 2529 PetscFunctionBeginUser; 2530 PetscCall(TSGetDM(ts, &sw)); 2531 PetscCall(DMGetApplicationContext(sw, &user)); 2532 PetscCall(DMGetDimension(sw, &dim)); 2533 if (useInitial) { 2534 PetscReal v0[2] = {1., 0.}; 2535 if (user->perturbed_weights) { 2536 PetscCall(InitializeParticles_PerturbedWeights(sw, user)); 2537 } else { 2538 PetscCall(DMSwarmComputeLocalSizeFromOptions(sw)); 2539 PetscCall(DMSwarmInitializeCoordinates(sw)); 2540 if (user->fake_1D) { 2541 PetscCall(InitializeVelocities_Fake1D(sw, user)); 2542 } else { 2543 PetscCall(DMSwarmInitializeVelocitiesFromOptions(sw, v0)); 2544 } 2545 } 2546 PetscCall(DMSwarmMigrate(sw, PETSC_TRUE)); 2547 PetscCall(DMSwarmTSRedistribute(ts)); 2548 } 2549 PetscCall(TSGetSolution(ts, &u)); 2550 PetscCall(TSRHSSplitGetIS(ts, "position", &isx)); 2551 PetscCall(TSRHSSplitGetIS(ts, "momentum", &isv)); 2552 PetscCall(DMSwarmCreateGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc)); 2553 PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "initCoordinates", &gc0)); 2554 PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "velocity", &gv)); 2555 PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "initVelocity", &gv0)); 2556 if (useInitial) { 2557 PetscCall(VecCopy(gc, gc0)); 2558 PetscCall(VecCopy(gv, gv0)); 2559 } 2560 PetscCall(VecISCopy(u, isx, SCATTER_FORWARD, gc)); 2561 PetscCall(VecISCopy(u, isv, SCATTER_FORWARD, gv)); 2562 PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc)); 2563 PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "initCoordinates", &gc0)); 2564 PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "velocity", &gv)); 2565 PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "initVelocity", &gv0)); 2566 PetscFunctionReturn(PETSC_SUCCESS); 2567 } 2568 2569 static PetscErrorCode InitializeSolve(TS ts, Vec u) 2570 { 2571 PetscFunctionBegin; 2572 PetscCall(TSSetSolution(ts, u)); 2573 PetscCall(InitializeSolveAndSwarm(ts, PETSC_TRUE)); 2574 PetscFunctionReturn(PETSC_SUCCESS); 2575 } 2576 2577 static PetscErrorCode ComputeError(TS ts, Vec U, Vec E) 2578 { 2579 MPI_Comm comm; 2580 DM sw; 2581 AppCtx *user; 2582 const PetscScalar *u; 2583 const PetscReal *coords, *vel; 2584 PetscScalar *e; 2585 PetscReal t; 2586 PetscInt dim, Np, p; 2587 2588 PetscFunctionBeginUser; 2589 PetscCall(PetscObjectGetComm((PetscObject)ts, &comm)); 2590 PetscCall(TSGetDM(ts, &sw)); 2591 PetscCall(DMGetApplicationContext(sw, &user)); 2592 PetscCall(DMGetDimension(sw, &dim)); 2593 PetscCall(TSGetSolveTime(ts, &t)); 2594 PetscCall(VecGetArray(E, &e)); 2595 PetscCall(VecGetArrayRead(U, &u)); 2596 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 2597 PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **)&coords)); 2598 PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **)&vel)); 2599 Np /= 2 * dim; 2600 for (p = 0; p < Np; ++p) { 2601 /* TODO generalize initial conditions and project into plane instead of assuming x-y */ 2602 const PetscReal r0 = DMPlex_NormD_Internal(dim, &coords[p * dim]); 2603 const PetscReal th0 = PetscAtan2Real(coords[p * dim + 1], coords[p * dim + 0]); 2604 const PetscReal v0 = DMPlex_NormD_Internal(dim, &vel[p * dim]); 2605 const PetscReal omega = v0 / r0; 2606 const PetscReal ct = PetscCosReal(omega * t + th0); 2607 const PetscReal st = PetscSinReal(omega * t + th0); 2608 const PetscScalar *x = &u[(p * 2 + 0) * dim]; 2609 const PetscScalar *v = &u[(p * 2 + 1) * dim]; 2610 const PetscReal xe[3] = {r0 * ct, r0 * st, 0.0}; 2611 const PetscReal ve[3] = {-v0 * st, v0 * ct, 0.0}; 2612 PetscInt d; 2613 2614 for (d = 0; d < dim; ++d) { 2615 e[(p * 2 + 0) * dim + d] = x[d] - xe[d]; 2616 e[(p * 2 + 1) * dim + d] = v[d] - ve[d]; 2617 } 2618 if (user->error) { 2619 const PetscReal en = 0.5 * DMPlex_DotRealD_Internal(dim, v, v); 2620 const PetscReal exen = 0.5 * PetscSqr(v0); 2621 PetscCall(PetscPrintf(comm, "t %.4g: p%" PetscInt_FMT " error [%.2g %.2g] sol [(%.6lf %.6lf) (%.6lf %.6lf)] exact [(%.6lf %.6lf) (%.6lf %.6lf)] energy/exact energy %g / %g (%.10lf%%)\n", (double)t, p, (double)DMPlex_NormD_Internal(dim, &e[(p * 2 + 0) * dim]), (double)DMPlex_NormD_Internal(dim, &e[(p * 2 + 1) * dim]), (double)x[0], (double)x[1], (double)v[0], (double)v[1], (double)xe[0], (double)xe[1], (double)ve[0], (double)ve[1], (double)en, (double)exen, (double)(PetscAbsReal(exen - en) * 100. / exen))); 2622 } 2623 } 2624 PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **)&coords)); 2625 PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **)&vel)); 2626 PetscCall(VecRestoreArrayRead(U, &u)); 2627 PetscCall(VecRestoreArray(E, &e)); 2628 PetscFunctionReturn(PETSC_SUCCESS); 2629 } 2630 2631 static PetscErrorCode MigrateParticles(TS ts) 2632 { 2633 DM sw, cdm; 2634 const PetscReal *L; 2635 AppCtx *ctx; 2636 2637 PetscFunctionBeginUser; 2638 PetscCall(TSGetDM(ts, &sw)); 2639 PetscCall(DMGetApplicationContext(sw, &ctx)); 2640 PetscCall(DMViewFromOptions(sw, NULL, "-migrate_view_pre")); 2641 { 2642 Vec u, gc, gv, position, momentum; 2643 IS isx, isv; 2644 PetscReal *pos, *mom; 2645 2646 PetscCall(TSGetSolution(ts, &u)); 2647 PetscCall(TSRHSSplitGetIS(ts, "position", &isx)); 2648 PetscCall(TSRHSSplitGetIS(ts, "momentum", &isv)); 2649 PetscCall(VecGetSubVector(u, isx, &position)); 2650 PetscCall(VecGetSubVector(u, isv, &momentum)); 2651 PetscCall(VecGetArray(position, &pos)); 2652 PetscCall(VecGetArray(momentum, &mom)); 2653 PetscCall(DMSwarmCreateGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc)); 2654 PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "velocity", &gv)); 2655 PetscCall(VecISCopy(u, isx, SCATTER_REVERSE, gc)); 2656 PetscCall(VecISCopy(u, isv, SCATTER_REVERSE, gv)); 2657 2658 PetscCall(DMSwarmGetCellDM(sw, &cdm)); 2659 PetscCall(DMGetPeriodicity(cdm, NULL, NULL, &L)); 2660 if ((L[0] || L[1]) >= 0.) { 2661 PetscReal *x, *v, upper[3], lower[3]; 2662 PetscInt Np, dim; 2663 2664 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 2665 PetscCall(DMGetDimension(cdm, &dim)); 2666 PetscCall(DMGetBoundingBox(cdm, lower, upper)); 2667 PetscCall(VecGetArray(gc, &x)); 2668 PetscCall(VecGetArray(gv, &v)); 2669 for (PetscInt p = 0; p < Np; ++p) { 2670 for (PetscInt d = 0; d < dim; ++d) { 2671 if (pos[p * dim + d] < lower[d]) { 2672 x[p * dim + d] = pos[p * dim + d] + (upper[d] - lower[d]); 2673 } else if (pos[p * dim + d] > upper[d]) { 2674 x[p * dim + d] = pos[p * dim + d] - (upper[d] - lower[d]); 2675 } else { 2676 x[p * dim + d] = pos[p * dim + d]; 2677 } 2678 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]); 2679 v[p * dim + d] = mom[p * dim + d]; 2680 } 2681 } 2682 PetscCall(VecRestoreArray(gc, &x)); 2683 PetscCall(VecRestoreArray(gv, &v)); 2684 } 2685 PetscCall(VecRestoreArray(position, &pos)); 2686 PetscCall(VecRestoreArray(momentum, &mom)); 2687 PetscCall(VecRestoreSubVector(u, isx, &position)); 2688 PetscCall(VecRestoreSubVector(u, isv, &momentum)); 2689 PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "velocity", &gv)); 2690 PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc)); 2691 } 2692 PetscCall(DMSwarmMigrate(sw, PETSC_TRUE)); 2693 PetscInt step; 2694 2695 PetscCall(TSGetStepNumber(ts, &step)); 2696 if (ctx->remap && !(step % ctx->remapFreq)) { 2697 // Monitor electric field before we destroy it 2698 PetscReal ptime; 2699 PetscInt step; 2700 2701 PetscCall(TSGetStepNumber(ts, &step)); 2702 PetscCall(TSGetTime(ts, &ptime)); 2703 if (ctx->efield_monitor) PetscCall(MonitorEField(ts, step, ptime, NULL, ctx)); 2704 if (ctx->poisson_monitor) PetscCall(MonitorPoisson(ts, step, ptime, NULL, ctx)); 2705 PetscCall(DMSwarmRemap(sw)); 2706 ctx->validE = PETSC_FALSE; 2707 } 2708 // This MUST come last, since it recreates the subswarms and they must DMClone() the new swarm 2709 PetscCall(DMSwarmTSRedistribute(ts)); 2710 PetscCall(InitializeSolveAndSwarm(ts, PETSC_FALSE)); 2711 PetscFunctionReturn(PETSC_SUCCESS); 2712 } 2713 2714 int main(int argc, char **argv) 2715 { 2716 DM dm, sw; 2717 TS ts; 2718 Vec u; 2719 PetscReal dt; 2720 PetscInt maxn; 2721 AppCtx user; 2722 2723 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 2724 PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 2725 PetscCall(PetscBagCreate(PETSC_COMM_SELF, sizeof(Parameter), &user.bag)); 2726 PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 2727 PetscCall(CreatePoisson(dm, &user)); 2728 PetscCall(CreateSwarm(dm, &user, &sw)); 2729 PetscCall(SetupParameters(PETSC_COMM_WORLD, &user)); 2730 PetscCall(InitializeConstants(sw, &user)); 2731 PetscCall(DMSetApplicationContext(sw, &user)); 2732 2733 PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 2734 PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 2735 PetscCall(TSSetDM(ts, sw)); 2736 PetscCall(TSSetMaxTime(ts, 0.1)); 2737 PetscCall(TSSetTimeStep(ts, 0.00001)); 2738 PetscCall(TSSetMaxSteps(ts, 100)); 2739 PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 2740 2741 if (user.efield_monitor) PetscCall(TSMonitorSet(ts, MonitorEField, &user, NULL)); 2742 if (user.moment_monitor) PetscCall(TSMonitorSet(ts, MonitorMoments, &user, NULL)); 2743 if (user.initial_monitor) PetscCall(TSMonitorSet(ts, MonitorInitialConditions, &user, NULL)); 2744 if (user.positions_monitor) PetscCall(TSMonitorSet(ts, MonitorPositions_2D, &user, NULL)); 2745 if (user.poisson_monitor) PetscCall(TSMonitorSet(ts, MonitorPoisson, &user, NULL)); 2746 2747 PetscCall(TSSetFromOptions(ts)); 2748 PetscCall(TSGetTimeStep(ts, &dt)); 2749 PetscCall(TSGetMaxSteps(ts, &maxn)); 2750 user.steps = maxn; 2751 user.stepSize = dt; 2752 PetscCall(SetupContext(dm, sw, &user)); 2753 PetscCall(TSSetComputeInitialCondition(ts, InitializeSolve)); 2754 PetscCall(TSSetComputeExactError(ts, ComputeError)); 2755 PetscCall(TSSetPostStep(ts, MigrateParticles)); 2756 PetscCall(CreateSolution(ts)); 2757 PetscCall(TSGetSolution(ts, &u)); 2758 PetscCall(TSComputeInitialCondition(ts, u)); 2759 PetscCall(CheckNonNegativeWeights(sw, &user)); 2760 PetscCall(TSSolve(ts, NULL)); 2761 2762 PetscCall(SNESDestroy(&user.snes)); 2763 PetscCall(MatDestroy(&user.M)); 2764 PetscCall(PetscFEGeomDestroy(&user.fegeom)); 2765 PetscCall(TSDestroy(&ts)); 2766 PetscCall(DMDestroy(&sw)); 2767 PetscCall(DMDestroy(&dm)); 2768 PetscCall(DestroyContext(&user)); 2769 PetscCall(PetscFinalize()); 2770 return 0; 2771 } 2772 2773 /*TEST 2774 2775 build: 2776 requires: !complex double 2777 2778 # This tests that we can put particles in a box and compute the Coulomb force 2779 # Recommend -draw_size 500,500 2780 testset: 2781 requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT) 2782 args: -dm_plex_dim 2 -fake_1D -dm_plex_simplex 0 -dm_plex_box_faces 20,1 \ 2783 -dm_plex_box_lower 0,-1 -dm_plex_box_upper 12.5664,1 \ 2784 -dm_plex_box_bd periodic,none \ 2785 -vdm_plex_simplex 0 \ 2786 -dm_swarm_coordinate_density constant -dm_swarm_num_particles 100 \ 2787 -sigma 1.0e-8 -timeScale 2.0e-14 \ 2788 -ts_type basicsymplectic -ts_basicsymplectic_type 1 \ 2789 -ts_monitor_sp_swarm -ts_monitor_sp_swarm_retain 0 -ts_monitor_sp_swarm_phase 0 \ 2790 -output_step 50 -check_vel_res 2791 test: 2792 suffix: none_1d 2793 requires: 2794 args: -em_type none -error 2795 test: 2796 suffix: coulomb_1d 2797 args: -em_type coulomb 2798 2799 # for viewers 2800 #-ts_monitor_sp_swarm_phase -ts_monitor_sp_swarm -em_snes_monitor -ts_monitor_sp_swarm_multi_species 0 -ts_monitor_sp_swarm_retain 0 2801 testset: 2802 nsize: {{1 2}} 2803 requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT) 2804 args: -dm_plex_dim 2 -fake_1D -dm_plex_simplex 0 -dm_plex_box_faces 36,1 \ 2805 -dm_plex_box_lower 0.,-0.5 -dm_plex_box_upper 12.5664,0.5 \ 2806 -dm_plex_box_bd periodic,none \ 2807 -vdm_plex_dim 1 -vdm_plex_simplex 0 -vdm_plex_box_faces 10 \ 2808 -vdm_plex_box_lower -3 -vdm_plex_box_upper 3 \ 2809 -dm_swarm_num_species 1 -twostream -charges -1.,1. -sigma 1.0e-8 \ 2810 -cosine_coefficients 0.01,0.5 -perturbed_weights -total_weight 1. \ 2811 -ts_type basicsymplectic -ts_basicsymplectic_type 2 \ 2812 -ts_dt 0.01 -ts_max_time 5 -ts_max_steps 10 \ 2813 -em_snes_atol 1.e-15 -em_snes_error_if_not_converged -em_ksp_error_if_not_converged \ 2814 -output_step 1 -check_vel_res -dm_swarm_print_coords 1 -dm_swarm_print_weights 1 2815 test: 2816 suffix: two_stream_c0 2817 args: -em_type primal -petscfe_default_quadrature_order 2 -petscspace_degree 2 -em_pc_type svd 2818 test: 2819 suffix: two_stream_rt 2820 requires: superlu_dist 2821 args: -em_type mixed \ 2822 -potential_petscspace_degree 0 \ 2823 -potential_petscdualspace_lagrange_use_moments \ 2824 -potential_petscdualspace_lagrange_moment_order 2 \ 2825 -field_petscspace_degree 2 -field_petscfe_default_quadrature_order 1 \ 2826 -field_petscspace_type sum \ 2827 -field_petscspace_variables 2 \ 2828 -field_petscspace_components 2 \ 2829 -field_petscspace_sum_spaces 2 \ 2830 -field_petscspace_sum_concatenate true \ 2831 -field_sumcomp_0_petscspace_variables 2 \ 2832 -field_sumcomp_0_petscspace_type tensor \ 2833 -field_sumcomp_0_petscspace_tensor_spaces 2 \ 2834 -field_sumcomp_0_petscspace_tensor_uniform false \ 2835 -field_sumcomp_0_tensorcomp_0_petscspace_degree 1 \ 2836 -field_sumcomp_0_tensorcomp_1_petscspace_degree 0 \ 2837 -field_sumcomp_1_petscspace_variables 2 \ 2838 -field_sumcomp_1_petscspace_type tensor \ 2839 -field_sumcomp_1_petscspace_tensor_spaces 2 \ 2840 -field_sumcomp_1_petscspace_tensor_uniform false \ 2841 -field_sumcomp_1_tensorcomp_0_petscspace_degree 0 \ 2842 -field_sumcomp_1_tensorcomp_1_petscspace_degree 1 \ 2843 -field_petscdualspace_form_degree -1 \ 2844 -field_petscdualspace_order 1 \ 2845 -field_petscdualspace_lagrange_trimmed true \ 2846 -em_snes_error_if_not_converged \ 2847 -em_ksp_type preonly -em_ksp_error_if_not_converged \ 2848 -em_pc_type fieldsplit -em_pc_fieldsplit_type schur \ 2849 -em_pc_fieldsplit_schur_fact_type full -em_pc_fieldsplit_schur_precondition full \ 2850 -em_fieldsplit_field_pc_type lu \ 2851 -em_fieldsplit_field_pc_factor_mat_solver_type superlu_dist \ 2852 -em_fieldsplit_potential_pc_type svd 2853 2854 # For an eyeball check, we use 2855 # -ts_max_steps 1000 -dm_plex_box_faces 10,1 -vdm_plex_box_faces 2000 -efield_monitor 2856 # For verification, we use 2857 # -ts_max_steps 1000 -dm_plex_box_faces 100,1 -vdm_plex_box_faces 8000 -dm_swarm_num_particles 800000 -monitor_efield 2858 # -ts_monitor_sp_swarm_multi_species 0 -ts_monitor_sp_swarm_retain 0 -ts_monitor_sp_swarm_phase 1 -draw_size 500,500 2859 testset: 2860 nsize: {{1 2}} 2861 requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT) 2862 args: -dm_plex_dim 2 -fake_1D -dm_plex_simplex 0 -dm_plex_box_faces 10,1 \ 2863 -dm_plex_box_lower 0.,-0.5 -dm_plex_box_upper 12.5664,0.5 \ 2864 -dm_plex_box_bd periodic,none \ 2865 -vdm_plex_dim 1 -vdm_plex_simplex 0 -vdm_plex_box_faces 10 \ 2866 -vdm_plex_box_lower -10 -vdm_plex_box_upper 10 \ 2867 -vpetscspace_degree 2 -vdm_plex_hash_location \ 2868 -dm_swarm_num_species 1 -charges -1.,1. \ 2869 -cosine_coefficients 0.01,0.5 -perturbed_weights -total_weight 1. \ 2870 -ts_type basicsymplectic -ts_basicsymplectic_type 1 \ 2871 -ts_dt 0.03 -ts_max_time 500 -ts_max_steps 1 \ 2872 -em_snes_atol 1.e-12 -em_snes_error_if_not_converged -em_ksp_error_if_not_converged \ 2873 -output_step 1 -check_vel_res -dm_swarm_print_coords 1 -dm_swarm_print_weights 1 2874 2875 test: 2876 suffix: uniform_equilibrium_1d 2877 args: -cosine_coefficients 0.0,0.5 -em_type primal -petscspace_degree 1 -em_pc_type svd 2878 test: 2879 suffix: uniform_equilibrium_1d_real 2880 args: -dm_plex_dim 1 -dm_plex_simplex 1 -fake_1D 0 -dm_plex_box_faces 10 \ 2881 -dm_plex_box_lower 0. -dm_plex_box_upper 12.5664 -dm_plex_box_bd periodic \ 2882 -cosine_coefficients 0.0 -em_type primal -petscspace_degree 1 -em_pc_type svd 2883 test: 2884 suffix: uniform_primal_1d 2885 args: -em_type primal -petscspace_degree 1 -em_pc_type svd 2886 test: 2887 suffix: uniform_primal_1d_real 2888 args: -dm_plex_dim 1 -dm_plex_simplex 1 -fake_1D 0 -dm_plex_box_faces 10 \ 2889 -dm_plex_box_lower 0. -dm_plex_box_upper 12.5664 -dm_plex_box_bd periodic \ 2890 -cosine_coefficients 0.01 -em_type primal -petscspace_degree 1 -em_pc_type svd 2891 # NOT WORKING -ftop_pc_type bjacobi -ftop_sub_pc_type lu -ftop_sub_pc_factor_shift_type nonzero 2892 test: 2893 suffix: uniform_primal_1d_real_pfak 2894 nsize: 1 2895 args: -dm_plex_dim 1 -dm_plex_simplex 1 -fake_1D 0 -dm_plex_box_faces 10 \ 2896 -dm_plex_box_lower 0. -dm_plex_box_upper 12.5664 -dm_plex_box_bd periodic \ 2897 -remap_dm_plex_dim 2 -remap_dm_plex_simplex 0 -remap_dm_plex_box_faces 10,10 -remap_dm_plex_box_bd periodic,none \ 2898 -remap_dm_plex_box_lower 0.,-10. -remap_dm_plex_box_upper 12.5664,10. \ 2899 -remap_petscspace_degree 2 -remap_dm_plex_hash_location \ 2900 -remap -remap_freq 1 -remap_type pfak \ 2901 -ftop_ksp_type lsqr -ftop_pc_type none -ftop_ksp_rtol 1.e-14 \ 2902 -ptof_pc_type lu \ 2903 -cosine_coefficients 0.01 -em_type primal -petscspace_degree 1 -em_pc_type svd -em_proj_pc_type lu 2904 test: 2905 requires: superlu_dist 2906 suffix: uniform_mixed_1d 2907 args: -em_type mixed \ 2908 -potential_petscspace_degree 0 \ 2909 -potential_petscdualspace_lagrange_use_moments \ 2910 -potential_petscdualspace_lagrange_moment_order 2 \ 2911 -field_petscspace_degree 2 -field_petscfe_default_quadrature_order 1 \ 2912 -field_petscspace_type sum \ 2913 -field_petscspace_variables 2 \ 2914 -field_petscspace_components 2 \ 2915 -field_petscspace_sum_spaces 2 \ 2916 -field_petscspace_sum_concatenate true \ 2917 -field_sumcomp_0_petscspace_variables 2 \ 2918 -field_sumcomp_0_petscspace_type tensor \ 2919 -field_sumcomp_0_petscspace_tensor_spaces 2 \ 2920 -field_sumcomp_0_petscspace_tensor_uniform false \ 2921 -field_sumcomp_0_tensorcomp_0_petscspace_degree 1 \ 2922 -field_sumcomp_0_tensorcomp_1_petscspace_degree 0 \ 2923 -field_sumcomp_1_petscspace_variables 2 \ 2924 -field_sumcomp_1_petscspace_type tensor \ 2925 -field_sumcomp_1_petscspace_tensor_spaces 2 \ 2926 -field_sumcomp_1_petscspace_tensor_uniform false \ 2927 -field_sumcomp_1_tensorcomp_0_petscspace_degree 0 \ 2928 -field_sumcomp_1_tensorcomp_1_petscspace_degree 1 \ 2929 -field_petscdualspace_form_degree -1 \ 2930 -field_petscdualspace_order 1 \ 2931 -field_petscdualspace_lagrange_trimmed true \ 2932 -em_snes_error_if_not_converged \ 2933 -em_ksp_type preonly -em_ksp_error_if_not_converged \ 2934 -em_pc_type fieldsplit -em_pc_fieldsplit_type schur \ 2935 -em_pc_fieldsplit_schur_fact_type full -em_pc_fieldsplit_schur_precondition full \ 2936 -em_fieldsplit_field_pc_type lu \ 2937 -em_fieldsplit_field_pc_factor_mat_solver_type superlu_dist \ 2938 -em_fieldsplit_potential_pc_type svd 2939 2940 TEST*/ 2941