1d0c080abSJoseph Pusztay static char help[] = "Particle Basis Landau Example using nonlinear solve + Implicit Midpoint-like time stepping."; 2d0c080abSJoseph Pusztay 3d0c080abSJoseph Pusztay /* TODO 4d0c080abSJoseph Pusztay 5d0c080abSJoseph Pusztay 1) SNES is sensitive to epsilon (but not to h). Should we do continuation in it? 6d0c080abSJoseph Pusztay 7d0c080abSJoseph Pusztay 2) Put this timestepper in library, maybe by changing DG 8d0c080abSJoseph Pusztay 9d0c080abSJoseph Pusztay 3) Add monitor to visualize distributions 10d0c080abSJoseph Pusztay 11d0c080abSJoseph Pusztay */ 12d0c080abSJoseph Pusztay 13d0c080abSJoseph Pusztay /* References 14d0c080abSJoseph Pusztay [1] https://arxiv.org/abs/1910.03080v2 15d0c080abSJoseph Pusztay */ 16d0c080abSJoseph Pusztay 17d0c080abSJoseph Pusztay #include <petscdmplex.h> 18d0c080abSJoseph Pusztay #include <petscdmswarm.h> 19d0c080abSJoseph Pusztay #include <petscts.h> 20d0c080abSJoseph Pusztay #include <petscviewer.h> 21d0c080abSJoseph Pusztay #include <petscmath.h> 22d0c080abSJoseph Pusztay 23d0c080abSJoseph Pusztay typedef struct { 24d0c080abSJoseph Pusztay /* Velocity space grid and functions */ 25d0c080abSJoseph Pusztay PetscInt N; /* The number of partices per spatial cell */ 26d0c080abSJoseph Pusztay PetscReal L; /* Velocity space is [-L, L]^d */ 27d0c080abSJoseph Pusztay PetscReal h; /* Spacing for grid 2L / N^{1/d} */ 28d0c080abSJoseph Pusztay PetscReal epsilon; /* gaussian regularization parameter */ 29d0c080abSJoseph Pusztay PetscReal momentTol; /* Tolerance for checking moment conservation */ 30d0c080abSJoseph Pusztay } AppCtx; 31d0c080abSJoseph Pusztay 329371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) { 33d0c080abSJoseph Pusztay PetscFunctionBeginUser; 34d0c080abSJoseph Pusztay options->N = 1; 35d0c080abSJoseph Pusztay options->momentTol = 100.0 * PETSC_MACHINE_EPSILON; 36d0c080abSJoseph Pusztay options->L = 1.0; 37d0c080abSJoseph Pusztay options->h = -1.0; 38d0c080abSJoseph Pusztay options->epsilon = -1.0; 39d0c080abSJoseph Pusztay 40d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Collision Options", "DMPLEX"); 419566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-N", "Number of particles per spatial cell", "ex27.c", options->N, &options->N, NULL)); 429566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-L", "Velocity-space extent", "ex27.c", options->L, &options->L, NULL)); 439566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-h", "Velocity-space resolution", "ex27.c", options->h, &options->h, NULL)); 449566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-epsilon", "Mollifier regularization parameter", "ex27.c", options->epsilon, &options->epsilon, NULL)); 45d0609cedSBarry Smith PetscOptionsEnd(); 46d0c080abSJoseph Pusztay PetscFunctionReturn(0); 47d0c080abSJoseph Pusztay } 48d0c080abSJoseph Pusztay 499371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user) { 50d0c080abSJoseph Pusztay PetscFunctionBeginUser; 519566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 529566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 539566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 549566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 55d0c080abSJoseph Pusztay PetscFunctionReturn(0); 56d0c080abSJoseph Pusztay } 57d0c080abSJoseph Pusztay 589371c9d4SSatish Balay static PetscErrorCode SetInitialCoordinates(DM sw) { 59d0c080abSJoseph Pusztay AppCtx *user; 60d0c080abSJoseph Pusztay PetscRandom rnd, rndv; 61d0c080abSJoseph Pusztay DM dm; 62d0c080abSJoseph Pusztay DMPolytopeType ct; 63d0c080abSJoseph Pusztay PetscBool simplex; 64d0c080abSJoseph Pusztay PetscReal *centroid, *coords, *velocity, *xi0, *v0, *J, *invJ, detJ, *vals; 65d0c080abSJoseph Pusztay PetscInt dim, d, cStart, cEnd, c, Np, p; 66d0c080abSJoseph Pusztay 67d0c080abSJoseph Pusztay PetscFunctionBeginUser; 68d0c080abSJoseph Pusztay /* Randomization for coordinates */ 699566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)sw), &rnd)); 709566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0)); 719566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rnd)); 729566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)sw), &rndv)); 739566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rndv, -1., 1.)); 749566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rndv)); 759566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(sw, &user)); 76d0c080abSJoseph Pusztay Np = user->N; 779566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 789566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(sw, &dm)); 799566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 809566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 81d0c080abSJoseph Pusztay simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; 829566063dSJacob Faibussowitsch PetscCall(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; 849566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 859566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&velocity)); 869566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&vals)); 87d0c080abSJoseph Pusztay for (c = cStart; c < cEnd; ++c) { 88d0c080abSJoseph Pusztay if (Np == 1) { 899566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL)); 909371c9d4SSatish Balay for (d = 0; d < dim; ++d) { coords[c * dim + d] = centroid[d]; } 91d0c080abSJoseph Pusztay vals[c] = 1.0; 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 } 1029371c9d4SSatish Balay if (simplex && sum > 0.0) 1039371c9d4SSatish 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 } 109d0c080abSJoseph Pusztay /* Random velocity IC */ 110d0c080abSJoseph Pusztay for (c = cStart; c < cEnd; ++c) { 111d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 112d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) { 113d0c080abSJoseph Pusztay PetscReal v_val; 114d0c080abSJoseph Pusztay 1159566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(rndv, &v_val)); 116d0c080abSJoseph Pusztay velocity[p * dim + d] = v_val; 117d0c080abSJoseph Pusztay } 118d0c080abSJoseph Pusztay } 119d0c080abSJoseph Pusztay } 1209566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 1219566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&velocity)); 1229566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&vals)); 1239566063dSJacob Faibussowitsch PetscCall(PetscFree5(centroid, xi0, v0, J, invJ)); 1249566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rnd)); 1259566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rndv)); 126d0c080abSJoseph Pusztay PetscFunctionReturn(0); 127d0c080abSJoseph Pusztay } 128d0c080abSJoseph Pusztay 129d0c080abSJoseph Pusztay /* Get velocities from swarm and place in solution vector */ 1309371c9d4SSatish Balay static PetscErrorCode SetInitialConditions(DM dmSw, Vec u) { 131d0c080abSJoseph Pusztay DM dm; 132d0c080abSJoseph Pusztay AppCtx *user; 133d0c080abSJoseph Pusztay PetscReal *velocity; 134d0c080abSJoseph Pusztay PetscScalar *initialConditions; 135d0c080abSJoseph Pusztay PetscInt dim, d, cStart, cEnd, c, Np, p, n; 136d0c080abSJoseph Pusztay 137d0c080abSJoseph Pusztay PetscFunctionBeginUser; 1389566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(u, &n)); 1399566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(dmSw, &user)); 140d0c080abSJoseph Pusztay Np = user->N; 1419566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dmSw, &dm)); 1429566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1439566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 1449566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dmSw, "velocity", NULL, NULL, (void **)&velocity)); 1459566063dSJacob Faibussowitsch PetscCall(VecGetArray(u, &initialConditions)); 146d0c080abSJoseph Pusztay for (c = cStart; c < cEnd; ++c) { 147d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 148d0c080abSJoseph Pusztay const PetscInt n = c * Np + p; 1499371c9d4SSatish Balay for (d = 0; d < dim; d++) { initialConditions[n * dim + d] = velocity[n * dim + d]; } 150d0c080abSJoseph Pusztay } 151d0c080abSJoseph Pusztay } 1529566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(u, &initialConditions)); 1539566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dmSw, "velocity", NULL, NULL, (void **)&velocity)); 154d0c080abSJoseph Pusztay PetscFunctionReturn(0); 155d0c080abSJoseph Pusztay } 156d0c080abSJoseph Pusztay 1579371c9d4SSatish Balay static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user) { 158d0c080abSJoseph Pusztay PetscInt *cellid; 159d0c080abSJoseph Pusztay PetscInt dim, cStart, cEnd, c, Np = user->N, p; 160d0c080abSJoseph Pusztay PetscBool view = PETSC_FALSE; 161d0c080abSJoseph Pusztay 162d0c080abSJoseph Pusztay PetscFunctionBeginUser; 1639566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1649566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw)); 1659566063dSJacob Faibussowitsch PetscCall(DMSetType(*sw, DMSWARM)); 1669566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*sw, dim)); 167d0c080abSJoseph Pusztay /* h = 2L/n and N = n^d */ 168d0c080abSJoseph Pusztay if (user->h < 0.) user->h = 2. * user->L / PetscPowReal(user->N, 1. / dim); 169d0c080abSJoseph Pusztay /* From Section 4 in [1], \epsilon = 0.64 h^.98 */ 170d0c080abSJoseph Pusztay if (user->epsilon < 0.) user->epsilon = 0.64 * pow(user->h, 1.98); 1719566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-param_view", &view, NULL)); 172*48a46eb9SPierre Jolivet if (view) PetscCall(PetscPrintf(PETSC_COMM_SELF, "N: %" PetscInt_FMT " L: %g h: %g eps: %g\n", user->N, (double)user->L, (double)user->h, (double)user->epsilon)); 1739566063dSJacob Faibussowitsch PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC)); 1749566063dSJacob Faibussowitsch PetscCall(DMSwarmSetCellDM(*sw, dm)); 1759566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "velocity", dim, PETSC_REAL)); 1769566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR)); 1779566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(*sw)); 1789566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 1799566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0)); 1809566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*sw)); 1819566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 182d0c080abSJoseph Pusztay for (c = cStart; c < cEnd; ++c) { 183d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 184d0c080abSJoseph Pusztay const PetscInt n = c * Np + p; 185d0c080abSJoseph Pusztay cellid[n] = c; 186d0c080abSJoseph Pusztay } 187d0c080abSJoseph Pusztay } 1889566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 1899566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles")); 1909566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view")); 191d0c080abSJoseph Pusztay PetscFunctionReturn(0); 192d0c080abSJoseph Pusztay } 193d0c080abSJoseph Pusztay 194d0c080abSJoseph Pusztay /* Internal dmplex function, same as found in dmpleximpl.h */ 1959371c9d4SSatish Balay static void DMPlex_WaxpyD_Internal(PetscInt dim, PetscReal a, const PetscReal *x, const PetscReal *y, PetscReal *w) { 196d0c080abSJoseph Pusztay PetscInt d; 197d0c080abSJoseph Pusztay 198d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) w[d] = a * x[d] + y[d]; 199d0c080abSJoseph Pusztay } 200d0c080abSJoseph Pusztay 201d0c080abSJoseph Pusztay /* Internal dmplex function, same as found in dmpleximpl.h */ 2029371c9d4SSatish Balay static PetscReal DMPlex_DotD_Internal(PetscInt dim, const PetscScalar *x, const PetscReal *y) { 203d0c080abSJoseph Pusztay PetscReal sum = 0.0; 204d0c080abSJoseph Pusztay PetscInt d; 205d0c080abSJoseph Pusztay 206d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) sum += PetscRealPart(x[d]) * y[d]; 207d0c080abSJoseph Pusztay return sum; 208d0c080abSJoseph Pusztay } 209d0c080abSJoseph Pusztay 210d0c080abSJoseph Pusztay /* Internal dmplex function, same as found in dmpleximpl.h */ 2119371c9d4SSatish Balay static void DMPlex_MultAdd2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[]) { 212d0c080abSJoseph Pusztay PetscScalar z[2]; 2139371c9d4SSatish Balay z[0] = x[0]; 2149371c9d4SSatish Balay z[1] = x[ldx]; 215d0c080abSJoseph Pusztay y[0] += A[0] * z[0] + A[1] * z[1]; 216d0c080abSJoseph Pusztay y[ldx] += A[2] * z[0] + A[3] * z[1]; 217d0c080abSJoseph Pusztay (void)PetscLogFlops(6.0); 218d0c080abSJoseph Pusztay } 219d0c080abSJoseph Pusztay 220d0c080abSJoseph Pusztay /* Internal dmplex function, same as found in dmpleximpl.h to avoid private includes. */ 2219371c9d4SSatish Balay static void DMPlex_MultAdd3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[]) { 222d0c080abSJoseph Pusztay PetscScalar z[3]; 2239371c9d4SSatish Balay z[0] = x[0]; 2249371c9d4SSatish Balay z[1] = x[ldx]; 2259371c9d4SSatish Balay z[2] = x[ldx * 2]; 226d0c080abSJoseph Pusztay y[0] += A[0] * z[0] + A[1] * z[1] + A[2] * z[2]; 227d0c080abSJoseph Pusztay y[ldx] += A[3] * z[0] + A[4] * z[1] + A[5] * z[2]; 228d0c080abSJoseph Pusztay y[ldx * 2] += A[6] * z[0] + A[7] * z[1] + A[8] * z[2]; 229d0c080abSJoseph Pusztay (void)PetscLogFlops(15.0); 230d0c080abSJoseph Pusztay } 231d0c080abSJoseph Pusztay 232d0c080abSJoseph Pusztay /* 233d0c080abSJoseph Pusztay Gaussian - The Gaussian function G(x) 234d0c080abSJoseph Pusztay 235d0c080abSJoseph Pusztay Input Parameters: 236d0c080abSJoseph Pusztay + dim - The number of dimensions, or size of x 237d0c080abSJoseph Pusztay . mu - The mean, or center 238d0c080abSJoseph Pusztay . sigma - The standard deviation, or width 239d0c080abSJoseph Pusztay - x - The evaluation point of the function 240d0c080abSJoseph Pusztay 241d0c080abSJoseph Pusztay Output Parameter: 242d0c080abSJoseph Pusztay . ret - The value G(x) 243d0c080abSJoseph Pusztay */ 2449371c9d4SSatish Balay static PetscReal Gaussian(PetscInt dim, const PetscReal mu[], PetscReal sigma, const PetscReal x[]) { 245d0c080abSJoseph Pusztay PetscReal arg = 0.0; 246d0c080abSJoseph Pusztay PetscInt d; 247d0c080abSJoseph Pusztay 248d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) arg += PetscSqr(x[d] - mu[d]); 249d0c080abSJoseph Pusztay return PetscPowReal(2.0 * PETSC_PI * sigma, -dim / 2.0) * PetscExpReal(-arg / (2.0 * sigma)); 250d0c080abSJoseph Pusztay } 251d0c080abSJoseph Pusztay 252d0c080abSJoseph Pusztay /* 253d0c080abSJoseph Pusztay ComputeGradS - Compute grad_v dS_eps/df 254d0c080abSJoseph Pusztay 255d0c080abSJoseph Pusztay Input Parameters: 256d0c080abSJoseph Pusztay + dim - The dimension 257d0c080abSJoseph Pusztay . Np - The number of particles 258d0c080abSJoseph Pusztay . vp - The velocity v_p of the particle at which we evaluate 259d0c080abSJoseph Pusztay . velocity - The velocity field for all particles 260d0c080abSJoseph Pusztay . epsilon - The regularization strength 261d0c080abSJoseph Pusztay 262d0c080abSJoseph Pusztay Output Parameter: 263d0c080abSJoseph Pusztay . integral - The output grad_v dS_eps/df (v_p) 264d0c080abSJoseph Pusztay 265d0c080abSJoseph Pusztay Note: 266d0c080abSJoseph Pusztay This comes from (3.6) in [1], and we are computing 267d0c080abSJoseph Pusztay $ \nabla_v S_p = \grad \psi_\epsilon(v_p - v) log \sum_q \psi_\epsilon(v - v_q) 268d0c080abSJoseph Pusztay which is discretized by using a one-point quadrature in each box l at its center v^c_l 269d0c080abSJoseph Pusztay $ \sum_l h^d \nabla\psi_\epsilon(v_p - v^c_l) \log\left( \sum_q w_q \psi_\epsilon(v^c_l - v_q) \right) 270d0c080abSJoseph Pusztay where h^d is the volume of each box. 271d0c080abSJoseph Pusztay */ 2729371c9d4SSatish Balay static PetscErrorCode ComputeGradS(PetscInt dim, PetscInt Np, const PetscReal vp[], const PetscReal velocity[], PetscReal integral[], AppCtx *ctx) { 273d0c080abSJoseph Pusztay PetscReal vc_l[3], L = ctx->L, h = ctx->h, epsilon = ctx->epsilon, init = 0.5 * h - L; 274d0c080abSJoseph Pusztay PetscInt nx = roundf(2. * L / h); 275d0c080abSJoseph Pusztay PetscInt ny = dim > 1 ? nx : 1; 276d0c080abSJoseph Pusztay PetscInt nz = dim > 2 ? nx : 1; 277d0c080abSJoseph Pusztay PetscInt i, j, k, d, q, dbg = 0; 278d0c080abSJoseph Pusztay 279d0c080abSJoseph Pusztay PetscFunctionBeginHot; 280d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) integral[d] = 0.0; 281d0c080abSJoseph Pusztay for (k = 0, vc_l[2] = init; k < nz; ++k, vc_l[2] += h) { 282d0c080abSJoseph Pusztay for (j = 0, vc_l[1] = init; j < ny; ++j, vc_l[1] += h) { 283d0c080abSJoseph Pusztay for (i = 0, vc_l[0] = init; i < nx; ++i, vc_l[0] += h) { 284d0c080abSJoseph Pusztay PetscReal sum = 0.0; 285d0c080abSJoseph Pusztay 28663a3b9bcSJacob Faibussowitsch if (dbg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "(%" PetscInt_FMT " %" PetscInt_FMT ") vc_l: %g %g\n", i, j, (double)vc_l[0], (double)vc_l[1])); 287d0c080abSJoseph Pusztay /* \log \sum_k \psi(v - v_k) */ 288d0c080abSJoseph Pusztay for (q = 0; q < Np; ++q) sum += Gaussian(dim, &velocity[q * dim], epsilon, vc_l); 289d0c080abSJoseph Pusztay sum = PetscLogReal(sum); 290d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) integral[d] += (-1. / (epsilon)) * PetscAbsReal(vp[d] - vc_l[d]) * (Gaussian(dim, vp, epsilon, vc_l)) * sum; 291d0c080abSJoseph Pusztay } 292d0c080abSJoseph Pusztay } 293d0c080abSJoseph Pusztay } 294d0c080abSJoseph Pusztay PetscFunctionReturn(0); 295d0c080abSJoseph Pusztay } 296d0c080abSJoseph Pusztay 297d0c080abSJoseph Pusztay /* Q = 1/|xi| (I - xi xi^T / |xi|^2), xi = vp - vq */ 2989371c9d4SSatish Balay static PetscErrorCode QCompute(PetscInt dim, const PetscReal vp[], const PetscReal vq[], PetscReal Q[]) { 299d0c080abSJoseph Pusztay PetscReal xi[3], xi2, xi3, mag; 300d0c080abSJoseph Pusztay PetscInt d, e; 301d0c080abSJoseph Pusztay 302d0c080abSJoseph Pusztay PetscFunctionBeginHot; 303d0c080abSJoseph Pusztay DMPlex_WaxpyD_Internal(dim, -1.0, vq, vp, xi); 304d0c080abSJoseph Pusztay xi2 = DMPlex_DotD_Internal(dim, xi, xi); 305d0c080abSJoseph Pusztay mag = PetscSqrtReal(xi2); 306d0c080abSJoseph Pusztay xi3 = xi2 * mag; 307d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) { 3089371c9d4SSatish Balay for (e = 0; e < dim; ++e) { Q[d * dim + e] = -xi[d] * xi[e] / xi3; } 309d0c080abSJoseph Pusztay Q[d * dim + d] += 1. / mag; 310d0c080abSJoseph Pusztay } 311d0c080abSJoseph Pusztay PetscFunctionReturn(0); 312d0c080abSJoseph Pusztay } 313d0c080abSJoseph Pusztay 3149371c9d4SSatish Balay static PetscErrorCode RHSFunctionParticles(TS ts, PetscReal t, Vec U, Vec R, void *ctx) { 315d0c080abSJoseph Pusztay AppCtx *user = (AppCtx *)ctx; 316d0c080abSJoseph Pusztay PetscInt dbg = 0; 317d0c080abSJoseph Pusztay DM sw; /* Particles */ 318d0c080abSJoseph Pusztay Vec sol; /* Solution vector at current time */ 319d0c080abSJoseph Pusztay const PetscScalar *u; /* input solution vector */ 320d0c080abSJoseph Pusztay PetscScalar *r; 321d0c080abSJoseph Pusztay PetscReal *velocity; 322d0c080abSJoseph Pusztay PetscInt dim, Np, p, q; 323d0c080abSJoseph Pusztay 324d0c080abSJoseph Pusztay PetscFunctionBeginUser; 3259566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(R)); 3269566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &sw)); 3279566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 3289566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(U, &Np)); 3299566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &sol)); 3309566063dSJacob Faibussowitsch PetscCall(VecGetArray(sol, &velocity)); 3319566063dSJacob Faibussowitsch PetscCall(VecGetArray(R, &r)); 3329566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 333d0c080abSJoseph Pusztay Np /= dim; 3349566063dSJacob Faibussowitsch if (dbg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Part ppr x y\n")); 335d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 336d0c080abSJoseph Pusztay PetscReal gradS_p[3] = {0., 0., 0.}; 337d0c080abSJoseph Pusztay 3389566063dSJacob Faibussowitsch PetscCall(ComputeGradS(dim, Np, &velocity[p * dim], velocity, gradS_p, user)); 339d0c080abSJoseph Pusztay for (q = 0; q < Np; ++q) { 340d0c080abSJoseph Pusztay PetscReal gradS_q[3] = {0., 0., 0.}, GammaS[3] = {0., 0., 0.}, Q[9]; 341d0c080abSJoseph Pusztay 342d0c080abSJoseph Pusztay if (q == p) continue; 3439566063dSJacob Faibussowitsch PetscCall(ComputeGradS(dim, Np, &velocity[q * dim], velocity, gradS_q, user)); 344d0c080abSJoseph Pusztay DMPlex_WaxpyD_Internal(dim, -1.0, gradS_q, gradS_p, GammaS); 3459566063dSJacob Faibussowitsch PetscCall(QCompute(dim, &u[p * dim], &u[q * dim], Q)); 346d0c080abSJoseph Pusztay switch (dim) { 347d0c080abSJoseph Pusztay case 2: DMPlex_MultAdd2DReal_Internal(Q, 1, GammaS, &r[p * dim]); break; 348d0c080abSJoseph Pusztay case 3: DMPlex_MultAdd3DReal_Internal(Q, 1, GammaS, &r[p * dim]); break; 34963a3b9bcSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Do not support dimension %" PetscInt_FMT, dim); 350d0c080abSJoseph Pusztay } 351d0c080abSJoseph Pusztay } 35263a3b9bcSJacob Faibussowitsch if (dbg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final %4" PetscInt_FMT " %10.8lf %10.8lf\n", p, r[p * dim + 0], r[p * dim + 1])); 353d0c080abSJoseph Pusztay } 3549566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 3559566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(R, &r)); 3569566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(sol, &velocity)); 3579566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(R, NULL, "-residual_view")); 358d0c080abSJoseph Pusztay PetscFunctionReturn(0); 359d0c080abSJoseph Pusztay } 360d0c080abSJoseph Pusztay 361d0c080abSJoseph Pusztay /* 362d0c080abSJoseph Pusztay TS Post Step Function. Copy the solution back into the swarm for migration. We may also need to reform 363d0c080abSJoseph Pusztay the solution vector in cases of particle migration, but we forgo that here since there is no velocity space grid 364d0c080abSJoseph Pusztay to migrate between. 365d0c080abSJoseph Pusztay */ 3669371c9d4SSatish Balay static PetscErrorCode UpdateSwarm(TS ts) { 367d0c080abSJoseph Pusztay PetscInt idx, n; 368d0c080abSJoseph Pusztay const PetscScalar *u; 369d0c080abSJoseph Pusztay PetscScalar *velocity; 370d0c080abSJoseph Pusztay DM sw; 371d0c080abSJoseph Pusztay Vec sol; 372d0c080abSJoseph Pusztay 373d0c080abSJoseph Pusztay PetscFunctionBeginUser; 3749566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &sw)); 3759566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&velocity)); 3769566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &sol)); 3779566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(sol, &u)); 3789566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(sol, &n)); 379d0c080abSJoseph Pusztay for (idx = 0; idx < n; ++idx) velocity[idx] = u[idx]; 3809566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(sol, &u)); 3819566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&velocity)); 382d0c080abSJoseph Pusztay PetscFunctionReturn(0); 383d0c080abSJoseph Pusztay } 384d0c080abSJoseph Pusztay 3859371c9d4SSatish Balay static PetscErrorCode InitializeSolve(TS ts, Vec u) { 386d0c080abSJoseph Pusztay DM dm; 387d0c080abSJoseph Pusztay AppCtx *user; 388d0c080abSJoseph Pusztay 389d0c080abSJoseph Pusztay PetscFunctionBeginUser; 3909566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 3919566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(dm, &user)); 3929566063dSJacob Faibussowitsch PetscCall(SetInitialCoordinates(dm)); 3939566063dSJacob Faibussowitsch PetscCall(SetInitialConditions(dm, u)); 394d0c080abSJoseph Pusztay PetscFunctionReturn(0); 395d0c080abSJoseph Pusztay } 396d0c080abSJoseph Pusztay 3979371c9d4SSatish Balay int main(int argc, char **argv) { 398d0c080abSJoseph Pusztay TS ts; /* nonlinear solver */ 399d0c080abSJoseph Pusztay DM dm, sw; /* Velocity space mesh and Particle Swarm */ 400d0c080abSJoseph Pusztay Vec u, v; /* problem vector */ 401d0c080abSJoseph Pusztay MPI_Comm comm; 402d0c080abSJoseph Pusztay AppCtx user; 403d0c080abSJoseph Pusztay 404327415f7SBarry Smith PetscFunctionBeginUser; 4059566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 406d0c080abSJoseph Pusztay comm = PETSC_COMM_WORLD; 4079566063dSJacob Faibussowitsch PetscCall(ProcessOptions(comm, &user)); 408d0c080abSJoseph Pusztay /* Initialize objects and set initial conditions */ 4099566063dSJacob Faibussowitsch PetscCall(CreateMesh(comm, &dm, &user)); 4109566063dSJacob Faibussowitsch PetscCall(CreateParticles(dm, &sw, &user)); 4119566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(sw, &user)); 4129566063dSJacob Faibussowitsch PetscCall(DMSwarmVectorDefineField(sw, "velocity")); 4139566063dSJacob Faibussowitsch PetscCall(TSCreate(comm, &ts)); 4149566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, sw)); 4159566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 10.0)); 4169566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 0.1)); 4179566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, 1)); 4189566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 4199566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunctionParticles, &user)); 4209566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 4219566063dSJacob Faibussowitsch PetscCall(TSSetComputeInitialCondition(ts, InitializeSolve)); 4229566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "velocity", &v)); 4239566063dSJacob Faibussowitsch PetscCall(VecDuplicate(v, &u)); 4249566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "velocity", &v)); 4259566063dSJacob Faibussowitsch PetscCall(TSComputeInitialCondition(ts, u)); 4269566063dSJacob Faibussowitsch PetscCall(TSSetPostStep(ts, UpdateSwarm)); 4279566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, u)); 4289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 4299566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 4309566063dSJacob Faibussowitsch PetscCall(DMDestroy(&sw)); 4319566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 4329566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 433b122ec5aSJacob Faibussowitsch return 0; 434d0c080abSJoseph Pusztay } 435d0c080abSJoseph Pusztay 436d0c080abSJoseph Pusztay /*TEST 437d0c080abSJoseph Pusztay build: 438d0c080abSJoseph Pusztay requires: triangle !single !complex 439d0c080abSJoseph Pusztay test: 440d0c080abSJoseph Pusztay suffix: midpoint 44130602db0SMatthew G. Knepley args: -N 3 -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1,-1 -dm_plex_box_upper 1,1 -dm_view \ 442d0c080abSJoseph Pusztay -ts_type theta -ts_theta_theta 0.5 -ts_dmswarm_monitor_moments -ts_monitor_frequency 1 -snes_fd 443d0c080abSJoseph Pusztay TEST*/ 444