1c4762a1bSJed Brown static char help[] = "Tests L2 projection with DMSwarm using delta function particles.\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown #include <petscfe.h> 5c4762a1bSJed Brown #include <petscdmswarm.h> 6c4762a1bSJed Brown #include <petscds.h> 7c4762a1bSJed Brown #include <petscksp.h> 8c4762a1bSJed Brown #include <petsc/private/petscfeimpl.h> 9c4762a1bSJed Brown typedef struct { 10832e2b15SMatthew G. Knepley PetscReal L[3]; /* Dimensions of the mesh bounding box */ 11c4762a1bSJed Brown PetscInt particlesPerCell; /* The number of partices per cell */ 12c4762a1bSJed Brown PetscReal particleRelDx; /* Relative particle position perturbation compared to average cell diameter h */ 13c4762a1bSJed Brown PetscReal meshRelDx; /* Relative vertex position perturbation compared to average cell diameter h */ 14c4762a1bSJed Brown PetscInt k; /* Mode number for test function */ 15c4762a1bSJed Brown PetscReal momentTol; /* Tolerance for checking moment conservation */ 16832e2b15SMatthew G. Knepley PetscBool useBlockDiagPrec; /* Use the block diagonal of the normal equations as a preconditioner */ 17c4762a1bSJed Brown PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *); 18c4762a1bSJed Brown } AppCtx; 19c4762a1bSJed Brown 20c4762a1bSJed Brown /* const char *const ex2FunctionTypes[] = {"linear","x2_x4","sin","ex2FunctionTypes","EX2_FUNCTION_",0}; */ 21c4762a1bSJed Brown 22c4762a1bSJed Brown static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx) 23c4762a1bSJed Brown { 24c4762a1bSJed Brown AppCtx *ctx = (AppCtx *) a_ctx; 25c4762a1bSJed Brown PetscInt d; 26c4762a1bSJed Brown 27c4762a1bSJed Brown u[0] = 0.0; 28832e2b15SMatthew G. Knepley for (d = 0; d < dim; ++d) u[0] += x[d]/(ctx->L[d]); 29c4762a1bSJed Brown return 0; 30c4762a1bSJed Brown } 31c4762a1bSJed Brown 32c4762a1bSJed Brown static PetscErrorCode x2_x4(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx) 33c4762a1bSJed Brown { 34c4762a1bSJed Brown AppCtx *ctx = (AppCtx *) a_ctx; 35c4762a1bSJed Brown PetscInt d; 36c4762a1bSJed Brown 37c4762a1bSJed Brown u[0] = 1; 38832e2b15SMatthew G. Knepley for (d = 0; d < dim; ++d) u[0] *= PetscSqr(x[d])*PetscSqr(ctx->L[d]) - PetscPowRealInt(x[d], 4); 39c4762a1bSJed Brown return 0; 40c4762a1bSJed Brown } 41c4762a1bSJed Brown 42c4762a1bSJed Brown static PetscErrorCode sinx(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx) 43c4762a1bSJed Brown { 44c4762a1bSJed Brown AppCtx *ctx = (AppCtx *) a_ctx; 45c4762a1bSJed Brown 46832e2b15SMatthew G. Knepley u[0] = PetscSinScalar(2*PETSC_PI*ctx->k*x[0]/(ctx->L[0])); 47c4762a1bSJed Brown return 0; 48c4762a1bSJed Brown } 49c4762a1bSJed Brown 50c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 51c4762a1bSJed Brown { 52c4762a1bSJed Brown char fstring[PETSC_MAX_PATH_LEN] = "linear"; 53c4762a1bSJed Brown PetscBool flag; 54c4762a1bSJed Brown PetscErrorCode ierr; 55c4762a1bSJed Brown 56c4762a1bSJed Brown PetscFunctionBeginUser; 57c4762a1bSJed Brown options->particlesPerCell = 1; 58c4762a1bSJed Brown options->k = 1; 59c4762a1bSJed Brown options->particleRelDx = 1.e-20; 60c4762a1bSJed Brown options->meshRelDx = 1.e-20; 61c4762a1bSJed Brown options->momentTol = 100.*PETSC_MACHINE_EPSILON; 62832e2b15SMatthew G. Knepley options->useBlockDiagPrec = PETSC_FALSE; 63c4762a1bSJed Brown 64c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "L2 Projection Options", "DMPLEX");CHKERRQ(ierr); 65c4762a1bSJed Brown ierr = PetscOptionsInt("-k", "Mode number of test", "ex2.c", options->k, &options->k, NULL);CHKERRQ(ierr); 66c4762a1bSJed Brown ierr = PetscOptionsInt("-particlesPerCell", "Number of particles per cell", "ex2.c", options->particlesPerCell, &options->particlesPerCell, NULL);CHKERRQ(ierr); 67c4762a1bSJed Brown ierr = PetscOptionsReal("-particle_perturbation", "Relative perturbation of particles (0,1)", "ex2.c", options->particleRelDx, &options->particleRelDx, NULL);CHKERRQ(ierr); 68c4762a1bSJed Brown ierr = PetscOptionsReal("-mesh_perturbation", "Relative perturbation of mesh points (0,1)", "ex2.c", options->meshRelDx, &options->meshRelDx, NULL);CHKERRQ(ierr); 69589a23caSBarry Smith ierr = PetscOptionsString("-function", "Name of test function", "ex2.c", fstring, fstring, sizeof(fstring), NULL);CHKERRQ(ierr); 70c4762a1bSJed Brown ierr = PetscStrcmp(fstring, "linear", &flag);CHKERRQ(ierr); 71c4762a1bSJed Brown if (flag) { 72c4762a1bSJed Brown options->func = linear; 73c4762a1bSJed Brown } else { 74c4762a1bSJed Brown ierr = PetscStrcmp(fstring, "sin", &flag);CHKERRQ(ierr); 75c4762a1bSJed Brown if (flag) { 76c4762a1bSJed Brown options->func = sinx; 77c4762a1bSJed Brown } else { 78c4762a1bSJed Brown ierr = PetscStrcmp(fstring, "x2_x4", &flag);CHKERRQ(ierr); 79c4762a1bSJed Brown options->func = x2_x4; 80c4762a1bSJed Brown if (!flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown function %s",fstring); 81c4762a1bSJed Brown } 82c4762a1bSJed Brown } 83ec817849SMatthew G. Knepley ierr = PetscOptionsReal("-moment_tol", "Tolerance for moment checks", "ex2.c", options->momentTol, &options->momentTol, NULL);CHKERRQ(ierr); 84832e2b15SMatthew G. Knepley ierr = PetscOptionsBool("-block_diag_prec", "Use the block diagonal of the normal equations to precondition the particle projection", "ex2.c", options->useBlockDiagPrec, &options->useBlockDiagPrec, NULL);CHKERRQ(ierr); 851e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 86c4762a1bSJed Brown 87c4762a1bSJed Brown PetscFunctionReturn(0); 88c4762a1bSJed Brown } 89c4762a1bSJed Brown 90c4762a1bSJed Brown static PetscErrorCode PerturbVertices(DM dm, AppCtx *user) 91c4762a1bSJed Brown { 92c4762a1bSJed Brown PetscRandom rnd; 93c4762a1bSJed Brown PetscReal interval = user->meshRelDx; 94c4762a1bSJed Brown Vec coordinates; 95c4762a1bSJed Brown PetscScalar *coords; 96832e2b15SMatthew G. Knepley PetscReal *hh, low[3], high[3]; 97832e2b15SMatthew G. Knepley PetscInt d, cdim, cEnd, N, p, bs; 98c4762a1bSJed Brown PetscErrorCode ierr; 99c4762a1bSJed Brown 100c4762a1bSJed Brown PetscFunctionBeginUser; 101832e2b15SMatthew G. Knepley ierr = DMGetBoundingBox(dm, low, high);CHKERRQ(ierr); 102832e2b15SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, NULL, &cEnd);CHKERRQ(ierr); 103c4762a1bSJed Brown ierr = PetscRandomCreate(PetscObjectComm((PetscObject) dm), &rnd);CHKERRQ(ierr); 104c4762a1bSJed Brown ierr = PetscRandomSetInterval(rnd, -interval, interval);CHKERRQ(ierr); 105c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rnd);CHKERRQ(ierr); 106c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 107c4762a1bSJed Brown ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr); 108832e2b15SMatthew G. Knepley ierr = PetscCalloc1(cdim,&hh);CHKERRQ(ierr); 109832e2b15SMatthew G. Knepley for (d = 0; d < cdim; ++d) hh[d] = (user->L[d])/PetscPowReal(cEnd, 1./cdim); 110c4762a1bSJed Brown ierr = VecGetLocalSize(coordinates, &N);CHKERRQ(ierr); 111c4762a1bSJed Brown ierr = VecGetBlockSize(coordinates, &bs);CHKERRQ(ierr); 112c4762a1bSJed Brown if (bs != cdim) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Coordinate vector has wrong block size %D != %D", bs, cdim); 113c4762a1bSJed Brown ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 114c4762a1bSJed Brown for (p = 0; p < N; p += cdim) { 115c4762a1bSJed Brown PetscScalar *coord = &coords[p], value; 116c4762a1bSJed Brown 117c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 118c4762a1bSJed Brown ierr = PetscRandomGetValue(rnd, &value);CHKERRQ(ierr); 119832e2b15SMatthew G. Knepley coord[d] = PetscMax(low[d], PetscMin(high[d], PetscRealPart(coord[d] + value*hh[d]))); 120c4762a1bSJed Brown } 121c4762a1bSJed Brown } 122c4762a1bSJed Brown ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 123c4762a1bSJed Brown ierr = PetscRandomDestroy(&rnd);CHKERRQ(ierr); 124c4762a1bSJed Brown ierr = PetscFree(hh);CHKERRQ(ierr); 125c4762a1bSJed Brown PetscFunctionReturn(0); 126c4762a1bSJed Brown } 127c4762a1bSJed Brown 128c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user) 129c4762a1bSJed Brown { 130832e2b15SMatthew G. Knepley PetscReal low[3], high[3]; 131832e2b15SMatthew G. Knepley PetscInt cdim, d; 132c4762a1bSJed Brown PetscErrorCode ierr; 133c4762a1bSJed Brown 134c4762a1bSJed Brown PetscFunctionBeginUser; 13530602db0SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 13630602db0SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 13730602db0SMatthew G. Knepley ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 13830602db0SMatthew G. Knepley ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 13930602db0SMatthew G. Knepley 140832e2b15SMatthew G. Knepley ierr = DMGetCoordinateDim(*dm, &cdim);CHKERRQ(ierr); 141832e2b15SMatthew G. Knepley ierr = DMGetBoundingBox(*dm, low, high);CHKERRQ(ierr); 142832e2b15SMatthew G. Knepley for (d = 0; d < cdim; ++d) user->L[d] = high[d] - low[d]; 143c4762a1bSJed Brown ierr = PerturbVertices(*dm, user);CHKERRQ(ierr); 144c4762a1bSJed Brown PetscFunctionReturn(0); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 147c4762a1bSJed Brown static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux, 148c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 149c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 150c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 151c4762a1bSJed Brown { 152c4762a1bSJed Brown g0[0] = 1.0; 153c4762a1bSJed Brown } 154c4762a1bSJed Brown 155c4762a1bSJed Brown static PetscErrorCode CreateFEM(DM dm, AppCtx *user) 156c4762a1bSJed Brown { 157c4762a1bSJed Brown PetscFE fe; 158c4762a1bSJed Brown PetscDS ds; 159832e2b15SMatthew G. Knepley DMPolytopeType ct; 160832e2b15SMatthew G. Knepley PetscBool simplex; 161832e2b15SMatthew G. Knepley PetscInt dim, cStart; 162c4762a1bSJed Brown PetscErrorCode ierr; 163c4762a1bSJed Brown 164c4762a1bSJed Brown PetscFunctionBeginUser; 165c4762a1bSJed Brown ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 166832e2b15SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr); 167832e2b15SMatthew G. Knepley ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr); 168832e2b15SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 169832e2b15SMatthew G. Knepley ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, simplex, NULL, -1, &fe);CHKERRQ(ierr); 170c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fe, "fe");CHKERRQ(ierr); 171c4762a1bSJed Brown ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 172c4762a1bSJed Brown ierr = DMCreateDS(dm);CHKERRQ(ierr); 173c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 174c4762a1bSJed Brown /* Setup to form mass matrix */ 175c4762a1bSJed Brown ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 176c4762a1bSJed Brown ierr = PetscDSSetJacobian(ds, 0, 0, identity, NULL, NULL, NULL);CHKERRQ(ierr); 177c4762a1bSJed Brown PetscFunctionReturn(0); 178c4762a1bSJed Brown } 179c4762a1bSJed Brown 180c4762a1bSJed Brown static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user) 181c4762a1bSJed Brown { 182c4762a1bSJed Brown PetscRandom rnd, rndp; 183c4762a1bSJed Brown PetscReal interval = user->particleRelDx; 184832e2b15SMatthew G. Knepley DMPolytopeType ct; 185832e2b15SMatthew G. Knepley PetscBool simplex; 186c4762a1bSJed Brown PetscScalar value, *vals; 187c4762a1bSJed Brown PetscReal *centroid, *coords, *xi0, *v0, *J, *invJ, detJ; 188c4762a1bSJed Brown PetscInt *cellid; 189832e2b15SMatthew G. Knepley PetscInt Ncell, Np = user->particlesPerCell, p, cStart, c, dim, d; 190c4762a1bSJed Brown PetscErrorCode ierr; 191c4762a1bSJed Brown 192c4762a1bSJed Brown PetscFunctionBeginUser; 193c4762a1bSJed Brown ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 194832e2b15SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr); 195832e2b15SMatthew G. Knepley ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr); 196832e2b15SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 197c4762a1bSJed Brown ierr = DMCreate(PetscObjectComm((PetscObject) dm), sw);CHKERRQ(ierr); 198c4762a1bSJed Brown ierr = DMSetType(*sw, DMSWARM);CHKERRQ(ierr); 199c4762a1bSJed Brown ierr = DMSetDimension(*sw, dim);CHKERRQ(ierr); 200c4762a1bSJed Brown 201c4762a1bSJed Brown ierr = PetscRandomCreate(PetscObjectComm((PetscObject) dm), &rnd);CHKERRQ(ierr); 202c4762a1bSJed Brown ierr = PetscRandomSetInterval(rnd, -1.0, 1.0);CHKERRQ(ierr); 203c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rnd);CHKERRQ(ierr); 204c4762a1bSJed Brown ierr = PetscRandomCreate(PetscObjectComm((PetscObject) dm), &rndp);CHKERRQ(ierr); 205c4762a1bSJed Brown ierr = PetscRandomSetInterval(rndp, -interval, interval);CHKERRQ(ierr); 206c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rndp);CHKERRQ(ierr); 207c4762a1bSJed Brown 208c4762a1bSJed Brown ierr = DMSwarmSetType(*sw, DMSWARM_PIC);CHKERRQ(ierr); 209c4762a1bSJed Brown ierr = DMSwarmSetCellDM(*sw, dm);CHKERRQ(ierr); 210c4762a1bSJed Brown ierr = DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR);CHKERRQ(ierr); 211c4762a1bSJed Brown ierr = DMSwarmFinalizeFieldRegister(*sw);CHKERRQ(ierr); 212c4762a1bSJed Brown ierr = DMPlexGetHeightStratum(dm, 0, NULL, &Ncell);CHKERRQ(ierr); 213c4762a1bSJed Brown ierr = DMSwarmSetLocalSizes(*sw, Ncell * Np, 0);CHKERRQ(ierr); 214c4762a1bSJed Brown ierr = DMSetFromOptions(*sw);CHKERRQ(ierr); 215c4762a1bSJed Brown ierr = DMSwarmGetField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 216c4762a1bSJed Brown ierr = DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr); 217c4762a1bSJed Brown ierr = DMSwarmGetField(*sw, "w_q", NULL, NULL, (void **) &vals);CHKERRQ(ierr); 218c4762a1bSJed Brown 219c4762a1bSJed Brown ierr = PetscMalloc5(dim, ¢roid, dim, &xi0, dim, &v0, dim*dim, &J, dim*dim, &invJ);CHKERRQ(ierr); 220c4762a1bSJed Brown for (c = 0; c < Ncell; ++c) { 221c4762a1bSJed Brown if (Np == 1) { 222c4762a1bSJed Brown ierr = DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL);CHKERRQ(ierr); 223c4762a1bSJed Brown cellid[c] = c; 224c4762a1bSJed Brown for (d = 0; d < dim; ++d) coords[c*dim+d] = centroid[d]; 225c4762a1bSJed Brown } else { 226c4762a1bSJed Brown ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); /* affine */ 227c4762a1bSJed Brown for (d = 0; d < dim; ++d) xi0[d] = -1.0; 228c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 229c4762a1bSJed Brown const PetscInt n = c*Np + p; 230c4762a1bSJed Brown PetscReal sum = 0.0, refcoords[3]; 231c4762a1bSJed Brown 232c4762a1bSJed Brown cellid[n] = c; 233c4762a1bSJed Brown for (d = 0; d < dim; ++d) {ierr = PetscRandomGetValue(rnd, &value);CHKERRQ(ierr); refcoords[d] = PetscRealPart(value); sum += refcoords[d];} 234832e2b15SMatthew G. Knepley if (simplex && sum > 0.0) for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim)*sum; 235c4762a1bSJed Brown CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n*dim]); 236c4762a1bSJed Brown } 237c4762a1bSJed Brown } 238c4762a1bSJed Brown } 239c4762a1bSJed Brown ierr = PetscFree5(centroid, xi0, v0, J, invJ);CHKERRQ(ierr); 240c4762a1bSJed Brown for (c = 0; c < Ncell; ++c) { 241c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 242c4762a1bSJed Brown const PetscInt n = c*Np + p; 243c4762a1bSJed Brown 244c4762a1bSJed Brown for (d = 0; d < dim; ++d) {ierr = PetscRandomGetValue(rndp, &value);CHKERRQ(ierr); coords[n*dim+d] += PetscRealPart(value);} 245c4762a1bSJed Brown user->func(dim, 0.0, &coords[n*dim], 1, &vals[c], user); 246c4762a1bSJed Brown } 247c4762a1bSJed Brown } 248c4762a1bSJed Brown ierr = DMSwarmRestoreField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 249c4762a1bSJed Brown ierr = DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr); 250c4762a1bSJed Brown ierr = DMSwarmRestoreField(*sw, "w_q", NULL, NULL, (void **) &vals);CHKERRQ(ierr); 251c4762a1bSJed Brown ierr = PetscRandomDestroy(&rnd);CHKERRQ(ierr); 252c4762a1bSJed Brown ierr = PetscRandomDestroy(&rndp);CHKERRQ(ierr); 253c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) *sw, "Particles");CHKERRQ(ierr); 254c4762a1bSJed Brown ierr = DMViewFromOptions(*sw, NULL, "-sw_view");CHKERRQ(ierr); 255c4762a1bSJed Brown PetscFunctionReturn(0); 256c4762a1bSJed Brown } 257c4762a1bSJed Brown 258c4762a1bSJed Brown static PetscErrorCode computeParticleMoments(DM sw, PetscReal moments[3], AppCtx *user) 259c4762a1bSJed Brown { 260c4762a1bSJed Brown DM dm; 261c4762a1bSJed Brown const PetscReal *coords; 262c4762a1bSJed Brown const PetscScalar *w; 263c4762a1bSJed Brown PetscReal mom[3] = {0.0, 0.0, 0.0}; 264c4762a1bSJed Brown PetscInt cell, cStart, cEnd, dim; 265c4762a1bSJed Brown PetscErrorCode ierr; 266c4762a1bSJed Brown 267c4762a1bSJed Brown PetscFunctionBeginUser; 268c4762a1bSJed Brown ierr = DMGetDimension(sw, &dim);CHKERRQ(ierr); 269c4762a1bSJed Brown ierr = DMSwarmGetCellDM(sw, &dm);CHKERRQ(ierr); 270c4762a1bSJed Brown ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 271c4762a1bSJed Brown ierr = DMSwarmSortGetAccess(sw);CHKERRQ(ierr); 272c4762a1bSJed Brown ierr = DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 273c4762a1bSJed Brown ierr = DMSwarmGetField(sw, "w_q", NULL, NULL, (void **) &w);CHKERRQ(ierr); 274c4762a1bSJed Brown for (cell = cStart; cell < cEnd; ++cell) { 275c4762a1bSJed Brown PetscInt *pidx; 276c4762a1bSJed Brown PetscInt Np, p, d; 277c4762a1bSJed Brown 278c4762a1bSJed Brown ierr = DMSwarmSortGetPointsPerCell(sw, cell, &Np, &pidx);CHKERRQ(ierr); 279c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 280c4762a1bSJed Brown const PetscInt idx = pidx[p]; 281c4762a1bSJed Brown const PetscReal *c = &coords[idx*dim]; 282c4762a1bSJed Brown 283c4762a1bSJed Brown mom[0] += PetscRealPart(w[idx]); 284c4762a1bSJed Brown mom[1] += PetscRealPart(w[idx]) * c[0]; 285c4762a1bSJed Brown for (d = 0; d < dim; ++d) mom[2] += PetscRealPart(w[idx]) * c[d]*c[d]; 286c4762a1bSJed Brown } 287c4762a1bSJed Brown ierr = PetscFree(pidx);CHKERRQ(ierr); 288c4762a1bSJed Brown } 289c4762a1bSJed Brown ierr = DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 290c4762a1bSJed Brown ierr = DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **) &w);CHKERRQ(ierr); 291c4762a1bSJed Brown ierr = DMSwarmSortRestoreAccess(sw);CHKERRQ(ierr); 292ffc4695bSBarry Smith ierr = MPI_Allreduce(mom, moments, 3, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject) sw));CHKERRMPI(ierr); 293c4762a1bSJed Brown PetscFunctionReturn(0); 294c4762a1bSJed Brown } 295c4762a1bSJed Brown 296c4762a1bSJed Brown static void f0_1(PetscInt dim, PetscInt Nf, PetscInt NfAux, 297c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 298c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 299c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 300c4762a1bSJed Brown { 301c4762a1bSJed Brown f0[0] = u[0]; 302c4762a1bSJed Brown } 303c4762a1bSJed Brown 304c4762a1bSJed Brown static void f0_x(PetscInt dim, PetscInt Nf, PetscInt NfAux, 305c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 306c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 307c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 308c4762a1bSJed Brown { 309c4762a1bSJed Brown f0[0] = x[0]*u[0]; 310c4762a1bSJed Brown } 311c4762a1bSJed Brown 312c4762a1bSJed Brown static void f0_r2(PetscInt dim, PetscInt Nf, PetscInt NfAux, 313c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 314c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 315c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 316c4762a1bSJed Brown { 317c4762a1bSJed Brown PetscInt d; 318c4762a1bSJed Brown 319c4762a1bSJed Brown f0[0] = 0.0; 320c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[0] += PetscSqr(x[d])*u[0]; 321c4762a1bSJed Brown } 322c4762a1bSJed Brown 323c4762a1bSJed Brown static PetscErrorCode computeFEMMoments(DM dm, Vec u, PetscReal moments[3], AppCtx *user) 324c4762a1bSJed Brown { 325c4762a1bSJed Brown PetscDS prob; 326c4762a1bSJed Brown PetscErrorCode ierr; 327c4762a1bSJed Brown PetscScalar mom; 328c4762a1bSJed Brown 329c4762a1bSJed Brown PetscFunctionBeginUser; 330c4762a1bSJed Brown ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 331c4762a1bSJed Brown ierr = PetscDSSetObjective(prob, 0, &f0_1);CHKERRQ(ierr); 332c4762a1bSJed Brown ierr = DMPlexComputeIntegralFEM(dm, u, &mom, user);CHKERRQ(ierr); 333c4762a1bSJed Brown moments[0] = PetscRealPart(mom); 334c4762a1bSJed Brown ierr = PetscDSSetObjective(prob, 0, &f0_x);CHKERRQ(ierr); 335c4762a1bSJed Brown ierr = DMPlexComputeIntegralFEM(dm, u, &mom, user);CHKERRQ(ierr); 336c4762a1bSJed Brown moments[1] = PetscRealPart(mom); 337c4762a1bSJed Brown ierr = PetscDSSetObjective(prob, 0, &f0_r2);CHKERRQ(ierr); 338c4762a1bSJed Brown ierr = DMPlexComputeIntegralFEM(dm, u, &mom, user);CHKERRQ(ierr); 339c4762a1bSJed Brown moments[2] = PetscRealPart(mom); 340c4762a1bSJed Brown PetscFunctionReturn(0); 341c4762a1bSJed Brown } 342c4762a1bSJed Brown 343c4762a1bSJed Brown static PetscErrorCode TestL2ProjectionParticlesToField(DM dm, DM sw, AppCtx *user) 344c4762a1bSJed Brown { 345c4762a1bSJed Brown MPI_Comm comm; 346c4762a1bSJed Brown KSP ksp; 347c4762a1bSJed Brown Mat M; /* FEM mass matrix */ 348c4762a1bSJed Brown Mat M_p; /* Particle mass matrix */ 349c4762a1bSJed Brown Vec f, rhs, fhat; /* Particle field f, \int phi_i f, FEM field */ 350c4762a1bSJed Brown PetscReal pmoments[3]; /* \int f, \int x f, \int r^2 f */ 351c4762a1bSJed Brown PetscReal fmoments[3]; /* \int \hat f, \int x \hat f, \int r^2 \hat f */ 352c4762a1bSJed Brown PetscInt m; 353c4762a1bSJed Brown PetscErrorCode ierr; 354c4762a1bSJed Brown 355c4762a1bSJed Brown PetscFunctionBeginUser; 356c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 357c4762a1bSJed Brown ierr = KSPCreate(comm, &ksp);CHKERRQ(ierr); 358c4762a1bSJed Brown ierr = KSPSetOptionsPrefix(ksp, "ptof_");CHKERRQ(ierr); 359c4762a1bSJed Brown ierr = DMGetGlobalVector(dm, &fhat);CHKERRQ(ierr); 360c4762a1bSJed Brown ierr = DMGetGlobalVector(dm, &rhs);CHKERRQ(ierr); 361c4762a1bSJed Brown 362c4762a1bSJed Brown ierr = DMCreateMassMatrix(sw, dm, &M_p);CHKERRQ(ierr); 363c4762a1bSJed Brown ierr = MatViewFromOptions(M_p, NULL, "-M_p_view");CHKERRQ(ierr); 364c4762a1bSJed Brown 365c4762a1bSJed Brown /* make particle weight vector */ 366c4762a1bSJed Brown ierr = DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f);CHKERRQ(ierr); 367c4762a1bSJed Brown 368c4762a1bSJed Brown /* create matrix RHS vector */ 369c4762a1bSJed Brown ierr = MatMultTranspose(M_p, f, rhs);CHKERRQ(ierr); 370c4762a1bSJed Brown ierr = DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f);CHKERRQ(ierr); 371c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) rhs,"rhs");CHKERRQ(ierr); 372c4762a1bSJed Brown ierr = VecViewFromOptions(rhs, NULL, "-rhs_view");CHKERRQ(ierr); 373c4762a1bSJed Brown 374c4762a1bSJed Brown ierr = DMCreateMatrix(dm, &M);CHKERRQ(ierr); 375c4762a1bSJed Brown ierr = DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user);CHKERRQ(ierr); 376c4762a1bSJed Brown ierr = MatViewFromOptions(M, NULL, "-M_view");CHKERRQ(ierr); 377c4762a1bSJed Brown ierr = KSPSetOperators(ksp, M, M);CHKERRQ(ierr); 378c4762a1bSJed Brown ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 379c4762a1bSJed Brown ierr = KSPSolve(ksp, rhs, fhat);CHKERRQ(ierr); 380c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fhat,"fhat");CHKERRQ(ierr); 381c4762a1bSJed Brown ierr = VecViewFromOptions(fhat, NULL, "-fhat_view");CHKERRQ(ierr); 382c4762a1bSJed Brown 383c4762a1bSJed Brown /* Check moments of field */ 384c4762a1bSJed Brown ierr = computeParticleMoments(sw, pmoments, user);CHKERRQ(ierr); 385c4762a1bSJed Brown ierr = computeFEMMoments(dm, fhat, fmoments, user);CHKERRQ(ierr); 386c4762a1bSJed Brown ierr = PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2]);CHKERRQ(ierr); 387c4762a1bSJed Brown for (m = 0; m < 3; ++m) { 388c4762a1bSJed Brown if (PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]) > user->momentTol) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "Moment %D error too large %g > %g", m, PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]), user->momentTol); 389c4762a1bSJed Brown } 390c4762a1bSJed Brown 391c4762a1bSJed Brown ierr = KSPDestroy(&ksp);CHKERRQ(ierr); 392c4762a1bSJed Brown ierr = MatDestroy(&M);CHKERRQ(ierr); 393c4762a1bSJed Brown ierr = MatDestroy(&M_p);CHKERRQ(ierr); 394c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dm, &fhat);CHKERRQ(ierr); 395c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dm, &rhs);CHKERRQ(ierr); 396c4762a1bSJed Brown 397c4762a1bSJed Brown PetscFunctionReturn(0); 398c4762a1bSJed Brown } 399c4762a1bSJed Brown 400c4762a1bSJed Brown static PetscErrorCode TestL2ProjectionFieldToParticles(DM dm, DM sw, AppCtx *user) 401c4762a1bSJed Brown { 402c4762a1bSJed Brown 403c4762a1bSJed Brown MPI_Comm comm; 404c4762a1bSJed Brown KSP ksp; 405c4762a1bSJed Brown Mat M; /* FEM mass matrix */ 406832e2b15SMatthew G. Knepley Mat M_p, PM_p; /* Particle mass matrix M_p, and the preconditioning matrix, e.g. M_p M^T_p */ 407c4762a1bSJed Brown Vec f, rhs, fhat; /* Particle field f, \int phi_i fhat, FEM field */ 408c4762a1bSJed Brown PetscReal pmoments[3]; /* \int f, \int x f, \int r^2 f */ 409c4762a1bSJed Brown PetscReal fmoments[3]; /* \int \hat f, \int x \hat f, \int r^2 \hat f */ 410c4762a1bSJed Brown PetscInt m; 411c4762a1bSJed Brown PetscErrorCode ierr; 412c4762a1bSJed Brown 413c4762a1bSJed Brown PetscFunctionBeginUser; 414c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 415c4762a1bSJed Brown 416c4762a1bSJed Brown ierr = KSPCreate(comm, &ksp);CHKERRQ(ierr); 417c4762a1bSJed Brown ierr = KSPSetOptionsPrefix(ksp, "ftop_");CHKERRQ(ierr); 418c4762a1bSJed Brown ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 419c4762a1bSJed Brown 420c4762a1bSJed Brown ierr = DMGetGlobalVector(dm, &fhat);CHKERRQ(ierr); 421c4762a1bSJed Brown ierr = DMGetGlobalVector(dm, &rhs);CHKERRQ(ierr); 422c4762a1bSJed Brown 423c4762a1bSJed Brown ierr = DMCreateMassMatrix(sw, dm, &M_p);CHKERRQ(ierr); 424c4762a1bSJed Brown ierr = MatViewFromOptions(M_p, NULL, "-M_p_view");CHKERRQ(ierr); 425c4762a1bSJed Brown 426c4762a1bSJed Brown /* make particle weight vector */ 427c4762a1bSJed Brown ierr = DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f);CHKERRQ(ierr); 428c4762a1bSJed Brown 429c4762a1bSJed Brown /* create matrix RHS vector, in this case the FEM field fhat with the coefficients vector #alpha */ 430c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) rhs,"rhs");CHKERRQ(ierr); 431c4762a1bSJed Brown ierr = VecViewFromOptions(rhs, NULL, "-rhs_view");CHKERRQ(ierr); 432c4762a1bSJed Brown ierr = DMCreateMatrix(dm, &M);CHKERRQ(ierr); 433c4762a1bSJed Brown ierr = DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user);CHKERRQ(ierr); 434c4762a1bSJed Brown ierr = MatViewFromOptions(M, NULL, "-M_view");CHKERRQ(ierr); 435c4762a1bSJed Brown ierr = MatMultTranspose(M, fhat, rhs);CHKERRQ(ierr); 436832e2b15SMatthew G. Knepley if (user->useBlockDiagPrec) {ierr = DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p);CHKERRQ(ierr);} 437832e2b15SMatthew G. Knepley else {ierr = PetscObjectReference((PetscObject) M_p);CHKERRQ(ierr); PM_p = M_p;} 438c4762a1bSJed Brown 439832e2b15SMatthew G. Knepley ierr = KSPSetOperators(ksp, M_p, PM_p);CHKERRQ(ierr); 440c4762a1bSJed Brown ierr = KSPSolveTranspose(ksp, rhs, f);CHKERRQ(ierr); 441c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fhat,"fhat");CHKERRQ(ierr); 442c4762a1bSJed Brown ierr = VecViewFromOptions(fhat, NULL, "-fhat_view");CHKERRQ(ierr); 443c4762a1bSJed Brown 444c4762a1bSJed Brown ierr = DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f);CHKERRQ(ierr); 445c4762a1bSJed Brown 446c4762a1bSJed Brown /* Check moments */ 447c4762a1bSJed Brown ierr = computeParticleMoments(sw, pmoments, user);CHKERRQ(ierr); 448c4762a1bSJed Brown ierr = computeFEMMoments(dm, fhat, fmoments, user);CHKERRQ(ierr); 449c4762a1bSJed Brown ierr = PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2]);CHKERRQ(ierr); 450c4762a1bSJed Brown for (m = 0; m < 3; ++m) { 451c4762a1bSJed Brown if (PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]) > user->momentTol) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "Moment %D error too large %g > %g", m, PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]), user->momentTol); 452c4762a1bSJed Brown } 453c4762a1bSJed Brown ierr = KSPDestroy(&ksp);CHKERRQ(ierr); 454c4762a1bSJed Brown ierr = MatDestroy(&M);CHKERRQ(ierr); 455c4762a1bSJed Brown ierr = MatDestroy(&M_p);CHKERRQ(ierr); 456832e2b15SMatthew G. Knepley ierr = MatDestroy(&PM_p);CHKERRQ(ierr); 457c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dm, &fhat);CHKERRQ(ierr); 458c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dm, &rhs);CHKERRQ(ierr); 459c4762a1bSJed Brown PetscFunctionReturn(0); 460c4762a1bSJed Brown } 461c4762a1bSJed Brown 462c4762a1bSJed Brown /* Interpolate the gradient of an FEM (FVM) field. Code repurposed from DMPlexComputeGradientClementInterpolant */ 46370a7d78aSStefano Zampini static PetscErrorCode InterpolateGradient(DM dm, Vec locX, Vec locC) 46470a7d78aSStefano Zampini { 465c4762a1bSJed Brown DM_Plex *mesh = (DM_Plex *) dm->data; 466c4762a1bSJed Brown PetscInt debug = mesh->printFEM; 467c4762a1bSJed Brown DM dmC; 468c4762a1bSJed Brown PetscSection section; 469c4762a1bSJed Brown PetscQuadrature quad = NULL; 470c4762a1bSJed Brown PetscScalar *interpolant, *gradsum; 471c4762a1bSJed Brown PetscFEGeom fegeom; 472c4762a1bSJed Brown PetscReal *coords; 473c4762a1bSJed Brown const PetscReal *quadPoints, *quadWeights; 474c4762a1bSJed Brown PetscInt dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, vStart, vEnd, v, field, fieldOffset; 475c4762a1bSJed Brown PetscErrorCode ierr; 476c4762a1bSJed Brown 477c4762a1bSJed Brown PetscFunctionBegin; 478c4762a1bSJed Brown ierr = VecGetDM(locC, &dmC);CHKERRQ(ierr); 479c4762a1bSJed Brown ierr = VecSet(locC, 0.0);CHKERRQ(ierr); 480c4762a1bSJed Brown ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 481c4762a1bSJed Brown ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr); 482c4762a1bSJed Brown fegeom.dimEmbed = coordDim; 483c4762a1bSJed Brown ierr = DMGetLocalSection(dm, §ion);CHKERRQ(ierr); 484c4762a1bSJed Brown ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 485c4762a1bSJed Brown for (field = 0; field < numFields; ++field) { 486c4762a1bSJed Brown PetscObject obj; 487c4762a1bSJed Brown PetscClassId id; 488c4762a1bSJed Brown PetscInt Nc; 489c4762a1bSJed Brown 490c4762a1bSJed Brown ierr = DMGetField(dm, field, NULL, &obj);CHKERRQ(ierr); 491c4762a1bSJed Brown ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 492c4762a1bSJed Brown if (id == PETSCFE_CLASSID) { 493c4762a1bSJed Brown PetscFE fe = (PetscFE) obj; 494c4762a1bSJed Brown 495c4762a1bSJed Brown ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 496c4762a1bSJed Brown ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 497c4762a1bSJed Brown } else if (id == PETSCFV_CLASSID) { 498c4762a1bSJed Brown PetscFV fv = (PetscFV) obj; 499c4762a1bSJed Brown 500c4762a1bSJed Brown ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr); 501c4762a1bSJed Brown ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 502c4762a1bSJed Brown } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 503c4762a1bSJed Brown numComponents += Nc; 504c4762a1bSJed Brown } 505c4762a1bSJed Brown ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 506c4762a1bSJed Brown if ((qNc != 1) && (qNc != numComponents)) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Quadrature components %D != %D field components", qNc, numComponents); 507c4762a1bSJed Brown ierr = PetscMalloc6(coordDim*numComponents*2,&gradsum,coordDim*numComponents,&interpolant,coordDim*Nq,&coords,Nq,&fegeom.detJ,coordDim*coordDim*Nq,&fegeom.J,coordDim*coordDim*Nq,&fegeom.invJ);CHKERRQ(ierr); 508c4762a1bSJed Brown ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 50954fcfd0cSMatthew G. Knepley ierr = DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 510c4762a1bSJed Brown for (v = vStart; v < vEnd; ++v) { 511c4762a1bSJed Brown PetscInt *star = NULL; 512c4762a1bSJed Brown PetscInt starSize, st, d, fc; 513c4762a1bSJed Brown 514c4762a1bSJed Brown ierr = PetscArrayzero(gradsum, coordDim*numComponents);CHKERRQ(ierr); 515c4762a1bSJed Brown ierr = DMPlexGetTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 516c4762a1bSJed Brown for (st = 0; st < starSize*2; st += 2) { 517c4762a1bSJed Brown const PetscInt cell = star[st]; 518c4762a1bSJed Brown PetscScalar *grad = &gradsum[coordDim*numComponents]; 519c4762a1bSJed Brown PetscScalar *x = NULL; 520c4762a1bSJed Brown 521c4762a1bSJed Brown if ((cell < cStart) || (cell >= cEnd)) continue; 522c4762a1bSJed Brown ierr = DMPlexComputeCellGeometryFEM(dm, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ);CHKERRQ(ierr); 523c4762a1bSJed Brown ierr = DMPlexVecGetClosure(dm, NULL, locX, cell, NULL, &x);CHKERRQ(ierr); 524c4762a1bSJed Brown for (field = 0, fieldOffset = 0; field < numFields; ++field) { 525c4762a1bSJed Brown PetscObject obj; 526c4762a1bSJed Brown PetscClassId id; 527c4762a1bSJed Brown PetscInt Nb, Nc, q, qc = 0; 528c4762a1bSJed Brown 529c4762a1bSJed Brown ierr = PetscArrayzero(grad, coordDim*numComponents);CHKERRQ(ierr); 530c4762a1bSJed Brown ierr = DMGetField(dm, field, NULL, &obj);CHKERRQ(ierr); 531c4762a1bSJed Brown ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 532c4762a1bSJed Brown if (id == PETSCFE_CLASSID) {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);} 533c4762a1bSJed Brown else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;} 534c4762a1bSJed Brown else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 535c4762a1bSJed Brown for (q = 0; q < Nq; ++q) { 536c4762a1bSJed Brown if (fegeom.detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, quadrature points %D", (double)fegeom.detJ[q], cell, q); 537c4762a1bSJed Brown if (ierr) { 538c4762a1bSJed Brown PetscErrorCode ierr2; 539c4762a1bSJed Brown ierr2 = DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x);CHKERRQ(ierr2); 540c4762a1bSJed Brown ierr2 = DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr2); 541c4762a1bSJed Brown ierr2 = PetscFree6(gradsum,interpolant,coords,fegeom.detJ,fegeom.J,fegeom.invJ);CHKERRQ(ierr2); 542c4762a1bSJed Brown CHKERRQ(ierr); 543c4762a1bSJed Brown } 544f9244615SMatthew G. Knepley if (id == PETSCFE_CLASSID) {ierr = PetscFEInterpolateGradient_Static((PetscFE) obj, 1, &x[fieldOffset], &fegeom, q, interpolant);CHKERRQ(ierr);} 545c4762a1bSJed Brown else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 546c4762a1bSJed Brown for (fc = 0; fc < Nc; ++fc) { 547c4762a1bSJed Brown const PetscReal wt = quadWeights[q*qNc+qc+fc]; 548c4762a1bSJed Brown 549c4762a1bSJed Brown for (d = 0; d < coordDim; ++d) grad[fc*coordDim+d] += interpolant[fc*dim+d]*wt*fegeom.detJ[q]; 550c4762a1bSJed Brown } 551c4762a1bSJed Brown } 552c4762a1bSJed Brown fieldOffset += Nb; 553c4762a1bSJed Brown } 554c4762a1bSJed Brown ierr = DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x);CHKERRQ(ierr); 555c4762a1bSJed Brown for (fc = 0; fc < numComponents; ++fc) { 556c4762a1bSJed Brown for (d = 0; d < coordDim; ++d) { 557c4762a1bSJed Brown gradsum[fc*coordDim+d] += grad[fc*coordDim+d]; 558c4762a1bSJed Brown } 559c4762a1bSJed Brown } 560c4762a1bSJed Brown if (debug) { 561c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D gradient: [", cell);CHKERRQ(ierr); 562c4762a1bSJed Brown for (fc = 0; fc < numComponents; ++fc) { 563c4762a1bSJed Brown for (d = 0; d < coordDim; ++d) { 564c4762a1bSJed Brown if (fc || d > 0) {ierr = PetscPrintf(PETSC_COMM_SELF, ", ");CHKERRQ(ierr);} 565c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(grad[fc*coordDim+d]));CHKERRQ(ierr); 566c4762a1bSJed Brown } 567c4762a1bSJed Brown } 568c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "]\n");CHKERRQ(ierr); 569c4762a1bSJed Brown } 570c4762a1bSJed Brown } 571c4762a1bSJed Brown ierr = DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 572c4762a1bSJed Brown ierr = DMPlexVecSetClosure(dmC, NULL, locC, v, gradsum, INSERT_VALUES);CHKERRQ(ierr); 573c4762a1bSJed Brown } 574c4762a1bSJed Brown ierr = PetscFree6(gradsum,interpolant,coords,fegeom.detJ,fegeom.J,fegeom.invJ);CHKERRQ(ierr); 575c4762a1bSJed Brown PetscFunctionReturn(0); 576c4762a1bSJed Brown } 577c4762a1bSJed Brown 578c4762a1bSJed Brown static PetscErrorCode TestFieldGradientProjection(DM dm, DM sw, AppCtx *user) 579c4762a1bSJed Brown { 580c4762a1bSJed Brown 581c4762a1bSJed Brown MPI_Comm comm; 582c4762a1bSJed Brown KSP ksp; 583c4762a1bSJed Brown Mat M; /* FEM mass matrix */ 584c4762a1bSJed Brown Mat M_p; /* Particle mass matrix */ 585c4762a1bSJed Brown Vec f, rhs, fhat, grad; /* Particle field f, \int phi_i f, FEM field */ 586c4762a1bSJed Brown PetscReal pmoments[3]; /* \int f, \int x f, \int r^2 f */ 587c4762a1bSJed Brown PetscReal fmoments[3]; /* \int \hat f, \int x \hat f, \int r^2 \hat f */ 588c4762a1bSJed Brown PetscInt m; 589c4762a1bSJed Brown PetscErrorCode ierr; 590c4762a1bSJed Brown 591c4762a1bSJed Brown PetscFunctionBeginUser; 592c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 593c4762a1bSJed Brown ierr = KSPCreate(comm, &ksp);CHKERRQ(ierr); 594c4762a1bSJed Brown ierr = KSPSetOptionsPrefix(ksp, "ptof_");CHKERRQ(ierr); 595c4762a1bSJed Brown ierr = DMGetGlobalVector(dm, &fhat);CHKERRQ(ierr); 596c4762a1bSJed Brown ierr = DMGetGlobalVector(dm, &rhs);CHKERRQ(ierr); 597c4762a1bSJed Brown ierr = DMGetGlobalVector(dm, &grad);CHKERRQ(ierr); 598c4762a1bSJed Brown 599c4762a1bSJed Brown ierr = DMCreateMassMatrix(sw, dm, &M_p);CHKERRQ(ierr); 600c4762a1bSJed Brown ierr = MatViewFromOptions(M_p, NULL, "-M_p_view");CHKERRQ(ierr); 601c4762a1bSJed Brown 602c4762a1bSJed Brown /* make particle weight vector */ 603c4762a1bSJed Brown ierr = DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f);CHKERRQ(ierr); 604c4762a1bSJed Brown 605c4762a1bSJed Brown /* create matrix RHS vector */ 606c4762a1bSJed Brown ierr = MatMultTranspose(M_p, f, rhs);CHKERRQ(ierr); 607c4762a1bSJed Brown ierr = DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f);CHKERRQ(ierr); 608c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) rhs,"rhs");CHKERRQ(ierr); 609c4762a1bSJed Brown ierr = VecViewFromOptions(rhs, NULL, "-rhs_view");CHKERRQ(ierr); 610c4762a1bSJed Brown 611c4762a1bSJed Brown ierr = DMCreateMatrix(dm, &M);CHKERRQ(ierr); 612c4762a1bSJed Brown ierr = DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user);CHKERRQ(ierr); 613c4762a1bSJed Brown 614c4762a1bSJed Brown ierr = InterpolateGradient(dm, fhat, grad);CHKERRQ(ierr); 615c4762a1bSJed Brown 616c4762a1bSJed Brown ierr = MatViewFromOptions(M, NULL, "-M_view");CHKERRQ(ierr); 617c4762a1bSJed Brown ierr = KSPSetOperators(ksp, M, M);CHKERRQ(ierr); 618c4762a1bSJed Brown ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 619c4762a1bSJed Brown ierr = KSPSolve(ksp, rhs, grad);CHKERRQ(ierr); 620c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fhat,"fhat");CHKERRQ(ierr); 621c4762a1bSJed Brown ierr = VecViewFromOptions(fhat, NULL, "-fhat_view");CHKERRQ(ierr); 622c4762a1bSJed Brown 623c4762a1bSJed Brown /* Check moments of field */ 624c4762a1bSJed Brown ierr = computeParticleMoments(sw, pmoments, user);CHKERRQ(ierr); 625c4762a1bSJed Brown ierr = computeFEMMoments(dm, grad, fmoments, user);CHKERRQ(ierr); 626c4762a1bSJed Brown ierr = PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2]);CHKERRQ(ierr); 627c4762a1bSJed Brown for (m = 0; m < 3; ++m) { 628c4762a1bSJed Brown if (PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]) > user->momentTol) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "Moment %D error too large %g > %g", m, PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]), user->momentTol); 629c4762a1bSJed Brown } 630c4762a1bSJed Brown 631c4762a1bSJed Brown ierr = KSPDestroy(&ksp);CHKERRQ(ierr); 632c4762a1bSJed Brown ierr = MatDestroy(&M);CHKERRQ(ierr); 633c4762a1bSJed Brown ierr = MatDestroy(&M_p);CHKERRQ(ierr); 634c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dm, &fhat);CHKERRQ(ierr); 635c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dm, &rhs);CHKERRQ(ierr); 636c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dm, &grad);CHKERRQ(ierr); 637c4762a1bSJed Brown 638c4762a1bSJed Brown PetscFunctionReturn(0); 639c4762a1bSJed Brown } 640c4762a1bSJed Brown 641*5627991aSBarry Smith int main (int argc, char *argv[]) 642*5627991aSBarry Smith { 643c4762a1bSJed Brown MPI_Comm comm; 644c4762a1bSJed Brown DM dm, sw; 645c4762a1bSJed Brown AppCtx user; 646c4762a1bSJed Brown PetscErrorCode ierr; 647c4762a1bSJed Brown 648c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 649c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 650c4762a1bSJed Brown ierr = ProcessOptions(comm, &user);CHKERRQ(ierr); 651c4762a1bSJed Brown ierr = CreateMesh(comm, &dm, &user);CHKERRQ(ierr); 652c4762a1bSJed Brown ierr = CreateFEM(dm, &user);CHKERRQ(ierr); 653c4762a1bSJed Brown ierr = CreateParticles(dm, &sw, &user);CHKERRQ(ierr); 654c4762a1bSJed Brown ierr = TestL2ProjectionParticlesToField(dm, sw, &user);CHKERRQ(ierr); 655c4762a1bSJed Brown ierr = TestL2ProjectionFieldToParticles(dm, sw, &user);CHKERRQ(ierr); 656c4762a1bSJed Brown ierr = TestFieldGradientProjection(dm, sw, &user);CHKERRQ(ierr); 657c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 658c4762a1bSJed Brown ierr = DMDestroy(&sw);CHKERRQ(ierr); 659c4762a1bSJed Brown ierr = PetscFinalize(); 660c4762a1bSJed Brown return ierr; 661c4762a1bSJed Brown } 662c4762a1bSJed Brown 663c4762a1bSJed Brown /*TEST 664c4762a1bSJed Brown 665832e2b15SMatthew G. Knepley # Swarm does not handle complex or quad 666832e2b15SMatthew G. Knepley build: 667832e2b15SMatthew G. Knepley requires: !complex double 668832e2b15SMatthew G. Knepley 669c4762a1bSJed Brown test: 670c4762a1bSJed Brown suffix: proj_tri_0 671832e2b15SMatthew G. Knepley requires: triangle 672832e2b15SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 673c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 674c4762a1bSJed Brown 675c4762a1bSJed Brown test: 676c4762a1bSJed Brown suffix: proj_tri_2_faces 677832e2b15SMatthew G. Knepley requires: triangle 678832e2b15SMatthew G. Knepley args: -dm_plex_box_faces 2,2 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 679c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 680c4762a1bSJed Brown 681c4762a1bSJed Brown test: 682c4762a1bSJed Brown suffix: proj_quad_0 683832e2b15SMatthew G. Knepley requires: triangle 68430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 685c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 686c4762a1bSJed Brown 687c4762a1bSJed Brown test: 688c4762a1bSJed Brown suffix: proj_quad_2_faces 689832e2b15SMatthew G. Knepley requires: triangle 69030602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 2,2 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 691c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 692c4762a1bSJed Brown 693c4762a1bSJed Brown test: 694c4762a1bSJed Brown suffix: proj_tri_5P 695832e2b15SMatthew G. Knepley requires: triangle 696832e2b15SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -particlesPerCell 5 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 697c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 698c4762a1bSJed Brown 699c4762a1bSJed Brown test: 700c4762a1bSJed Brown suffix: proj_quad_5P 701832e2b15SMatthew G. Knepley requires: triangle 70230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -particlesPerCell 5 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 703c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 704c4762a1bSJed Brown 705c4762a1bSJed Brown test: 706c4762a1bSJed Brown suffix: proj_tri_mdx 707832e2b15SMatthew G. Knepley requires: triangle 708832e2b15SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 709c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 710c4762a1bSJed Brown 711c4762a1bSJed Brown test: 712c4762a1bSJed Brown suffix: proj_tri_mdx_5P 713832e2b15SMatthew G. Knepley requires: triangle 714832e2b15SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -particlesPerCell 5 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 715c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 716c4762a1bSJed Brown 717c4762a1bSJed Brown test: 718c4762a1bSJed Brown suffix: proj_tri_3d 719832e2b15SMatthew G. Knepley requires: ctetgen 72030602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 1,1,1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 721c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 722c4762a1bSJed Brown 723c4762a1bSJed Brown test: 724c4762a1bSJed Brown suffix: proj_tri_3d_2_faces 725832e2b15SMatthew G. Knepley requires: ctetgen 72630602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 727c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 728c4762a1bSJed Brown 729c4762a1bSJed Brown test: 730c4762a1bSJed Brown suffix: proj_tri_3d_5P 731832e2b15SMatthew G. Knepley requires: ctetgen 73230602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 1,1,1 -particlesPerCell 5 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 733c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 734c4762a1bSJed Brown 735c4762a1bSJed Brown test: 736c4762a1bSJed Brown suffix: proj_tri_3d_mdx 737832e2b15SMatthew G. Knepley requires: ctetgen 73830602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 1,1,1 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 739c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 740c4762a1bSJed Brown 741c4762a1bSJed Brown test: 742c4762a1bSJed Brown suffix: proj_tri_3d_mdx_5P 743832e2b15SMatthew G. Knepley requires: ctetgen 74430602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 1,1,1 -particlesPerCell 5 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 745c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 746c4762a1bSJed Brown 747c4762a1bSJed Brown test: 748c4762a1bSJed Brown suffix: proj_tri_3d_mdx_2_faces 749832e2b15SMatthew G. Knepley requires: ctetgen 75030602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 751c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 752c4762a1bSJed Brown 753c4762a1bSJed Brown test: 754c4762a1bSJed Brown suffix: proj_tri_3d_mdx_5P_2_faces 755832e2b15SMatthew G. Knepley requires: ctetgen 75630602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -particlesPerCell 5 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 757c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 758c4762a1bSJed Brown 759832e2b15SMatthew G. Knepley test: 760832e2b15SMatthew G. Knepley suffix: proj_quad_lsqr_scale 761832e2b15SMatthew G. Knepley requires: !complex 76230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 4,4 \ 763832e2b15SMatthew G. Knepley -petscspace_degree 2 -petscfe_default_quadrature_order 3 \ 764832e2b15SMatthew G. Knepley -particlesPerCell 200 \ 765832e2b15SMatthew G. Knepley -ptof_pc_type lu \ 766832e2b15SMatthew G. Knepley -ftop_ksp_rtol 1e-17 -ftop_ksp_type lsqr -ftop_pc_type none 767832e2b15SMatthew G. Knepley 768832e2b15SMatthew G. Knepley test: 769832e2b15SMatthew G. Knepley suffix: proj_quad_lsqr_prec_scale 770832e2b15SMatthew G. Knepley requires: !complex 77130602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 4,4 \ 772832e2b15SMatthew G. Knepley -petscspace_degree 2 -petscfe_default_quadrature_order 3 \ 773832e2b15SMatthew G. Knepley -particlesPerCell 200 \ 774832e2b15SMatthew G. Knepley -ptof_pc_type lu \ 775832e2b15SMatthew G. Knepley -ftop_ksp_rtol 1e-17 -ftop_ksp_type lsqr -ftop_pc_type lu -ftop_pc_factor_shift_type nonzero -block_diag_prec 776832e2b15SMatthew G. Knepley 777c4762a1bSJed Brown TEST*/ 778