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