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 25d0c080abSJoseph Pusztay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 26d0c080abSJoseph Pusztay { 27d0c080abSJoseph Pusztay PetscErrorCode ierr; 28d0c080abSJoseph Pusztay 29d0c080abSJoseph Pusztay PetscFunctionBeginUser; 30d0c080abSJoseph Pusztay options->monitorhg = PETSC_FALSE; 31d0c080abSJoseph Pusztay options->monitorsp = PETSC_FALSE; 32d0c080abSJoseph Pusztay options->monitorks = PETSC_FALSE; 33d0c080abSJoseph Pusztay options->particlesPerCell = 1; 34d0c080abSJoseph Pusztay options->momentTol = 100.0*PETSC_MACHINE_EPSILON; 35d0c080abSJoseph Pusztay options->ostep = 100; 36d0c080abSJoseph Pusztay 37d0c080abSJoseph Pusztay ierr = PetscOptionsBegin(comm, "", "Collision Options", "DMPLEX");CHKERRQ(ierr); 38*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-monitorhg", "Flag to use the TS histogram monitor", "ex28.c", options->monitorhg, &options->monitorhg, NULL)); 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-monitorsp", "Flag to use the TS scatter plot monitor", "ex28.c", options->monitorsp, &options->monitorsp, NULL)); 40*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-monitorks", "Flag to plot KS test results", "ex28.c", options->monitorks, &options->monitorks, NULL)); 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-particles_per_cell", "Number of particles per cell", "ex28.c", options->particlesPerCell, &options->particlesPerCell, NULL)); 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-output_step", "Number of time steps between output", "ex28.c", options->ostep, &options->ostep, PETSC_NULL)); 43d0c080abSJoseph Pusztay ierr = PetscOptionsEnd();CHKERRQ(ierr); 44d0c080abSJoseph Pusztay 45d0c080abSJoseph Pusztay PetscFunctionReturn(0); 46d0c080abSJoseph Pusztay } 47d0c080abSJoseph Pusztay 48d0c080abSJoseph Pusztay /* Create the mesh for velocity space */ 49d0c080abSJoseph Pusztay static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user) 50d0c080abSJoseph Pusztay { 51d0c080abSJoseph Pusztay PetscFunctionBeginUser; 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 54*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 56d0c080abSJoseph Pusztay PetscFunctionReturn(0); 57d0c080abSJoseph Pusztay } 58d0c080abSJoseph Pusztay 59d0c080abSJoseph Pusztay /* Since we are putting the same number of particles in each cell, this amounts to a uniform distribution of v */ 60d0c080abSJoseph Pusztay static PetscErrorCode SetInitialCoordinates(DM sw) 61d0c080abSJoseph Pusztay { 62d0c080abSJoseph Pusztay AppCtx *user; 63d0c080abSJoseph Pusztay PetscRandom rnd; 64d0c080abSJoseph Pusztay DM dm; 65d0c080abSJoseph Pusztay DMPolytopeType ct; 66d0c080abSJoseph Pusztay PetscBool simplex; 67d0c080abSJoseph Pusztay PetscReal *centroid, *coords, *xi0, *v0, *J, *invJ, detJ, *vals; 68d0c080abSJoseph Pusztay PetscInt dim, d, cStart, cEnd, c, Np, p; 69d0c080abSJoseph Pusztay 70d0c080abSJoseph Pusztay PetscFunctionBeginUser; 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PetscObjectComm((PetscObject) sw), &rnd)); 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetInterval(rnd, -1.0, 1.0)); 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rnd)); 74d0c080abSJoseph Pusztay 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetApplicationContext(sw, &user)); 76d0c080abSJoseph Pusztay Np = user->particlesPerCell; 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(sw, &dim)); 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetCellDM(sw, &dm)); 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCellType(dm, cStart, &ct)); 81d0c080abSJoseph Pusztay simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc5(dim, ¢roid, dim, &xi0, dim, &v0, dim*dim, &J, dim*dim, &invJ)); 83d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) xi0[d] = -1.0; 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords)); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **) &vals)); 86d0c080abSJoseph Pusztay for (c = cStart; c < cEnd; ++c) { 87d0c080abSJoseph Pusztay if (Np == 1) { 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL)); 89d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) { 90d0c080abSJoseph Pusztay coords[c*dim+d] = centroid[d]; 91d0c080abSJoseph Pusztay if ((coords[c*dim+d] >= -1) && (coords[c*dim+d] <= 1)) { 92d0c080abSJoseph Pusztay vals[c] = 1.0; 93d0c080abSJoseph Pusztay } else { 94d0c080abSJoseph Pusztay vals[c] = 0.; 95d0c080abSJoseph Pusztay } 96d0c080abSJoseph Pusztay } 97d0c080abSJoseph Pusztay } else { 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); /* affine */ 99d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 100d0c080abSJoseph Pusztay const PetscInt n = c*Np + p; 101d0c080abSJoseph Pusztay PetscReal sum = 0.0, refcoords[3]; 102d0c080abSJoseph Pusztay 103d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) { 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValueReal(rnd, &refcoords[d])); 105d0c080abSJoseph Pusztay sum += refcoords[d]; 106d0c080abSJoseph Pusztay } 107d0c080abSJoseph Pusztay if (simplex && sum > 0.0) for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim)*sum; 108d0c080abSJoseph Pusztay vals[n] = 1.0; 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexReferenceToCoordinates(dm, c, 1, refcoords, &coords[n*dim])); 110d0c080abSJoseph Pusztay } 111d0c080abSJoseph Pusztay } 112d0c080abSJoseph Pusztay } 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords)); 114*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **) &vals)); 115*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree5(centroid, xi0, v0, J, invJ)); 116*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rnd)); 117d0c080abSJoseph Pusztay PetscFunctionReturn(0); 118d0c080abSJoseph Pusztay } 119d0c080abSJoseph Pusztay 120d0c080abSJoseph Pusztay /* The intiial conditions are just the initial particle weights */ 121d0c080abSJoseph Pusztay static PetscErrorCode SetInitialConditions(DM dmSw, Vec u) 122d0c080abSJoseph Pusztay { 123d0c080abSJoseph Pusztay DM dm; 124d0c080abSJoseph Pusztay AppCtx *user; 125d0c080abSJoseph Pusztay PetscReal *vals; 126d0c080abSJoseph Pusztay PetscScalar *initialConditions; 127d0c080abSJoseph Pusztay PetscInt dim, d, cStart, cEnd, c, Np, p, n; 128d0c080abSJoseph Pusztay 129d0c080abSJoseph Pusztay PetscFunctionBeginUser; 130*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(u, &n)); 131*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetApplicationContext(dmSw, &user)); 132d0c080abSJoseph Pusztay Np = user->particlesPerCell; 133*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetCellDM(dmSw, &dm)); 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 135*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 1363c633725SBarry Smith PetscCheck(n == (cEnd-cStart)*Np,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "TS solution local size %D != %D nm particles", n, (cEnd-cStart)*Np); 137*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(dmSw, "w_q", NULL, NULL, (void **) &vals)); 138*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(u, &initialConditions)); 139d0c080abSJoseph Pusztay for (c = cStart; c < cEnd; ++c) { 140d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 141d0c080abSJoseph Pusztay const PetscInt n = c*Np + p; 142d0c080abSJoseph Pusztay for (d = 0; d < dim; d++) { 143d0c080abSJoseph Pusztay initialConditions[n] = vals[n]; 144d0c080abSJoseph Pusztay } 145d0c080abSJoseph Pusztay } 146d0c080abSJoseph Pusztay } 147*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(u, &initialConditions)); 148*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(dmSw, "w_q", NULL, NULL, (void **) &vals)); 149d0c080abSJoseph Pusztay PetscFunctionReturn(0); 150d0c080abSJoseph Pusztay } 151d0c080abSJoseph Pusztay 152d0c080abSJoseph Pusztay static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user) 153d0c080abSJoseph Pusztay { 154d0c080abSJoseph Pusztay PetscInt *cellid; 155d0c080abSJoseph Pusztay PetscInt dim, cStart, cEnd, c, Np = user->particlesPerCell, p; 156d0c080abSJoseph Pusztay 157d0c080abSJoseph Pusztay PetscFunctionBeginUser; 158*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 159*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(PetscObjectComm((PetscObject) dm), sw)); 160*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*sw, DMSWARM)); 161*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetDimension(*sw, dim)); 162*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetType(*sw, DMSWARM_PIC)); 163*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetCellDM(*sw, dm)); 164*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRegisterPetscDatatypeField(*sw, "kinematics", dim, PETSC_REAL)); 165*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR)); 166*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmFinalizeFieldRegister(*sw)); 167*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 168*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0)); 169*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*sw)); 170*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid)); 171d0c080abSJoseph Pusztay for (c = cStart; c < cEnd; ++c) { 172d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 173d0c080abSJoseph Pusztay const PetscInt n = c*Np + p; 174d0c080abSJoseph Pusztay cellid[n] = c; 175d0c080abSJoseph Pusztay } 176d0c080abSJoseph Pusztay } 177*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid)); 178*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) *sw, "Particles")); 179*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*sw, NULL, "-sw_view")); 180d0c080abSJoseph Pusztay PetscFunctionReturn(0); 181d0c080abSJoseph Pusztay } 182d0c080abSJoseph Pusztay 183d0c080abSJoseph Pusztay /* 184d0c080abSJoseph Pusztay f_t = 1/\tau \left( f_eq - f \right) 185d0c080abSJoseph Pusztay n_t = 1/\tau \left( \int f_eq - \int f \right) = 1/\tau (n - n) = 0 186d0c080abSJoseph Pusztay v_t = 1/\tau \left( \int v f_eq - \int v f \right) = 1/\tau (v - v) = 0 187d0c080abSJoseph Pusztay E_t = 1/\tau \left( \int v^2 f_eq - \int v^2 f \right) = 1/\tau (T - T) = 0 188d0c080abSJoseph Pusztay 189d0c080abSJoseph Pusztay Let's look at a single cell: 190d0c080abSJoseph Pusztay 191d0c080abSJoseph Pusztay \int_C f_t = 1/\tau \left( \int_C f_eq - \int_C f \right) 192d0c080abSJoseph Pusztay \sum_{x_i \in C} w^i_t = 1/\tau \left( F_eq(C) - \sum_{x_i \in C} w_i \right) 193d0c080abSJoseph Pusztay */ 194d0c080abSJoseph Pusztay 195d0c080abSJoseph Pusztay /* This computes the 1D Maxwellian distribution for given mass n, velocity v, and temperature T */ 196d0c080abSJoseph Pusztay static PetscReal ComputePDF(PetscReal m, PetscReal n, PetscReal T, PetscReal v[]) 197d0c080abSJoseph Pusztay { 198d0c080abSJoseph Pusztay return (n/PetscSqrtReal(2.0*PETSC_PI*T/m)) * PetscExpReal(-0.5*m*PetscSqr(v[0])/T); 199d0c080abSJoseph Pusztay } 200d0c080abSJoseph Pusztay 201d0c080abSJoseph Pusztay /* 202d0c080abSJoseph Pusztay erf z = \frac{2}{\sqrt\pi} \int^z_0 dt e^{-t^2} and erf \infty = 1 203d0c080abSJoseph Pusztay 204d0c080abSJoseph Pusztay We begin with our distribution 205d0c080abSJoseph Pusztay 206d0c080abSJoseph Pusztay \sqrt{\frac{m}{2 \pi T}} e^{-m v^2/2T} 207d0c080abSJoseph Pusztay 208d0c080abSJoseph Pusztay Let t = \sqrt{m/2T} v, z = \sqrt{m/2T} w, so that we now have 209d0c080abSJoseph Pusztay 210d0c080abSJoseph Pusztay \sqrt{\frac{m}{2 \pi T}} \int^w_0 dv e^{-m v^2/2T} 211d0c080abSJoseph Pusztay = \sqrt{\frac{m}{2 \pi T}} \int^{\sqrt{m/2T} w}_0 \sqrt{2T/m} dt e^{-t^2} 212d0c080abSJoseph Pusztay = 1/\sqrt{\pi} \int^{\sqrt{m/2T} w}_0 dt e^{-t^2} 213d0c080abSJoseph Pusztay = 1/2 erf(\sqrt{m/2T} w) 214d0c080abSJoseph Pusztay */ 215d0c080abSJoseph Pusztay static PetscReal ComputeCDF(PetscReal m, PetscReal n, PetscReal T, PetscReal va, PetscReal vb) 216d0c080abSJoseph Pusztay { 217d0c080abSJoseph Pusztay PetscReal alpha = PetscSqrtReal(0.5*m/T); 218d0c080abSJoseph Pusztay PetscReal za = alpha*va; 219d0c080abSJoseph Pusztay PetscReal zb = alpha*vb; 220d0c080abSJoseph Pusztay PetscReal sum = 0.0; 221d0c080abSJoseph Pusztay 222d0c080abSJoseph Pusztay sum += zb >= 0 ? erf(zb) : -erf(-zb); 223d0c080abSJoseph Pusztay sum -= za >= 0 ? erf(za) : -erf(-za); 224d0c080abSJoseph Pusztay return 0.5 * n * sum; 225d0c080abSJoseph Pusztay } 226d0c080abSJoseph Pusztay 227d0c080abSJoseph Pusztay static PetscErrorCode CheckDistribution(DM dm, PetscReal m, PetscReal n, PetscReal T, PetscReal v[]) 228d0c080abSJoseph Pusztay { 229d0c080abSJoseph Pusztay PetscSection coordSection; 230d0c080abSJoseph Pusztay Vec coordsLocal; 231d0c080abSJoseph Pusztay PetscReal *xq, *wq; 232d0c080abSJoseph Pusztay PetscReal vmin, vmax, neq, veq, Teq; 233d0c080abSJoseph Pusztay PetscInt Nq = 100, q, cStart, cEnd, c; 234d0c080abSJoseph Pusztay 235d0c080abSJoseph Pusztay PetscFunctionBegin; 236*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetBoundingBox(dm, &vmin, &vmax)); 237d0c080abSJoseph Pusztay /* Check analytic over entire line */ 238d0c080abSJoseph Pusztay neq = ComputeCDF(m, n, T, vmin, vmax); 2393c633725SBarry Smith PetscCheck(PetscAbsReal(neq - n) <= PETSC_SMALL,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Int f %g != %g mass (%g)", neq, n, neq-n); 240d0c080abSJoseph Pusztay /* Check analytic over cells */ 241*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocal(dm, &coordsLocal)); 242*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateSection(dm, &coordSection)); 243*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 244d0c080abSJoseph Pusztay neq = 0.0; 245d0c080abSJoseph Pusztay for (c = cStart; c < cEnd; ++c) { 246d0c080abSJoseph Pusztay PetscScalar *vcoords = NULL; 247d0c080abSJoseph Pusztay 248*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &vcoords)); 249d0c080abSJoseph Pusztay neq += ComputeCDF(m, n, T, vcoords[0], vcoords[1]); 250*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &vcoords)); 251d0c080abSJoseph Pusztay } 2523c633725SBarry Smith PetscCheck(PetscAbsReal(neq - n) <= PETSC_SMALL,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell Int f %g != %g mass (%g)", neq, n, neq-n); 253d0c080abSJoseph Pusztay /* Check quadrature over entire line */ 254*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(Nq, &xq, Nq, &wq)); 255*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTGaussQuadrature(100, vmin, vmax, xq, wq)); 256d0c080abSJoseph Pusztay neq = 0.0; 257d0c080abSJoseph Pusztay for (q = 0; q < Nq; ++q) { 258d0c080abSJoseph Pusztay neq += ComputePDF(m, n, T, &xq[q])*wq[q]; 259d0c080abSJoseph Pusztay } 2603c633725SBarry Smith PetscCheck(PetscAbsReal(neq - n) <= PETSC_SMALL,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Int f %g != %g mass (%g)", neq, n, neq-n); 261d0c080abSJoseph Pusztay /* Check omemnts with quadrature */ 262d0c080abSJoseph Pusztay veq = 0.0; 263d0c080abSJoseph Pusztay for (q = 0; q < Nq; ++q) { 264d0c080abSJoseph Pusztay veq += xq[q]*ComputePDF(m, n, T, &xq[q])*wq[q]; 265d0c080abSJoseph Pusztay } 266d0c080abSJoseph Pusztay veq /= neq; 2673c633725SBarry Smith PetscCheck(PetscAbsReal(veq - v[0]) <= PETSC_SMALL,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Int v f %g != %g velocity (%g)", veq, v[0], veq-v[0]); 268d0c080abSJoseph Pusztay Teq = 0.0; 269d0c080abSJoseph Pusztay for (q = 0; q < Nq; ++q) { 270d0c080abSJoseph Pusztay Teq += PetscSqr(xq[q])*ComputePDF(m, n, T, &xq[q])*wq[q]; 271d0c080abSJoseph Pusztay } 272d0c080abSJoseph Pusztay Teq = Teq * m/neq - PetscSqr(veq); 2733c633725SBarry Smith PetscCheck(PetscAbsReal(Teq - T) <= PETSC_SMALL,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Int v^2 f %g != %g temperature (%g)", Teq, T, Teq-T); 274*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(xq, wq)); 275d0c080abSJoseph Pusztay PetscFunctionReturn(0); 276d0c080abSJoseph Pusztay } 277d0c080abSJoseph Pusztay 278d0c080abSJoseph Pusztay static PetscErrorCode RHSFunctionParticles(TS ts, PetscReal t, Vec U, Vec R, void *ctx) 279d0c080abSJoseph Pusztay { 280d0c080abSJoseph Pusztay const PetscScalar *u; 281d0c080abSJoseph Pusztay PetscSection coordSection; 282d0c080abSJoseph Pusztay Vec coordsLocal; 283d0c080abSJoseph Pusztay PetscScalar *r, *coords; 284d0c080abSJoseph 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; 285d0c080abSJoseph Pusztay PetscInt dim, d, Np, Ncp, p, cStart, cEnd, c; 286d0c080abSJoseph Pusztay DM dmSw, plex; 287d0c080abSJoseph Pusztay 288d0c080abSJoseph Pusztay PetscFunctionBeginUser; 289*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(U, &Np)); 290*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(U, &u)); 291*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(R, &r)); 292*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts, &dmSw)); 293*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetCellDM(dmSw, &plex)); 294*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dmSw, &dim)); 295*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocal(plex, &coordsLocal)); 296*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateSection(plex, &coordSection)); 297*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd)); 298d0c080abSJoseph Pusztay Np /= dim; 299d0c080abSJoseph Pusztay Ncp = Np / (cEnd - cStart); 300d0c080abSJoseph Pusztay /* Calculate moments of particle distribution, note that velocity is in the coordinate */ 301*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords)); 302d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 303d0c080abSJoseph Pusztay PetscReal m1 = 0.0, m2 = 0.0; 304d0c080abSJoseph Pusztay 305d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) {m1 += PetscRealPart(coords[p*dim+d]); m2 += PetscSqr(PetscRealPart(coords[p*dim+d]));} 306d0c080abSJoseph Pusztay n += u[p]; 307d0c080abSJoseph Pusztay v += u[p]*m1; 308d0c080abSJoseph Pusztay E += u[p]*m2; 309d0c080abSJoseph Pusztay } 310d0c080abSJoseph Pusztay v /= n; 311d0c080abSJoseph Pusztay T = E*m/n - v*v; 312*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(ts, "Time %.2f: mass %.4f velocity: %+.4f temperature: %.4f\n", t, n, v, T)); 313*5f80ce2aSJacob Faibussowitsch CHKERRQ(CheckDistribution(plex, m, n, T, &v)); 314d0c080abSJoseph Pusztay /* 315d0c080abSJoseph Pusztay Begin cellwise evaluation of the collision operator. Essentially, penalize the weights of the particles 316d0c080abSJoseph Pusztay in that cell against the maxwellian for the number of particles expected to be in that cell 317d0c080abSJoseph Pusztay */ 318d0c080abSJoseph Pusztay for (c = cStart; c < cEnd; ++c) { 319d0c080abSJoseph Pusztay PetscScalar *vcoords = NULL; 320d0c080abSJoseph Pusztay PetscReal relaxation = 1.0, neq; 321d0c080abSJoseph Pusztay PetscInt sp = c*Ncp, q; 322d0c080abSJoseph Pusztay 323d0c080abSJoseph Pusztay /* Calculate equilibrium occupation for this velocity cell */ 324*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecGetClosure(plex, coordSection, coordsLocal, c, NULL, &vcoords)); 325d0c080abSJoseph Pusztay neq = ComputeCDF(m, n, T, vcoords[0], vcoords[1]); 326*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecRestoreClosure(plex, coordSection, coordsLocal, c, NULL, &vcoords)); 327d0c080abSJoseph Pusztay for (q = 0; q < Ncp; ++q) r[sp+q] = (1.0/relaxation)*(neq - u[sp+q]); 328d0c080abSJoseph Pusztay } 329d0c080abSJoseph Pusztay /* Check update */ 330d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 331d0c080abSJoseph Pusztay PetscReal m1 = 0.0, m2 = 0.0; 332d0c080abSJoseph Pusztay PetscScalar *vcoords = NULL; 333d0c080abSJoseph Pusztay 334d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) {m1 += PetscRealPart(coords[p*dim+d]); m2 += PetscSqr(PetscRealPart(coords[p*dim+d]));} 335d0c080abSJoseph Pusztay cn += r[p]; 336d0c080abSJoseph Pusztay cv += r[p]*m1; 337d0c080abSJoseph Pusztay cE += r[p]*m2; 338d0c080abSJoseph Pusztay pE += u[p]*m2; 339*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecGetClosure(plex, coordSection, coordsLocal, p, NULL, &vcoords)); 340d0c080abSJoseph Pusztay eqE += ComputeCDF(m, n, T, vcoords[0], vcoords[1])*m2; 341*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecRestoreClosure(plex, coordSection, coordsLocal, p, NULL, &vcoords)); 342d0c080abSJoseph Pusztay } 343*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(ts, "Time %.2f: mass update %.8f velocity update: %+.8f energy update: %.8f (%.8f, %.8f)\n", t, cn, cv, cE, pE, eqE)); 344*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords)); 345*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(U, &u)); 346*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(R, &r)); 347d0c080abSJoseph Pusztay PetscFunctionReturn(0); 348d0c080abSJoseph Pusztay } 349d0c080abSJoseph Pusztay 350d0c080abSJoseph Pusztay static PetscErrorCode HGMonitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) 351d0c080abSJoseph Pusztay { 352d0c080abSJoseph Pusztay AppCtx *user = (AppCtx *) ctx; 353d0c080abSJoseph Pusztay const PetscScalar *u; 354d0c080abSJoseph Pusztay DM sw, dm; 355d0c080abSJoseph Pusztay PetscInt dim, Np, p; 356d0c080abSJoseph Pusztay 357d0c080abSJoseph Pusztay PetscFunctionBeginUser; 358d0c080abSJoseph Pusztay if (step < 0) PetscFunctionReturn(0); 359d0c080abSJoseph Pusztay if (((user->ostep > 0) && (!(step % user->ostep)))) { 360d0c080abSJoseph Pusztay PetscDrawAxis axis; 361d0c080abSJoseph Pusztay 362*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts, &sw)); 363*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetCellDM(sw, &dm)); 364*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 365*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawHGReset(user->drawhg)); 366*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawHGGetAxis(user->drawhg,&axis)); 367*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawAxisSetLabels(axis,"Particles","V","f(V)")); 368*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawAxisSetLimits(axis, -3, 3, 0, 100)); 369*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawHGSetLimits(user->drawhg,-3.0, 3.0, 0, 10)); 370*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(U, &Np)); 371d0c080abSJoseph Pusztay Np /= dim; 372*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(U,&u)); 373d0c080abSJoseph Pusztay /* get points from solution vector */ 374*5f80ce2aSJacob Faibussowitsch for (p = 0; p < Np; ++p) CHKERRQ(PetscDrawHGAddValue(user->drawhg,u[p])); 375*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawHGDraw(user->drawhg)); 376*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(U,&u)); 377d0c080abSJoseph Pusztay } 378d0c080abSJoseph Pusztay PetscFunctionReturn(0); 379d0c080abSJoseph Pusztay } 380d0c080abSJoseph Pusztay 381d0c080abSJoseph Pusztay static PetscErrorCode SPMonitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) 382d0c080abSJoseph Pusztay { 383d0c080abSJoseph Pusztay AppCtx *user = (AppCtx *) ctx; 384d0c080abSJoseph Pusztay const PetscScalar *u; 385d0c080abSJoseph Pusztay PetscReal *v, *coords; 386d0c080abSJoseph Pusztay PetscInt Np, p; 387d0c080abSJoseph Pusztay DM dmSw; 388d0c080abSJoseph Pusztay 389d0c080abSJoseph Pusztay PetscFunctionBeginUser; 390d0c080abSJoseph Pusztay 391d0c080abSJoseph Pusztay if (step < 0) PetscFunctionReturn(0); 392d0c080abSJoseph Pusztay if (((user->ostep > 0) && (!(step % user->ostep)))) { 393d0c080abSJoseph Pusztay PetscDrawAxis axis; 394d0c080abSJoseph Pusztay 395*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts, &dmSw)); 396*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSPReset(user->drawsp)); 397*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSPGetAxis(user->drawsp,&axis)); 398*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawAxisSetLabels(axis,"Particles","V","w")); 399*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(U, &Np)); 400*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(U,&u)); 401d0c080abSJoseph Pusztay /* get points from solution vector */ 402*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(Np, &v)); 403d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) v[p] = PetscRealPart(u[p]); 404*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords)); 405*5f80ce2aSJacob Faibussowitsch for (p = 0; p < Np-1; ++p) CHKERRQ(PetscDrawSPAddPoint(user->drawsp, &coords[p], &v[p])); 406*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSPDraw(user->drawsp, PETSC_TRUE)); 407*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(U,&u)); 408*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords)); 409*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(v)); 410d0c080abSJoseph Pusztay } 411d0c080abSJoseph Pusztay PetscFunctionReturn(0); 412d0c080abSJoseph Pusztay } 413d0c080abSJoseph Pusztay 414d0c080abSJoseph Pusztay static PetscErrorCode KSConv(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) 415d0c080abSJoseph Pusztay { 416d0c080abSJoseph Pusztay AppCtx *user = (AppCtx *) ctx; 417d0c080abSJoseph Pusztay const PetscScalar *u; 418d0c080abSJoseph Pusztay PetscScalar sup; 419d0c080abSJoseph Pusztay PetscReal *v, *coords, T=0., vel=0., step_cast, w_sum; 420d0c080abSJoseph Pusztay PetscInt dim, Np, p, cStart, cEnd; 421d0c080abSJoseph Pusztay DM sw, plex; 422d0c080abSJoseph Pusztay 423d0c080abSJoseph Pusztay PetscFunctionBeginUser; 424d0c080abSJoseph Pusztay if (step < 0) PetscFunctionReturn(0); 425d0c080abSJoseph Pusztay if (((user->ostep > 0) && (!(step % user->ostep)))) { 426d0c080abSJoseph Pusztay PetscDrawAxis axis; 427*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSPGetAxis(user->drawks,&axis)); 428*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawAxisSetLabels(axis,"Particles","t","D_n")); 429*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSPSetLimits(user->drawks,0.,100,0.,3.5)); 430*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts, &sw)); 431*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(U, &Np)); 432*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(U,&u)); 433d0c080abSJoseph Pusztay /* get points from solution vector */ 434*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(Np, &v)); 435*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetCellDM(sw, &plex)); 436*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(sw, &dim)); 437*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd)); 438d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) v[p] = PetscRealPart(u[p]); 439*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords)); 440d0c080abSJoseph Pusztay w_sum = 0.; 441d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 442d0c080abSJoseph Pusztay w_sum += u[p]; 443d0c080abSJoseph Pusztay T += u[p]*coords[p]*coords[p]; 444d0c080abSJoseph Pusztay vel += u[p]*coords[p]; 445d0c080abSJoseph Pusztay } 446d0c080abSJoseph Pusztay vel /= w_sum; 447d0c080abSJoseph Pusztay T = T/w_sum - vel*vel; 448d0c080abSJoseph Pusztay sup = 0.0; 449d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 450d0c080abSJoseph Pusztay PetscReal tmp = 0.; 451d0c080abSJoseph Pusztay tmp = PetscAbs(u[p]-ComputePDF(1.0, w_sum, T, &coords[p*dim])); 452d0c080abSJoseph Pusztay if (tmp > sup) sup = tmp; 453d0c080abSJoseph Pusztay } 454d0c080abSJoseph Pusztay step_cast = (PetscReal)step; 455*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSPAddPoint(user->drawks, &step_cast, &sup)); 456*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSPDraw(user->drawks, PETSC_TRUE)); 457*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(U,&u)); 458*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords)); 459*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(v)); 460d0c080abSJoseph Pusztay } 461d0c080abSJoseph Pusztay PetscFunctionReturn(0); 462d0c080abSJoseph Pusztay } 463d0c080abSJoseph Pusztay 464d0c080abSJoseph Pusztay static PetscErrorCode InitializeSolve(TS ts, Vec u) 465d0c080abSJoseph Pusztay { 466d0c080abSJoseph Pusztay DM dm; 467d0c080abSJoseph Pusztay AppCtx *user; 468d0c080abSJoseph Pusztay 469d0c080abSJoseph Pusztay PetscFunctionBeginUser; 470*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts, &dm)); 471*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetApplicationContext(dm, &user)); 472*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetInitialCoordinates(dm)); 473*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetInitialConditions(dm, u)); 474d0c080abSJoseph Pusztay PetscFunctionReturn(0); 475d0c080abSJoseph Pusztay } 476d0c080abSJoseph Pusztay /* 477d0c080abSJoseph Pusztay A single particle is generated for each velocity space cell using the dmswarmpicfield_coor and is used to evaluate collisions in that cell. 478d0c080abSJoseph Pusztay 0 weight ghost particles are initialized outside of a small velocity domain to ensure the tails of the amxwellian are resolved. 479d0c080abSJoseph Pusztay */ 480d0c080abSJoseph Pusztay int main(int argc,char **argv) 481d0c080abSJoseph Pusztay { 482d0c080abSJoseph Pusztay TS ts; /* nonlinear solver */ 483d0c080abSJoseph Pusztay DM dm, sw; /* Velocity space mesh and Particle Swarm */ 484d0c080abSJoseph Pusztay Vec u, w; /* swarm vector */ 485d0c080abSJoseph Pusztay MPI_Comm comm; 486d0c080abSJoseph Pusztay AppCtx user; 487d0c080abSJoseph Pusztay PetscErrorCode ierr; 488d0c080abSJoseph Pusztay 489d0c080abSJoseph Pusztay ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 490d0c080abSJoseph Pusztay comm = PETSC_COMM_WORLD; 491*5f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(comm, &user)); 492*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(comm, &dm, &user)); 493*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateParticles(dm, &sw, &user)); 494*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetApplicationContext(sw, &user)); 495*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(comm, &ts)); 496*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetDM(ts, sw)); 497*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts, 10.0)); 498*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts, 0.01)); 499*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxSteps(ts, 100000)); 500*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 501d0c080abSJoseph Pusztay if (user.monitorhg) { 502*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawCreate(comm, NULL, "monitor", 0,0,400,300, &user.draw)); 503*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSetFromOptions(user.draw)); 504*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawHGCreate(user.draw, 20, &user.drawhg)); 505*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawHGSetColor(user.drawhg,3)); 506*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSMonitorSet(ts, HGMonitor, &user, NULL)); 507d0c080abSJoseph Pusztay } 508d0c080abSJoseph Pusztay else if (user.monitorsp) { 509*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawCreate(comm, NULL, "monitor", 0,0,400,300, &user.draw)); 510*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSetFromOptions(user.draw)); 511*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSPCreate(user.draw, 1, &user.drawsp)); 512*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSMonitorSet(ts, SPMonitor, &user, NULL)); 513d0c080abSJoseph Pusztay } 514d0c080abSJoseph Pusztay else if (user.monitorks) { 515*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawCreate(comm, NULL, "monitor", 0,0,400,300, &user.draw)); 516*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSetFromOptions(user.draw)); 517*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSPCreate(user.draw, 1, &user.drawks)); 518*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSMonitorSet(ts, KSConv, &user, NULL)); 519d0c080abSJoseph Pusztay } 520*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts, NULL, RHSFunctionParticles, &user)); 521*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 522*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetComputeInitialCondition(ts, InitializeSolve)); 523*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &w)); 524*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(w, &u)); 525*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(w, u)); 526*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &w)); 527*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSComputeInitialCondition(ts, u)); 528*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts, u)); 529d0c080abSJoseph Pusztay if (user.monitorhg) { 530*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSave(user.draw)); 531*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawHGDestroy(&user.drawhg)); 532*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawDestroy(&user.draw)); 533d0c080abSJoseph Pusztay } 534d0c080abSJoseph Pusztay if (user.monitorsp) { 535*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSave(user.draw)); 536*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSPDestroy(&user.drawsp)); 537*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawDestroy(&user.draw)); 538d0c080abSJoseph Pusztay } 539d0c080abSJoseph Pusztay if (user.monitorks) { 540*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSave(user.draw)); 541*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSPDestroy(&user.drawks)); 542*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawDestroy(&user.draw)); 543d0c080abSJoseph Pusztay } 544*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 545*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 546*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&sw)); 547*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 548d0c080abSJoseph Pusztay ierr = PetscFinalize(); 549d0c080abSJoseph Pusztay return ierr; 550d0c080abSJoseph Pusztay } 551d0c080abSJoseph Pusztay 552d0c080abSJoseph Pusztay /*TEST 553d0c080abSJoseph Pusztay build: 55430602db0SMatthew G. Knepley requires: double !complex 555d0c080abSJoseph Pusztay test: 556d0c080abSJoseph Pusztay suffix: 1 55730602db0SMatthew 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 558d0c080abSJoseph Pusztay test: 559d0c080abSJoseph Pusztay suffix: 2 56030602db0SMatthew 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 561d0c080abSJoseph Pusztay TEST*/ 562