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