1d0c080abSJoseph Pusztay static char help[] = "Example application of the Bhatnagar-Gross-Krook (BGK) collision operator.\n\ 2d0c080abSJoseph Pusztay This example is a 0D-1V setting for the kinetic equation\n\ 3d0c080abSJoseph Pusztay https://en.wikipedia.org/wiki/Bhatnagar%E2%80%93Gross%E2%80%93Krook_operator\n"; 4d0c080abSJoseph Pusztay 5d0c080abSJoseph Pusztay #include <petscdmplex.h> 6d0c080abSJoseph Pusztay #include <petscdmswarm.h> 7d0c080abSJoseph Pusztay #include <petscts.h> 8d0c080abSJoseph Pusztay #include <petscdraw.h> 9d0c080abSJoseph Pusztay #include <petscviewer.h> 10d0c080abSJoseph Pusztay 11d0c080abSJoseph Pusztay typedef struct { 12d0c080abSJoseph Pusztay PetscInt particlesPerCell; /* The number of partices per cell */ 13d0c080abSJoseph Pusztay PetscReal momentTol; /* Tolerance for checking moment conservation */ 14d0c080abSJoseph Pusztay PetscBool monitorhg; /* Flag for using the TS histogram monitor */ 15d0c080abSJoseph Pusztay PetscBool monitorsp; /* Flag for using the TS scatter monitor */ 16d0c080abSJoseph Pusztay PetscBool monitorks; /* Monitor to perform KS test to the maxwellian */ 17d0c080abSJoseph Pusztay PetscBool error; /* Flag for printing the error */ 18d0c080abSJoseph Pusztay PetscInt ostep; /* print the energy at each ostep time steps */ 19d0c080abSJoseph Pusztay PetscDraw draw; /* The draw object for histogram monitoring */ 20d0c080abSJoseph Pusztay PetscDrawHG drawhg; /* The histogram draw context for monitoring */ 21d0c080abSJoseph Pusztay PetscDrawSP drawsp; /* The scatter plot draw context for the monitor */ 22d0c080abSJoseph Pusztay PetscDrawSP drawks; /* Scatterplot draw context for KS test */ 23d0c080abSJoseph Pusztay } AppCtx; 24d0c080abSJoseph Pusztay 25*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) { 26d0c080abSJoseph Pusztay PetscFunctionBeginUser; 27d0c080abSJoseph Pusztay options->monitorhg = PETSC_FALSE; 28d0c080abSJoseph Pusztay options->monitorsp = PETSC_FALSE; 29d0c080abSJoseph Pusztay options->monitorks = PETSC_FALSE; 30d0c080abSJoseph Pusztay options->particlesPerCell = 1; 31d0c080abSJoseph Pusztay options->momentTol = 100.0 * PETSC_MACHINE_EPSILON; 32d0c080abSJoseph Pusztay options->ostep = 100; 33d0c080abSJoseph Pusztay 34d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Collision Options", "DMPLEX"); 359566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-monitorhg", "Flag to use the TS histogram monitor", "ex28.c", options->monitorhg, &options->monitorhg, NULL)); 369566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-monitorsp", "Flag to use the TS scatter plot monitor", "ex28.c", options->monitorsp, &options->monitorsp, NULL)); 379566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-monitorks", "Flag to plot KS test results", "ex28.c", options->monitorks, &options->monitorks, NULL)); 389566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-particles_per_cell", "Number of particles per cell", "ex28.c", options->particlesPerCell, &options->particlesPerCell, NULL)); 399566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-output_step", "Number of time steps between output", "ex28.c", options->ostep, &options->ostep, PETSC_NULL)); 40d0609cedSBarry Smith PetscOptionsEnd(); 41d0c080abSJoseph Pusztay 42d0c080abSJoseph Pusztay PetscFunctionReturn(0); 43d0c080abSJoseph Pusztay } 44d0c080abSJoseph Pusztay 45d0c080abSJoseph Pusztay /* Create the mesh for velocity space */ 46*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user) { 47d0c080abSJoseph Pusztay PetscFunctionBeginUser; 489566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 499566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 509566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 519566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 52d0c080abSJoseph Pusztay PetscFunctionReturn(0); 53d0c080abSJoseph Pusztay } 54d0c080abSJoseph Pusztay 55d0c080abSJoseph Pusztay /* Since we are putting the same number of particles in each cell, this amounts to a uniform distribution of v */ 56*9371c9d4SSatish Balay static PetscErrorCode SetInitialCoordinates(DM sw) { 57d0c080abSJoseph Pusztay AppCtx *user; 58d0c080abSJoseph Pusztay PetscRandom rnd; 59d0c080abSJoseph Pusztay DM dm; 60d0c080abSJoseph Pusztay DMPolytopeType ct; 61d0c080abSJoseph Pusztay PetscBool simplex; 62d0c080abSJoseph Pusztay PetscReal *centroid, *coords, *xi0, *v0, *J, *invJ, detJ, *vals; 63d0c080abSJoseph Pusztay PetscInt dim, d, cStart, cEnd, c, Np, p; 64d0c080abSJoseph Pusztay 65d0c080abSJoseph Pusztay PetscFunctionBeginUser; 669566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)sw), &rnd)); 679566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0)); 689566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rnd)); 69d0c080abSJoseph Pusztay 709566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(sw, &user)); 71d0c080abSJoseph Pusztay Np = user->particlesPerCell; 729566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 739566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(sw, &dm)); 749566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 759566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 76d0c080abSJoseph Pusztay simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; 779566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(dim, ¢roid, dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ)); 78d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) xi0[d] = -1.0; 799566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 809566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&vals)); 81d0c080abSJoseph Pusztay for (c = cStart; c < cEnd; ++c) { 82d0c080abSJoseph Pusztay if (Np == 1) { 839566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL)); 84d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) { 85d0c080abSJoseph Pusztay coords[c * dim + d] = centroid[d]; 86d0c080abSJoseph Pusztay if ((coords[c * dim + d] >= -1) && (coords[c * dim + d] <= 1)) { 87d0c080abSJoseph Pusztay vals[c] = 1.0; 88d0c080abSJoseph Pusztay } else { 89d0c080abSJoseph Pusztay vals[c] = 0.; 90d0c080abSJoseph Pusztay } 91d0c080abSJoseph Pusztay } 92d0c080abSJoseph Pusztay } else { 939566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); /* affine */ 94d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 95d0c080abSJoseph Pusztay const PetscInt n = c * Np + p; 96d0c080abSJoseph Pusztay PetscReal sum = 0.0, refcoords[3]; 97d0c080abSJoseph Pusztay 98d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) { 999566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(rnd, &refcoords[d])); 100d0c080abSJoseph Pusztay sum += refcoords[d]; 101d0c080abSJoseph Pusztay } 102*9371c9d4SSatish Balay if (simplex && sum > 0.0) 103*9371c9d4SSatish Balay for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim) * sum; 104d0c080abSJoseph Pusztay vals[n] = 1.0; 1059566063dSJacob Faibussowitsch PetscCall(DMPlexReferenceToCoordinates(dm, c, 1, refcoords, &coords[n * dim])); 106d0c080abSJoseph Pusztay } 107d0c080abSJoseph Pusztay } 108d0c080abSJoseph Pusztay } 1099566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 1109566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&vals)); 1119566063dSJacob Faibussowitsch PetscCall(PetscFree5(centroid, xi0, v0, J, invJ)); 1129566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rnd)); 113d0c080abSJoseph Pusztay PetscFunctionReturn(0); 114d0c080abSJoseph Pusztay } 115d0c080abSJoseph Pusztay 116d0c080abSJoseph Pusztay /* The intiial conditions are just the initial particle weights */ 117*9371c9d4SSatish Balay static PetscErrorCode SetInitialConditions(DM dmSw, Vec u) { 118d0c080abSJoseph Pusztay DM dm; 119d0c080abSJoseph Pusztay AppCtx *user; 120d0c080abSJoseph Pusztay PetscReal *vals; 121d0c080abSJoseph Pusztay PetscScalar *initialConditions; 122d0c080abSJoseph Pusztay PetscInt dim, d, cStart, cEnd, c, Np, p, n; 123d0c080abSJoseph Pusztay 124d0c080abSJoseph Pusztay PetscFunctionBeginUser; 1259566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &n)); 1269566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(dmSw, &user)); 127d0c080abSJoseph Pusztay Np = user->particlesPerCell; 1289566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dmSw, &dm)); 1299566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1309566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 13163a3b9bcSJacob Faibussowitsch PetscCheck(n == (cEnd - cStart) * Np, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "TS solution local size %" PetscInt_FMT " != %" PetscInt_FMT " nm particles", n, (cEnd - cStart) * Np); 1329566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dmSw, "w_q", NULL, NULL, (void **)&vals)); 1339566063dSJacob Faibussowitsch PetscCall(VecGetArray(u, &initialConditions)); 134d0c080abSJoseph Pusztay for (c = cStart; c < cEnd; ++c) { 135d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 136d0c080abSJoseph Pusztay const PetscInt n = c * Np + p; 137*9371c9d4SSatish Balay for (d = 0; d < dim; d++) { initialConditions[n] = vals[n]; } 138d0c080abSJoseph Pusztay } 139d0c080abSJoseph Pusztay } 1409566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(u, &initialConditions)); 1419566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dmSw, "w_q", NULL, NULL, (void **)&vals)); 142d0c080abSJoseph Pusztay PetscFunctionReturn(0); 143d0c080abSJoseph Pusztay } 144d0c080abSJoseph Pusztay 145*9371c9d4SSatish Balay static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user) { 146d0c080abSJoseph Pusztay PetscInt *cellid; 147d0c080abSJoseph Pusztay PetscInt dim, cStart, cEnd, c, Np = user->particlesPerCell, p; 148d0c080abSJoseph Pusztay 149d0c080abSJoseph Pusztay PetscFunctionBeginUser; 1509566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1519566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw)); 1529566063dSJacob Faibussowitsch PetscCall(DMSetType(*sw, DMSWARM)); 1539566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*sw, dim)); 1549566063dSJacob Faibussowitsch PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC)); 1559566063dSJacob Faibussowitsch PetscCall(DMSwarmSetCellDM(*sw, dm)); 1569566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "kinematics", dim, PETSC_REAL)); 1579566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR)); 1589566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(*sw)); 1599566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 1609566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0)); 1619566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*sw)); 1629566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 163d0c080abSJoseph Pusztay for (c = cStart; c < cEnd; ++c) { 164d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 165d0c080abSJoseph Pusztay const PetscInt n = c * Np + p; 166d0c080abSJoseph Pusztay cellid[n] = c; 167d0c080abSJoseph Pusztay } 168d0c080abSJoseph Pusztay } 1699566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 1709566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles")); 1719566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view")); 172d0c080abSJoseph Pusztay PetscFunctionReturn(0); 173d0c080abSJoseph Pusztay } 174d0c080abSJoseph Pusztay 175d0c080abSJoseph Pusztay /* 176d0c080abSJoseph Pusztay f_t = 1/\tau \left( f_eq - f \right) 177d0c080abSJoseph Pusztay n_t = 1/\tau \left( \int f_eq - \int f \right) = 1/\tau (n - n) = 0 178d0c080abSJoseph Pusztay v_t = 1/\tau \left( \int v f_eq - \int v f \right) = 1/\tau (v - v) = 0 179d0c080abSJoseph Pusztay E_t = 1/\tau \left( \int v^2 f_eq - \int v^2 f \right) = 1/\tau (T - T) = 0 180d0c080abSJoseph Pusztay 181d0c080abSJoseph Pusztay Let's look at a single cell: 182d0c080abSJoseph Pusztay 183d0c080abSJoseph Pusztay \int_C f_t = 1/\tau \left( \int_C f_eq - \int_C f \right) 184d0c080abSJoseph Pusztay \sum_{x_i \in C} w^i_t = 1/\tau \left( F_eq(C) - \sum_{x_i \in C} w_i \right) 185d0c080abSJoseph Pusztay */ 186d0c080abSJoseph Pusztay 187d0c080abSJoseph Pusztay /* This computes the 1D Maxwellian distribution for given mass n, velocity v, and temperature T */ 188*9371c9d4SSatish Balay static PetscReal ComputePDF(PetscReal m, PetscReal n, PetscReal T, PetscReal v[]) { 189d0c080abSJoseph Pusztay return (n / PetscSqrtReal(2.0 * PETSC_PI * T / m)) * PetscExpReal(-0.5 * m * PetscSqr(v[0]) / T); 190d0c080abSJoseph Pusztay } 191d0c080abSJoseph Pusztay 192d0c080abSJoseph Pusztay /* 193d0c080abSJoseph Pusztay erf z = \frac{2}{\sqrt\pi} \int^z_0 dt e^{-t^2} and erf \infty = 1 194d0c080abSJoseph Pusztay 195d0c080abSJoseph Pusztay We begin with our distribution 196d0c080abSJoseph Pusztay 197d0c080abSJoseph Pusztay \sqrt{\frac{m}{2 \pi T}} e^{-m v^2/2T} 198d0c080abSJoseph Pusztay 199d0c080abSJoseph Pusztay Let t = \sqrt{m/2T} v, z = \sqrt{m/2T} w, so that we now have 200d0c080abSJoseph Pusztay 201d0c080abSJoseph Pusztay \sqrt{\frac{m}{2 \pi T}} \int^w_0 dv e^{-m v^2/2T} 202d0c080abSJoseph Pusztay = \sqrt{\frac{m}{2 \pi T}} \int^{\sqrt{m/2T} w}_0 \sqrt{2T/m} dt e^{-t^2} 203d0c080abSJoseph Pusztay = 1/\sqrt{\pi} \int^{\sqrt{m/2T} w}_0 dt e^{-t^2} 204d0c080abSJoseph Pusztay = 1/2 erf(\sqrt{m/2T} w) 205d0c080abSJoseph Pusztay */ 206*9371c9d4SSatish Balay static PetscReal ComputeCDF(PetscReal m, PetscReal n, PetscReal T, PetscReal va, PetscReal vb) { 207d0c080abSJoseph Pusztay PetscReal alpha = PetscSqrtReal(0.5 * m / T); 208d0c080abSJoseph Pusztay PetscReal za = alpha * va; 209d0c080abSJoseph Pusztay PetscReal zb = alpha * vb; 210d0c080abSJoseph Pusztay PetscReal sum = 0.0; 211d0c080abSJoseph Pusztay 212d0c080abSJoseph Pusztay sum += zb >= 0 ? erf(zb) : -erf(-zb); 213d0c080abSJoseph Pusztay sum -= za >= 0 ? erf(za) : -erf(-za); 214d0c080abSJoseph Pusztay return 0.5 * n * sum; 215d0c080abSJoseph Pusztay } 216d0c080abSJoseph Pusztay 217*9371c9d4SSatish Balay static PetscErrorCode CheckDistribution(DM dm, PetscReal m, PetscReal n, PetscReal T, PetscReal v[]) { 218d0c080abSJoseph Pusztay PetscSection coordSection; 219d0c080abSJoseph Pusztay Vec coordsLocal; 220d0c080abSJoseph Pusztay PetscReal *xq, *wq; 221d0c080abSJoseph Pusztay PetscReal vmin, vmax, neq, veq, Teq; 222d0c080abSJoseph Pusztay PetscInt Nq = 100, q, cStart, cEnd, c; 223d0c080abSJoseph Pusztay 2247510d9b0SBarry Smith PetscFunctionBeginUser; 2259566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(dm, &vmin, &vmax)); 226d0c080abSJoseph Pusztay /* Check analytic over entire line */ 227d0c080abSJoseph Pusztay neq = ComputeCDF(m, n, T, vmin, vmax); 22863a3b9bcSJacob Faibussowitsch PetscCheck(PetscAbsReal(neq - n) <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Int f %g != %g mass (%g)", (double)neq, (double)n, (double)(neq - n)); 229d0c080abSJoseph Pusztay /* Check analytic over cells */ 2309566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal)); 2319566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 2329566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 233d0c080abSJoseph Pusztay neq = 0.0; 234d0c080abSJoseph Pusztay for (c = cStart; c < cEnd; ++c) { 235d0c080abSJoseph Pusztay PetscScalar *vcoords = NULL; 236d0c080abSJoseph Pusztay 2379566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &vcoords)); 238d0c080abSJoseph Pusztay neq += ComputeCDF(m, n, T, vcoords[0], vcoords[1]); 2399566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &vcoords)); 240d0c080abSJoseph Pusztay } 24163a3b9bcSJacob Faibussowitsch PetscCheck(PetscAbsReal(neq - n) <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell Int f %g != %g mass (%g)", (double)neq, (double)n, (double)(neq - n)); 242d0c080abSJoseph Pusztay /* Check quadrature over entire line */ 2439566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Nq, &xq, Nq, &wq)); 2449566063dSJacob Faibussowitsch PetscCall(PetscDTGaussQuadrature(100, vmin, vmax, xq, wq)); 245d0c080abSJoseph Pusztay neq = 0.0; 246*9371c9d4SSatish Balay for (q = 0; q < Nq; ++q) { neq += ComputePDF(m, n, T, &xq[q]) * wq[q]; } 24763a3b9bcSJacob Faibussowitsch PetscCheck(PetscAbsReal(neq - n) <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Int f %g != %g mass (%g)", (double)neq, (double)n, (double)(neq - n)); 248d0c080abSJoseph Pusztay /* Check omemnts with quadrature */ 249d0c080abSJoseph Pusztay veq = 0.0; 250*9371c9d4SSatish Balay for (q = 0; q < Nq; ++q) { veq += xq[q] * ComputePDF(m, n, T, &xq[q]) * wq[q]; } 251d0c080abSJoseph Pusztay veq /= neq; 25263a3b9bcSJacob Faibussowitsch PetscCheck(PetscAbsReal(veq - v[0]) <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Int v f %g != %g velocity (%g)", (double)veq, (double)v[0], (double)(veq - v[0])); 253d0c080abSJoseph Pusztay Teq = 0.0; 254*9371c9d4SSatish Balay for (q = 0; q < Nq; ++q) { Teq += PetscSqr(xq[q]) * ComputePDF(m, n, T, &xq[q]) * wq[q]; } 255d0c080abSJoseph Pusztay Teq = Teq * m / neq - PetscSqr(veq); 25663a3b9bcSJacob Faibussowitsch PetscCheck(PetscAbsReal(Teq - T) <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Int v^2 f %g != %g temperature (%g)", (double)Teq, (double)T, (double)(Teq - T)); 2579566063dSJacob Faibussowitsch PetscCall(PetscFree2(xq, wq)); 258d0c080abSJoseph Pusztay PetscFunctionReturn(0); 259d0c080abSJoseph Pusztay } 260d0c080abSJoseph Pusztay 261*9371c9d4SSatish Balay static PetscErrorCode RHSFunctionParticles(TS ts, PetscReal t, Vec U, Vec R, void *ctx) { 262d0c080abSJoseph Pusztay const PetscScalar *u; 263d0c080abSJoseph Pusztay PetscSection coordSection; 264d0c080abSJoseph Pusztay Vec coordsLocal; 265d0c080abSJoseph Pusztay PetscScalar *r, *coords; 266d0c080abSJoseph Pusztay PetscReal n = 0.0, v = 0.0, E = 0.0, T = 0.0, m = 1.0, cn = 0.0, cv = 0.0, cE = 0.0, pE = 0.0, eqE = 0.0; 267d0c080abSJoseph Pusztay PetscInt dim, d, Np, Ncp, p, cStart, cEnd, c; 268d0c080abSJoseph Pusztay DM dmSw, plex; 269d0c080abSJoseph Pusztay 270d0c080abSJoseph Pusztay PetscFunctionBeginUser; 2719566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(U, &Np)); 2729566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 2739566063dSJacob Faibussowitsch PetscCall(VecGetArray(R, &r)); 2749566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dmSw)); 2759566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dmSw, &plex)); 2769566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmSw, &dim)); 2779566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(plex, &coordsLocal)); 2789566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(plex, &coordSection)); 2799566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd)); 280d0c080abSJoseph Pusztay Np /= dim; 281d0c080abSJoseph Pusztay Ncp = Np / (cEnd - cStart); 282d0c080abSJoseph Pusztay /* Calculate moments of particle distribution, note that velocity is in the coordinate */ 2839566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 284d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 285d0c080abSJoseph Pusztay PetscReal m1 = 0.0, m2 = 0.0; 286d0c080abSJoseph Pusztay 287*9371c9d4SSatish Balay for (d = 0; d < dim; ++d) { 288*9371c9d4SSatish Balay m1 += PetscRealPart(coords[p * dim + d]); 289*9371c9d4SSatish Balay m2 += PetscSqr(PetscRealPart(coords[p * dim + d])); 290*9371c9d4SSatish Balay } 291d0c080abSJoseph Pusztay n += u[p]; 292d0c080abSJoseph Pusztay v += u[p] * m1; 293d0c080abSJoseph Pusztay E += u[p] * m2; 294d0c080abSJoseph Pusztay } 295d0c080abSJoseph Pusztay v /= n; 296d0c080abSJoseph Pusztay T = E * m / n - v * v; 29763a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Time %.2f: mass %.4f velocity: %+.4f temperature: %.4f\n", (double)t, (double)n, (double)v, (double)T)); 2989566063dSJacob Faibussowitsch PetscCall(CheckDistribution(plex, m, n, T, &v)); 299d0c080abSJoseph Pusztay /* 300d0c080abSJoseph Pusztay Begin cellwise evaluation of the collision operator. Essentially, penalize the weights of the particles 301d0c080abSJoseph Pusztay in that cell against the maxwellian for the number of particles expected to be in that cell 302d0c080abSJoseph Pusztay */ 303d0c080abSJoseph Pusztay for (c = cStart; c < cEnd; ++c) { 304d0c080abSJoseph Pusztay PetscScalar *vcoords = NULL; 305d0c080abSJoseph Pusztay PetscReal relaxation = 1.0, neq; 306d0c080abSJoseph Pusztay PetscInt sp = c * Ncp, q; 307d0c080abSJoseph Pusztay 308d0c080abSJoseph Pusztay /* Calculate equilibrium occupation for this velocity cell */ 3099566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(plex, coordSection, coordsLocal, c, NULL, &vcoords)); 310d0c080abSJoseph Pusztay neq = ComputeCDF(m, n, T, vcoords[0], vcoords[1]); 3119566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(plex, coordSection, coordsLocal, c, NULL, &vcoords)); 312d0c080abSJoseph Pusztay for (q = 0; q < Ncp; ++q) r[sp + q] = (1.0 / relaxation) * (neq - u[sp + q]); 313d0c080abSJoseph Pusztay } 314d0c080abSJoseph Pusztay /* Check update */ 315d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 316d0c080abSJoseph Pusztay PetscReal m1 = 0.0, m2 = 0.0; 317d0c080abSJoseph Pusztay PetscScalar *vcoords = NULL; 318d0c080abSJoseph Pusztay 319*9371c9d4SSatish Balay for (d = 0; d < dim; ++d) { 320*9371c9d4SSatish Balay m1 += PetscRealPart(coords[p * dim + d]); 321*9371c9d4SSatish Balay m2 += PetscSqr(PetscRealPart(coords[p * dim + d])); 322*9371c9d4SSatish Balay } 323d0c080abSJoseph Pusztay cn += r[p]; 324d0c080abSJoseph Pusztay cv += r[p] * m1; 325d0c080abSJoseph Pusztay cE += r[p] * m2; 326d0c080abSJoseph Pusztay pE += u[p] * m2; 3279566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(plex, coordSection, coordsLocal, p, NULL, &vcoords)); 328d0c080abSJoseph Pusztay eqE += ComputeCDF(m, n, T, vcoords[0], vcoords[1]) * m2; 3299566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(plex, coordSection, coordsLocal, p, NULL, &vcoords)); 330d0c080abSJoseph Pusztay } 33163a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Time %.2f: mass update %.8f velocity update: %+.8f energy update: %.8f (%.8f, %.8f)\n", (double)t, (double)cn, (double)cv, (double)cE, (double)pE, (double)eqE)); 3329566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 3339566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 3349566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(R, &r)); 335d0c080abSJoseph Pusztay PetscFunctionReturn(0); 336d0c080abSJoseph Pusztay } 337d0c080abSJoseph Pusztay 338*9371c9d4SSatish Balay static PetscErrorCode HGMonitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) { 339d0c080abSJoseph Pusztay AppCtx *user = (AppCtx *)ctx; 340d0c080abSJoseph Pusztay const PetscScalar *u; 341d0c080abSJoseph Pusztay DM sw, dm; 342d0c080abSJoseph Pusztay PetscInt dim, Np, p; 343d0c080abSJoseph Pusztay 344d0c080abSJoseph Pusztay PetscFunctionBeginUser; 345d0c080abSJoseph Pusztay if (step < 0) PetscFunctionReturn(0); 346d0c080abSJoseph Pusztay if (((user->ostep > 0) && (!(step % user->ostep)))) { 347d0c080abSJoseph Pusztay PetscDrawAxis axis; 348d0c080abSJoseph Pusztay 3499566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &sw)); 3509566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(sw, &dm)); 3519566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 3529566063dSJacob Faibussowitsch PetscCall(PetscDrawHGReset(user->drawhg)); 3539566063dSJacob Faibussowitsch PetscCall(PetscDrawHGGetAxis(user->drawhg, &axis)); 3549566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "V", "f(V)")); 3559566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLimits(axis, -3, 3, 0, 100)); 3569566063dSJacob Faibussowitsch PetscCall(PetscDrawHGSetLimits(user->drawhg, -3.0, 3.0, 0, 10)); 3579566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(U, &Np)); 358d0c080abSJoseph Pusztay Np /= dim; 3599566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 360d0c080abSJoseph Pusztay /* get points from solution vector */ 3619566063dSJacob Faibussowitsch for (p = 0; p < Np; ++p) PetscCall(PetscDrawHGAddValue(user->drawhg, u[p])); 3629566063dSJacob Faibussowitsch PetscCall(PetscDrawHGDraw(user->drawhg)); 3639566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 364d0c080abSJoseph Pusztay } 365d0c080abSJoseph Pusztay PetscFunctionReturn(0); 366d0c080abSJoseph Pusztay } 367d0c080abSJoseph Pusztay 368*9371c9d4SSatish Balay static PetscErrorCode SPMonitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) { 369d0c080abSJoseph Pusztay AppCtx *user = (AppCtx *)ctx; 370d0c080abSJoseph Pusztay const PetscScalar *u; 371d0c080abSJoseph Pusztay PetscReal *v, *coords; 372d0c080abSJoseph Pusztay PetscInt Np, p; 373d0c080abSJoseph Pusztay DM dmSw; 374d0c080abSJoseph Pusztay 375d0c080abSJoseph Pusztay PetscFunctionBeginUser; 376d0c080abSJoseph Pusztay 377d0c080abSJoseph Pusztay if (step < 0) PetscFunctionReturn(0); 378d0c080abSJoseph Pusztay if (((user->ostep > 0) && (!(step % user->ostep)))) { 379d0c080abSJoseph Pusztay PetscDrawAxis axis; 380d0c080abSJoseph Pusztay 3819566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dmSw)); 3829566063dSJacob Faibussowitsch PetscCall(PetscDrawSPReset(user->drawsp)); 3839566063dSJacob Faibussowitsch PetscCall(PetscDrawSPGetAxis(user->drawsp, &axis)); 3849566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "V", "w")); 3859566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(U, &Np)); 3869566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 387d0c080abSJoseph Pusztay /* get points from solution vector */ 3889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Np, &v)); 389d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) v[p] = PetscRealPart(u[p]); 3909566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 3919566063dSJacob Faibussowitsch for (p = 0; p < Np - 1; ++p) PetscCall(PetscDrawSPAddPoint(user->drawsp, &coords[p], &v[p])); 3929566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDraw(user->drawsp, PETSC_TRUE)); 3939566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 3949566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 3959566063dSJacob Faibussowitsch PetscCall(PetscFree(v)); 396d0c080abSJoseph Pusztay } 397d0c080abSJoseph Pusztay PetscFunctionReturn(0); 398d0c080abSJoseph Pusztay } 399d0c080abSJoseph Pusztay 400*9371c9d4SSatish Balay static PetscErrorCode KSConv(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) { 401d0c080abSJoseph Pusztay AppCtx *user = (AppCtx *)ctx; 402d0c080abSJoseph Pusztay const PetscScalar *u; 403d0c080abSJoseph Pusztay PetscScalar sup; 404d0c080abSJoseph Pusztay PetscReal *v, *coords, T = 0., vel = 0., step_cast, w_sum; 405d0c080abSJoseph Pusztay PetscInt dim, Np, p, cStart, cEnd; 406d0c080abSJoseph Pusztay DM sw, plex; 407d0c080abSJoseph Pusztay 408d0c080abSJoseph Pusztay PetscFunctionBeginUser; 409d0c080abSJoseph Pusztay if (step < 0) PetscFunctionReturn(0); 410d0c080abSJoseph Pusztay if (((user->ostep > 0) && (!(step % user->ostep)))) { 411d0c080abSJoseph Pusztay PetscDrawAxis axis; 4129566063dSJacob Faibussowitsch PetscCall(PetscDrawSPGetAxis(user->drawks, &axis)); 4139566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "t", "D_n")); 4149566063dSJacob Faibussowitsch PetscCall(PetscDrawSPSetLimits(user->drawks, 0., 100, 0., 3.5)); 4159566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &sw)); 4169566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(U, &Np)); 4179566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 418d0c080abSJoseph Pusztay /* get points from solution vector */ 4199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Np, &v)); 4209566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(sw, &plex)); 4219566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 4229566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd)); 423d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) v[p] = PetscRealPart(u[p]); 4249566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 425d0c080abSJoseph Pusztay w_sum = 0.; 426d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 427d0c080abSJoseph Pusztay w_sum += u[p]; 428d0c080abSJoseph Pusztay T += u[p] * coords[p] * coords[p]; 429d0c080abSJoseph Pusztay vel += u[p] * coords[p]; 430d0c080abSJoseph Pusztay } 431d0c080abSJoseph Pusztay vel /= w_sum; 432d0c080abSJoseph Pusztay T = T / w_sum - vel * vel; 433d0c080abSJoseph Pusztay sup = 0.0; 434d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 435d0c080abSJoseph Pusztay PetscReal tmp = 0.; 436d0c080abSJoseph Pusztay tmp = PetscAbs(u[p] - ComputePDF(1.0, w_sum, T, &coords[p * dim])); 437d0c080abSJoseph Pusztay if (tmp > sup) sup = tmp; 438d0c080abSJoseph Pusztay } 439d0c080abSJoseph Pusztay step_cast = (PetscReal)step; 4409566063dSJacob Faibussowitsch PetscCall(PetscDrawSPAddPoint(user->drawks, &step_cast, &sup)); 4419566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDraw(user->drawks, PETSC_TRUE)); 4429566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 4439566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 4449566063dSJacob Faibussowitsch PetscCall(PetscFree(v)); 445d0c080abSJoseph Pusztay } 446d0c080abSJoseph Pusztay PetscFunctionReturn(0); 447d0c080abSJoseph Pusztay } 448d0c080abSJoseph Pusztay 449*9371c9d4SSatish Balay static PetscErrorCode InitializeSolve(TS ts, Vec u) { 450d0c080abSJoseph Pusztay DM dm; 451d0c080abSJoseph Pusztay AppCtx *user; 452d0c080abSJoseph Pusztay 453d0c080abSJoseph Pusztay PetscFunctionBeginUser; 4549566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 4559566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(dm, &user)); 4569566063dSJacob Faibussowitsch PetscCall(SetInitialCoordinates(dm)); 4579566063dSJacob Faibussowitsch PetscCall(SetInitialConditions(dm, u)); 458d0c080abSJoseph Pusztay PetscFunctionReturn(0); 459d0c080abSJoseph Pusztay } 460d0c080abSJoseph Pusztay /* 461d0c080abSJoseph Pusztay A single particle is generated for each velocity space cell using the dmswarmpicfield_coor and is used to evaluate collisions in that cell. 462d0c080abSJoseph Pusztay 0 weight ghost particles are initialized outside of a small velocity domain to ensure the tails of the amxwellian are resolved. 463d0c080abSJoseph Pusztay */ 464*9371c9d4SSatish Balay int main(int argc, char **argv) { 465d0c080abSJoseph Pusztay TS ts; /* nonlinear solver */ 466d0c080abSJoseph Pusztay DM dm, sw; /* Velocity space mesh and Particle Swarm */ 467d0c080abSJoseph Pusztay Vec u, w; /* swarm vector */ 468d0c080abSJoseph Pusztay MPI_Comm comm; 469d0c080abSJoseph Pusztay AppCtx user; 470d0c080abSJoseph Pusztay 471327415f7SBarry Smith PetscFunctionBeginUser; 4729566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 473d0c080abSJoseph Pusztay comm = PETSC_COMM_WORLD; 4749566063dSJacob Faibussowitsch PetscCall(ProcessOptions(comm, &user)); 4759566063dSJacob Faibussowitsch PetscCall(CreateMesh(comm, &dm, &user)); 4769566063dSJacob Faibussowitsch PetscCall(CreateParticles(dm, &sw, &user)); 4779566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(sw, &user)); 4789566063dSJacob Faibussowitsch PetscCall(TSCreate(comm, &ts)); 4799566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, sw)); 4809566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 10.0)); 4819566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 0.01)); 4829566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, 100000)); 4839566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 484d0c080abSJoseph Pusztay if (user.monitorhg) { 4859566063dSJacob Faibussowitsch PetscCall(PetscDrawCreate(comm, NULL, "monitor", 0, 0, 400, 300, &user.draw)); 4869566063dSJacob Faibussowitsch PetscCall(PetscDrawSetFromOptions(user.draw)); 4879566063dSJacob Faibussowitsch PetscCall(PetscDrawHGCreate(user.draw, 20, &user.drawhg)); 4889566063dSJacob Faibussowitsch PetscCall(PetscDrawHGSetColor(user.drawhg, 3)); 4899566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts, HGMonitor, &user, NULL)); 490*9371c9d4SSatish Balay } else if (user.monitorsp) { 4919566063dSJacob Faibussowitsch PetscCall(PetscDrawCreate(comm, NULL, "monitor", 0, 0, 400, 300, &user.draw)); 4929566063dSJacob Faibussowitsch PetscCall(PetscDrawSetFromOptions(user.draw)); 4939566063dSJacob Faibussowitsch PetscCall(PetscDrawSPCreate(user.draw, 1, &user.drawsp)); 4949566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts, SPMonitor, &user, NULL)); 495*9371c9d4SSatish Balay } else if (user.monitorks) { 4969566063dSJacob Faibussowitsch PetscCall(PetscDrawCreate(comm, NULL, "monitor", 0, 0, 400, 300, &user.draw)); 4979566063dSJacob Faibussowitsch PetscCall(PetscDrawSetFromOptions(user.draw)); 4989566063dSJacob Faibussowitsch PetscCall(PetscDrawSPCreate(user.draw, 1, &user.drawks)); 4999566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts, KSConv, &user, NULL)); 500d0c080abSJoseph Pusztay } 5019566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunctionParticles, &user)); 5029566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 5039566063dSJacob Faibussowitsch PetscCall(TSSetComputeInitialCondition(ts, InitializeSolve)); 5049566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &w)); 5059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(w, &u)); 5069566063dSJacob Faibussowitsch PetscCall(VecCopy(w, u)); 5079566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &w)); 5089566063dSJacob Faibussowitsch PetscCall(TSComputeInitialCondition(ts, u)); 5099566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, u)); 510d0c080abSJoseph Pusztay if (user.monitorhg) { 5119566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(user.draw)); 5129566063dSJacob Faibussowitsch PetscCall(PetscDrawHGDestroy(&user.drawhg)); 5139566063dSJacob Faibussowitsch PetscCall(PetscDrawDestroy(&user.draw)); 514d0c080abSJoseph Pusztay } 515d0c080abSJoseph Pusztay if (user.monitorsp) { 5169566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(user.draw)); 5179566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDestroy(&user.drawsp)); 5189566063dSJacob Faibussowitsch PetscCall(PetscDrawDestroy(&user.draw)); 519d0c080abSJoseph Pusztay } 520d0c080abSJoseph Pusztay if (user.monitorks) { 5219566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(user.draw)); 5229566063dSJacob Faibussowitsch PetscCall(PetscDrawSPDestroy(&user.drawks)); 5239566063dSJacob Faibussowitsch PetscCall(PetscDrawDestroy(&user.draw)); 524d0c080abSJoseph Pusztay } 5259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 5269566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 5279566063dSJacob Faibussowitsch PetscCall(DMDestroy(&sw)); 5289566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 5299566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 530b122ec5aSJacob Faibussowitsch return 0; 531d0c080abSJoseph Pusztay } 532d0c080abSJoseph Pusztay 533d0c080abSJoseph Pusztay /*TEST 534d0c080abSJoseph Pusztay build: 53530602db0SMatthew G. Knepley requires: double !complex 536d0c080abSJoseph Pusztay test: 537d0c080abSJoseph Pusztay suffix: 1 53830602db0SMatthew G. Knepley args: -particles_per_cell 1 -output_step 10 -ts_type euler -dm_plex_dim 1 -dm_plex_box_faces 200 -dm_plex_box_lower -10 -dm_plex_box_upper 10 -dm_view -monitorsp 539d0c080abSJoseph Pusztay test: 540d0c080abSJoseph Pusztay suffix: 2 54130602db0SMatthew G. Knepley args: -particles_per_cell 1 -output_step 50 -ts_type euler -dm_plex_dim 1 -dm_plex_box_faces 200 -dm_plex_box_lower -10 -dm_plex_box_upper 10 -dm_view -monitorks 542d0c080abSJoseph Pusztay TEST*/ 543