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 32d0c080abSJoseph Pusztay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 33d0c080abSJoseph Pusztay { 34d0c080abSJoseph Pusztay PetscErrorCode ierr; 35d0c080abSJoseph Pusztay 36d0c080abSJoseph Pusztay PetscFunctionBeginUser; 37d0c080abSJoseph Pusztay options->N = 1; 38d0c080abSJoseph Pusztay options->momentTol = 100.0*PETSC_MACHINE_EPSILON; 39d0c080abSJoseph Pusztay options->L = 1.0; 40d0c080abSJoseph Pusztay options->h = -1.0; 41d0c080abSJoseph Pusztay options->epsilon = -1.0; 42d0c080abSJoseph Pusztay 43d0c080abSJoseph Pusztay ierr = PetscOptionsBegin(comm, "", "Collision Options", "DMPLEX");CHKERRQ(ierr); 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-N", "Number of particles per spatial cell", "ex27.c", options->N, &options->N, NULL)); 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-L", "Velocity-space extent", "ex27.c", options->L, &options->L, NULL)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-h", "Velocity-space resolution", "ex27.c", options->h, &options->h, NULL)); 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-epsilon", "Mollifier regularization parameter", "ex27.c", options->epsilon, &options->epsilon, NULL)); 48d0c080abSJoseph Pusztay ierr = PetscOptionsEnd();CHKERRQ(ierr); 49d0c080abSJoseph Pusztay PetscFunctionReturn(0); 50d0c080abSJoseph Pusztay } 51d0c080abSJoseph Pusztay 52d0c080abSJoseph Pusztay static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user) 53d0c080abSJoseph Pusztay { 54d0c080abSJoseph Pusztay PetscFunctionBeginUser; 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 56*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 59d0c080abSJoseph Pusztay PetscFunctionReturn(0); 60d0c080abSJoseph Pusztay } 61d0c080abSJoseph Pusztay 62d0c080abSJoseph Pusztay static PetscErrorCode SetInitialCoordinates(DM sw) 63d0c080abSJoseph Pusztay { 64d0c080abSJoseph Pusztay AppCtx *user; 65d0c080abSJoseph Pusztay PetscRandom rnd, rndv; 66d0c080abSJoseph Pusztay DM dm; 67d0c080abSJoseph Pusztay DMPolytopeType ct; 68d0c080abSJoseph Pusztay PetscBool simplex; 69d0c080abSJoseph Pusztay PetscReal *centroid, *coords, *velocity, *xi0, *v0, *J, *invJ, detJ, *vals; 70d0c080abSJoseph Pusztay PetscInt dim, d, cStart, cEnd, c, Np, p; 71d0c080abSJoseph Pusztay 72d0c080abSJoseph Pusztay PetscFunctionBeginUser; 73d0c080abSJoseph Pusztay /* Randomization for coordinates */ 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PetscObjectComm((PetscObject) sw), &rnd)); 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetInterval(rnd, -1.0, 1.0)); 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rnd)); 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PetscObjectComm((PetscObject) sw), &rndv)); 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetInterval(rndv, -1., 1.)); 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rndv)); 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetApplicationContext(sw, &user)); 81d0c080abSJoseph Pusztay Np = user->N; 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(sw, &dim)); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetCellDM(sw, &dm)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCellType(dm, cStart, &ct)); 86d0c080abSJoseph Pusztay simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc5(dim, ¢roid, dim, &xi0, dim, &v0, dim*dim, &J, dim*dim, &invJ)); 88d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) xi0[d] = -1.0; 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords)); 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **) &velocity)); 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **) &vals)); 92d0c080abSJoseph Pusztay for (c = cStart; c < cEnd; ++c) { 93d0c080abSJoseph Pusztay if (Np == 1) { 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL)); 95d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) { 96d0c080abSJoseph Pusztay coords[c*dim+d] = centroid[d]; 97d0c080abSJoseph Pusztay } 98d0c080abSJoseph Pusztay vals[c] = 1.0; 99d0c080abSJoseph Pusztay } else { 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); /* affine */ 101d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 102d0c080abSJoseph Pusztay const PetscInt n = c*Np + p; 103d0c080abSJoseph Pusztay PetscReal sum = 0.0, refcoords[3]; 104d0c080abSJoseph Pusztay 105d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) { 106*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValueReal(rnd, &refcoords[d])); 107d0c080abSJoseph Pusztay sum += refcoords[d]; 108d0c080abSJoseph Pusztay } 109d0c080abSJoseph Pusztay if (simplex && sum > 0.0) for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim)*sum; 110d0c080abSJoseph Pusztay vals[n] = 1.0; 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexReferenceToCoordinates(dm, c, 1, refcoords, &coords[n*dim])); 112d0c080abSJoseph Pusztay } 113d0c080abSJoseph Pusztay } 114d0c080abSJoseph Pusztay } 115d0c080abSJoseph Pusztay /* Random velocity IC */ 116d0c080abSJoseph Pusztay for (c = cStart; c < cEnd; ++c) { 117d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 118d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) { 119d0c080abSJoseph Pusztay PetscReal v_val; 120d0c080abSJoseph Pusztay 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValueReal(rndv, &v_val)); 122d0c080abSJoseph Pusztay velocity[p*dim+d] = v_val; 123d0c080abSJoseph Pusztay } 124d0c080abSJoseph Pusztay } 125d0c080abSJoseph Pusztay } 126*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords)); 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **) &velocity)); 128*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **) &vals)); 129*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree5(centroid, xi0, v0, J, invJ)); 130*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rnd)); 131*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rndv)); 132d0c080abSJoseph Pusztay PetscFunctionReturn(0); 133d0c080abSJoseph Pusztay } 134d0c080abSJoseph Pusztay 135d0c080abSJoseph Pusztay /* Get velocities from swarm and place in solution vector */ 136d0c080abSJoseph Pusztay static PetscErrorCode SetInitialConditions(DM dmSw, Vec u) 137d0c080abSJoseph Pusztay { 138d0c080abSJoseph Pusztay DM dm; 139d0c080abSJoseph Pusztay AppCtx *user; 140d0c080abSJoseph Pusztay PetscReal *velocity; 141d0c080abSJoseph Pusztay PetscScalar *initialConditions; 142d0c080abSJoseph Pusztay PetscInt dim, d, cStart, cEnd, c, Np, p, n; 143d0c080abSJoseph Pusztay 144d0c080abSJoseph Pusztay PetscFunctionBeginUser; 145*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(u, &n)); 146*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetApplicationContext(dmSw, &user)); 147d0c080abSJoseph Pusztay Np = user->N; 148*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetCellDM(dmSw, &dm)); 149*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 150*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 151*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(dmSw, "velocity", NULL, NULL, (void **) &velocity)); 152*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(u, &initialConditions)); 153d0c080abSJoseph Pusztay for (c = cStart; c < cEnd; ++c) { 154d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 155d0c080abSJoseph Pusztay const PetscInt n = c*Np + p; 156d0c080abSJoseph Pusztay for (d = 0; d < dim; d++) { 157d0c080abSJoseph Pusztay initialConditions[n*dim+d] = velocity[n*dim+d]; 158d0c080abSJoseph Pusztay } 159d0c080abSJoseph Pusztay } 160d0c080abSJoseph Pusztay } 161*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(u, &initialConditions)); 162*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(dmSw, "velocity", NULL, NULL, (void **) &velocity)); 163d0c080abSJoseph Pusztay PetscFunctionReturn(0); 164d0c080abSJoseph Pusztay } 165d0c080abSJoseph Pusztay 166d0c080abSJoseph Pusztay static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user) 167d0c080abSJoseph Pusztay { 168d0c080abSJoseph Pusztay PetscInt *cellid; 169d0c080abSJoseph Pusztay PetscInt dim, cStart, cEnd, c, Np = user->N, p; 170d0c080abSJoseph Pusztay PetscBool view = PETSC_FALSE; 171d0c080abSJoseph Pusztay 172d0c080abSJoseph Pusztay PetscFunctionBeginUser; 173*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 174*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(PetscObjectComm((PetscObject) dm), sw)); 175*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*sw, DMSWARM)); 176*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetDimension(*sw, dim)); 177d0c080abSJoseph Pusztay /* h = 2L/n and N = n^d */ 178d0c080abSJoseph Pusztay if (user->h < 0.) user->h = 2.*user->L / PetscPowReal(user->N, 1./dim); 179d0c080abSJoseph Pusztay /* From Section 4 in [1], \epsilon = 0.64 h^.98 */ 180d0c080abSJoseph Pusztay if (user->epsilon < 0.) user->epsilon = 0.64*pow(user->h, 1.98); 181*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL, NULL, "-param_view", &view, NULL)); 182d0c080abSJoseph Pusztay if (view) { 183*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "N: %D L: %g h: %g eps: %g\n", user->N, user->L, user->h, user->epsilon)); 184d0c080abSJoseph Pusztay } 185*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetType(*sw, DMSWARM_PIC)); 186*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetCellDM(*sw, dm)); 187*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRegisterPetscDatatypeField(*sw, "velocity", dim, PETSC_REAL)); 188*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR)); 189*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmFinalizeFieldRegister(*sw)); 190*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 191*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0)); 192*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*sw)); 193*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid)); 194d0c080abSJoseph Pusztay for (c = cStart; c < cEnd; ++c) { 195d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 196d0c080abSJoseph Pusztay const PetscInt n = c*Np + p; 197d0c080abSJoseph Pusztay cellid[n] = c; 198d0c080abSJoseph Pusztay } 199d0c080abSJoseph Pusztay } 200*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid)); 201*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) *sw, "Particles")); 202*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*sw, NULL, "-sw_view")); 203d0c080abSJoseph Pusztay PetscFunctionReturn(0); 204d0c080abSJoseph Pusztay } 205d0c080abSJoseph Pusztay 206d0c080abSJoseph Pusztay /* Internal dmplex function, same as found in dmpleximpl.h */ 207d0c080abSJoseph Pusztay static void DMPlex_WaxpyD_Internal(PetscInt dim, PetscReal a, const PetscReal *x, const PetscReal *y, PetscReal *w) 208d0c080abSJoseph Pusztay { 209d0c080abSJoseph Pusztay PetscInt d; 210d0c080abSJoseph Pusztay 211d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) w[d] = a*x[d] + y[d]; 212d0c080abSJoseph Pusztay } 213d0c080abSJoseph Pusztay 214d0c080abSJoseph Pusztay /* Internal dmplex function, same as found in dmpleximpl.h */ 215d0c080abSJoseph Pusztay static PetscReal DMPlex_DotD_Internal(PetscInt dim, const PetscScalar *x, const PetscReal *y) 216d0c080abSJoseph Pusztay { 217d0c080abSJoseph Pusztay PetscReal sum = 0.0; 218d0c080abSJoseph Pusztay PetscInt d; 219d0c080abSJoseph Pusztay 220d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) sum += PetscRealPart(x[d])*y[d]; 221d0c080abSJoseph Pusztay return sum; 222d0c080abSJoseph Pusztay } 223d0c080abSJoseph Pusztay 224d0c080abSJoseph Pusztay /* Internal dmplex function, same as found in dmpleximpl.h */ 225d0c080abSJoseph Pusztay static void DMPlex_MultAdd2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[]) 226d0c080abSJoseph Pusztay { 227d0c080abSJoseph Pusztay PetscScalar z[2]; 228d0c080abSJoseph Pusztay z[0] = x[0]; z[1] = x[ldx]; 229d0c080abSJoseph Pusztay y[0] += A[0]*z[0] + A[1]*z[1]; 230d0c080abSJoseph Pusztay y[ldx] += A[2]*z[0] + A[3]*z[1]; 231d0c080abSJoseph Pusztay (void)PetscLogFlops(6.0); 232d0c080abSJoseph Pusztay } 233d0c080abSJoseph Pusztay 234d0c080abSJoseph Pusztay /* Internal dmplex function, same as found in dmpleximpl.h to avoid private includes. */ 235d0c080abSJoseph Pusztay static void DMPlex_MultAdd3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[]) 236d0c080abSJoseph Pusztay { 237d0c080abSJoseph Pusztay PetscScalar z[3]; 238d0c080abSJoseph Pusztay z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2]; 239d0c080abSJoseph Pusztay y[0] += A[0]*z[0] + A[1]*z[1] + A[2]*z[2]; 240d0c080abSJoseph Pusztay y[ldx] += A[3]*z[0] + A[4]*z[1] + A[5]*z[2]; 241d0c080abSJoseph Pusztay y[ldx*2] += A[6]*z[0] + A[7]*z[1] + A[8]*z[2]; 242d0c080abSJoseph Pusztay (void)PetscLogFlops(15.0); 243d0c080abSJoseph Pusztay } 244d0c080abSJoseph Pusztay 245d0c080abSJoseph Pusztay /* 246d0c080abSJoseph Pusztay Gaussian - The Gaussian function G(x) 247d0c080abSJoseph Pusztay 248d0c080abSJoseph Pusztay Input Parameters: 249d0c080abSJoseph Pusztay + dim - The number of dimensions, or size of x 250d0c080abSJoseph Pusztay . mu - The mean, or center 251d0c080abSJoseph Pusztay . sigma - The standard deviation, or width 252d0c080abSJoseph Pusztay - x - The evaluation point of the function 253d0c080abSJoseph Pusztay 254d0c080abSJoseph Pusztay Output Parameter: 255d0c080abSJoseph Pusztay . ret - The value G(x) 256d0c080abSJoseph Pusztay */ 257d0c080abSJoseph Pusztay static PetscReal Gaussian(PetscInt dim, const PetscReal mu[], PetscReal sigma, const PetscReal x[]) 258d0c080abSJoseph Pusztay { 259d0c080abSJoseph Pusztay PetscReal arg = 0.0; 260d0c080abSJoseph Pusztay PetscInt d; 261d0c080abSJoseph Pusztay 262d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) arg += PetscSqr(x[d] - mu[d]); 263d0c080abSJoseph Pusztay return PetscPowReal(2.0*PETSC_PI*sigma, -dim/2.0) * PetscExpReal(-arg/(2.0*sigma)); 264d0c080abSJoseph Pusztay } 265d0c080abSJoseph Pusztay 266d0c080abSJoseph Pusztay /* 267d0c080abSJoseph Pusztay ComputeGradS - Compute grad_v dS_eps/df 268d0c080abSJoseph Pusztay 269d0c080abSJoseph Pusztay Input Parameters: 270d0c080abSJoseph Pusztay + dim - The dimension 271d0c080abSJoseph Pusztay . Np - The number of particles 272d0c080abSJoseph Pusztay . vp - The velocity v_p of the particle at which we evaluate 273d0c080abSJoseph Pusztay . velocity - The velocity field for all particles 274d0c080abSJoseph Pusztay . epsilon - The regularization strength 275d0c080abSJoseph Pusztay 276d0c080abSJoseph Pusztay Output Parameter: 277d0c080abSJoseph Pusztay . integral - The output grad_v dS_eps/df (v_p) 278d0c080abSJoseph Pusztay 279d0c080abSJoseph Pusztay Note: 280d0c080abSJoseph Pusztay This comes from (3.6) in [1], and we are computing 281d0c080abSJoseph Pusztay $ \nabla_v S_p = \grad \psi_\epsilon(v_p - v) log \sum_q \psi_\epsilon(v - v_q) 282d0c080abSJoseph Pusztay which is discretized by using a one-point quadrature in each box l at its center v^c_l 283d0c080abSJoseph 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) 284d0c080abSJoseph Pusztay where h^d is the volume of each box. 285d0c080abSJoseph Pusztay */ 286d0c080abSJoseph Pusztay static PetscErrorCode ComputeGradS(PetscInt dim, PetscInt Np, const PetscReal vp[], const PetscReal velocity[], PetscReal integral[], AppCtx *ctx) 287d0c080abSJoseph Pusztay { 288d0c080abSJoseph Pusztay PetscReal vc_l[3], L = ctx->L, h = ctx->h, epsilon = ctx->epsilon, init = 0.5*h - L; 289d0c080abSJoseph Pusztay PetscInt nx = roundf(2.*L / h); 290d0c080abSJoseph Pusztay PetscInt ny = dim > 1 ? nx : 1; 291d0c080abSJoseph Pusztay PetscInt nz = dim > 2 ? nx : 1; 292d0c080abSJoseph Pusztay PetscInt i, j, k, d, q, dbg = 0; 293d0c080abSJoseph Pusztay 294d0c080abSJoseph Pusztay PetscFunctionBeginHot; 295d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) integral[d] = 0.0; 296d0c080abSJoseph Pusztay for (k = 0, vc_l[2] = init; k < nz; ++k, vc_l[2] += h) { 297d0c080abSJoseph Pusztay for (j = 0, vc_l[1] = init; j < ny; ++j, vc_l[1] += h) { 298d0c080abSJoseph Pusztay for (i = 0, vc_l[0] = init; i < nx; ++i, vc_l[0] += h) { 299d0c080abSJoseph Pusztay PetscReal sum = 0.0; 300d0c080abSJoseph Pusztay 301*5f80ce2aSJacob Faibussowitsch if (dbg) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "(%D %D) vc_l: %g %g\n", i, j, vc_l[0], vc_l[1])); 302d0c080abSJoseph Pusztay /* \log \sum_k \psi(v - v_k) */ 303d0c080abSJoseph Pusztay for (q = 0; q < Np; ++q) sum += Gaussian(dim, &velocity[q*dim], epsilon, vc_l); 304d0c080abSJoseph Pusztay sum = PetscLogReal(sum); 305d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) integral[d] += (-1./(epsilon))*PetscAbsReal(vp[d] - vc_l[d])*(Gaussian(dim, vp, epsilon, vc_l)) * sum; 306d0c080abSJoseph Pusztay } 307d0c080abSJoseph Pusztay } 308d0c080abSJoseph Pusztay } 309d0c080abSJoseph Pusztay PetscFunctionReturn(0); 310d0c080abSJoseph Pusztay } 311d0c080abSJoseph Pusztay 312d0c080abSJoseph Pusztay /* Q = 1/|xi| (I - xi xi^T / |xi|^2), xi = vp - vq */ 313d0c080abSJoseph Pusztay static PetscErrorCode QCompute(PetscInt dim, const PetscReal vp[], const PetscReal vq[], PetscReal Q[]) 314d0c080abSJoseph Pusztay { 315d0c080abSJoseph Pusztay PetscReal xi[3], xi2, xi3, mag; 316d0c080abSJoseph Pusztay PetscInt d, e; 317d0c080abSJoseph Pusztay 318d0c080abSJoseph Pusztay PetscFunctionBeginHot; 319d0c080abSJoseph Pusztay DMPlex_WaxpyD_Internal(dim, -1.0, vq, vp, xi); 320d0c080abSJoseph Pusztay xi2 = DMPlex_DotD_Internal(dim, xi, xi); 321d0c080abSJoseph Pusztay mag = PetscSqrtReal(xi2); 322d0c080abSJoseph Pusztay xi3 = xi2 * mag; 323d0c080abSJoseph Pusztay for (d = 0; d < dim; ++d) { 324d0c080abSJoseph Pusztay for (e = 0; e < dim; ++e) { 325d0c080abSJoseph Pusztay Q[d*dim+e] = -xi[d]*xi[e] / xi3; 326d0c080abSJoseph Pusztay } 327d0c080abSJoseph Pusztay Q[d*dim+d] += 1. / mag; 328d0c080abSJoseph Pusztay } 329d0c080abSJoseph Pusztay PetscFunctionReturn(0); 330d0c080abSJoseph Pusztay } 331d0c080abSJoseph Pusztay 332d0c080abSJoseph Pusztay static PetscErrorCode RHSFunctionParticles(TS ts, PetscReal t, Vec U, Vec R, void *ctx) 333d0c080abSJoseph Pusztay { 334d0c080abSJoseph Pusztay AppCtx *user = (AppCtx *) ctx; 335d0c080abSJoseph Pusztay PetscInt dbg = 0; 336d0c080abSJoseph Pusztay DM sw; /* Particles */ 337d0c080abSJoseph Pusztay Vec sol; /* Solution vector at current time */ 338d0c080abSJoseph Pusztay const PetscScalar *u; /* input solution vector */ 339d0c080abSJoseph Pusztay PetscScalar *r; 340d0c080abSJoseph Pusztay PetscReal *velocity; 341d0c080abSJoseph Pusztay PetscInt dim, Np, p, q; 342d0c080abSJoseph Pusztay 343d0c080abSJoseph Pusztay PetscFunctionBeginUser; 344*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(R)); 345*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts, &sw)); 346*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(sw, &dim)); 347*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(U, &Np)); 348*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSolution(ts, &sol)); 349*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(sol, &velocity)); 350*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(R, &r)); 351*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(U, &u)); 352d0c080abSJoseph Pusztay Np /= dim; 353*5f80ce2aSJacob Faibussowitsch if (dbg) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Part ppr x y\n")); 354d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 355d0c080abSJoseph Pusztay PetscReal gradS_p[3] = {0., 0., 0.}; 356d0c080abSJoseph Pusztay 357*5f80ce2aSJacob Faibussowitsch CHKERRQ(ComputeGradS(dim, Np, &velocity[p*dim], velocity, gradS_p, user)); 358d0c080abSJoseph Pusztay for (q = 0; q < Np; ++q) { 359d0c080abSJoseph Pusztay PetscReal gradS_q[3] = {0., 0., 0.}, GammaS[3] = {0., 0., 0.}, Q[9]; 360d0c080abSJoseph Pusztay 361d0c080abSJoseph Pusztay if (q == p) continue; 362*5f80ce2aSJacob Faibussowitsch CHKERRQ(ComputeGradS(dim, Np, &velocity[q*dim], velocity, gradS_q, user)); 363d0c080abSJoseph Pusztay DMPlex_WaxpyD_Internal(dim, -1.0, gradS_q, gradS_p, GammaS); 364*5f80ce2aSJacob Faibussowitsch CHKERRQ(QCompute(dim, &u[p*dim], &u[q*dim], Q)); 365d0c080abSJoseph Pusztay switch (dim) { 366d0c080abSJoseph Pusztay case 2: DMPlex_MultAdd2DReal_Internal(Q, 1, GammaS, &r[p*dim]);break; 367d0c080abSJoseph Pusztay case 3: DMPlex_MultAdd3DReal_Internal(Q, 1, GammaS, &r[p*dim]);break; 36898921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Do not support dimension %D", dim); 369d0c080abSJoseph Pusztay } 370d0c080abSJoseph Pusztay } 371*5f80ce2aSJacob Faibussowitsch if (dbg) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Final %4D %10.8lf %10.8lf\n", p, r[p*dim+0], r[p*dim+1])); 372d0c080abSJoseph Pusztay } 373*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(U, &u)); 374*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(R, &r)); 375*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(sol, &velocity)); 376*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(R, NULL, "-residual_view")); 377d0c080abSJoseph Pusztay PetscFunctionReturn(0); 378d0c080abSJoseph Pusztay } 379d0c080abSJoseph Pusztay 380d0c080abSJoseph Pusztay /* 381d0c080abSJoseph Pusztay TS Post Step Function. Copy the solution back into the swarm for migration. We may also need to reform 382d0c080abSJoseph Pusztay the solution vector in cases of particle migration, but we forgo that here since there is no velocity space grid 383d0c080abSJoseph Pusztay to migrate between. 384d0c080abSJoseph Pusztay */ 385d0c080abSJoseph Pusztay static PetscErrorCode UpdateSwarm(TS ts) 386d0c080abSJoseph Pusztay { 387d0c080abSJoseph Pusztay PetscInt idx, n; 388d0c080abSJoseph Pusztay const PetscScalar *u; 389d0c080abSJoseph Pusztay PetscScalar *velocity; 390d0c080abSJoseph Pusztay DM sw; 391d0c080abSJoseph Pusztay Vec sol; 392d0c080abSJoseph Pusztay 393d0c080abSJoseph Pusztay PetscFunctionBeginUser; 394*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts, &sw)); 395*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **) &velocity)); 396*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSolution(ts, &sol)); 397*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(sol, &u)); 398*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(sol, &n)); 399d0c080abSJoseph Pusztay for (idx = 0; idx < n; ++idx) velocity[idx] = u[idx]; 400*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(sol, &u)); 401*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **) &velocity)); 402d0c080abSJoseph Pusztay PetscFunctionReturn(0); 403d0c080abSJoseph Pusztay } 404d0c080abSJoseph Pusztay 405d0c080abSJoseph Pusztay static PetscErrorCode InitializeSolve(TS ts, Vec u) 406d0c080abSJoseph Pusztay { 407d0c080abSJoseph Pusztay DM dm; 408d0c080abSJoseph Pusztay AppCtx *user; 409d0c080abSJoseph Pusztay 410d0c080abSJoseph Pusztay PetscFunctionBeginUser; 411*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts, &dm)); 412*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetApplicationContext(dm, &user)); 413*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetInitialCoordinates(dm)); 414*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetInitialConditions(dm, u)); 415d0c080abSJoseph Pusztay PetscFunctionReturn(0); 416d0c080abSJoseph Pusztay } 417d0c080abSJoseph Pusztay 418d0c080abSJoseph Pusztay int main(int argc,char **argv) 419d0c080abSJoseph Pusztay { 420d0c080abSJoseph Pusztay TS ts; /* nonlinear solver */ 421d0c080abSJoseph Pusztay DM dm, sw; /* Velocity space mesh and Particle Swarm */ 422d0c080abSJoseph Pusztay Vec u, v; /* problem vector */ 423d0c080abSJoseph Pusztay MPI_Comm comm; 424d0c080abSJoseph Pusztay AppCtx user; 425d0c080abSJoseph Pusztay PetscErrorCode ierr; 426d0c080abSJoseph Pusztay 427d0c080abSJoseph Pusztay ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 428d0c080abSJoseph Pusztay comm = PETSC_COMM_WORLD; 429*5f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(comm, &user)); 430d0c080abSJoseph Pusztay /* Initialize objects and set initial conditions */ 431*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(comm, &dm, &user)); 432*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateParticles(dm, &sw, &user)); 433*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetApplicationContext(sw, &user)); 434*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmVectorDefineField(sw, "velocity")); 435*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(comm, &ts)); 436*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetDM(ts, sw)); 437*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts, 10.0)); 438*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts, 0.1)); 439*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxSteps(ts, 1)); 440*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 441*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts, NULL, RHSFunctionParticles, &user)); 442*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 443*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetComputeInitialCondition(ts, InitializeSolve)); 444*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmCreateGlobalVectorFromField(sw, "velocity", &v)); 445*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(v, &u)); 446*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmDestroyGlobalVectorFromField(sw, "velocity", &v)); 447*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSComputeInitialCondition(ts, u)); 448*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetPostStep(ts, UpdateSwarm)); 449*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts, u)); 450*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 451*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 452*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&sw)); 453*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 454d0c080abSJoseph Pusztay ierr = PetscFinalize(); 455d0c080abSJoseph Pusztay return ierr; 456d0c080abSJoseph Pusztay } 457d0c080abSJoseph Pusztay 458d0c080abSJoseph Pusztay /*TEST 459d0c080abSJoseph Pusztay build: 460d0c080abSJoseph Pusztay requires: triangle !single !complex 461d0c080abSJoseph Pusztay test: 462d0c080abSJoseph Pusztay suffix: midpoint 46330602db0SMatthew 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 \ 464d0c080abSJoseph Pusztay -ts_type theta -ts_theta_theta 0.5 -ts_dmswarm_monitor_moments -ts_monitor_frequency 1 -snes_fd 465d0c080abSJoseph Pusztay TEST*/ 466