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