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