1*98285ab4SZach Atkins // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2*98285ab4SZach Atkins // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3*98285ab4SZach Atkins // 4*98285ab4SZach Atkins // SPDX-License-Identifier: BSD-2-Clause 5*98285ab4SZach Atkins // 6*98285ab4SZach Atkins // This file is part of CEED: http://github.com/ceed 7*98285ab4SZach Atkins 8b6972d74SZach Atkins #include "../include/swarmutils.h" 9b6972d74SZach Atkins #include "../qfunctions/swarm/swarmmass.h" 10b6972d74SZach Atkins 11b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 12b6972d74SZach Atkins // Context utilities 13b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 14b6972d74SZach Atkins PetscErrorCode DMSwarmCeedContextCreate(DM dm_swarm, const char *ceed_resource, DMSwarmCeedContext *ctx) { 15b6972d74SZach Atkins DM dm_mesh, dm_coord; 16b6972d74SZach Atkins CeedElemRestriction elem_restr_u_mesh, elem_restr_x_mesh, elem_restr_x_points, elem_restr_u_points, elem_restr_q_data_points; 17b6972d74SZach Atkins CeedBasis basis_u, basis_x; 18b6972d74SZach Atkins CeedVector x_ref_points, q_data_points; 19b6972d74SZach Atkins CeedInt num_comp; 20b6972d74SZach Atkins 21b6972d74SZach Atkins PetscFunctionBeginUser; 22b6972d74SZach Atkins PetscCall(PetscNew(ctx)); 23b6972d74SZach Atkins PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh)); 24b6972d74SZach Atkins PetscCall(DMGetCoordinateDM(dm_mesh, &dm_coord)); 25b6972d74SZach Atkins 26b6972d74SZach Atkins CeedInit(ceed_resource, &(*ctx)->ceed); 27b6972d74SZach Atkins // Background mesh objects 28b6972d74SZach Atkins { 29b6972d74SZach Atkins BPData bp_data = {.q_mode = CEED_GAUSS}; 30b6972d74SZach Atkins 31b6972d74SZach Atkins PetscCall(CreateBasisFromPlex((*ctx)->ceed, dm_mesh, NULL, 0, 0, 0, bp_data, &basis_u)); 32b6972d74SZach Atkins PetscCall(CreateBasisFromPlex((*ctx)->ceed, dm_coord, NULL, 0, 0, 0, bp_data, &basis_x)); 33b6972d74SZach Atkins PetscCall(CreateRestrictionFromPlex((*ctx)->ceed, dm_mesh, 0, NULL, 0, &elem_restr_u_mesh)); 34b6972d74SZach Atkins PetscCall(CreateRestrictionFromPlex((*ctx)->ceed, dm_coord, 0, NULL, 0, &elem_restr_x_mesh)); 35b6972d74SZach Atkins 36b6972d74SZach Atkins // -- Mesh vectors 37b6972d74SZach Atkins CeedElemRestrictionCreateVector(elem_restr_u_mesh, &(*ctx)->u_mesh, NULL); 38b6972d74SZach Atkins CeedElemRestrictionCreateVector(elem_restr_u_mesh, &(*ctx)->v_mesh, NULL); 39b6972d74SZach Atkins } 40b6972d74SZach Atkins // Swarm objects 41b6972d74SZach Atkins { 42b6972d74SZach Atkins PetscInt dim; 43b6972d74SZach Atkins const PetscInt *cell_points; 44b6972d74SZach Atkins IS is_points; 45b6972d74SZach Atkins Vec X_ref; 46b6972d74SZach Atkins CeedInt num_elem; 47b6972d74SZach Atkins 48b6972d74SZach Atkins PetscCall(DMSwarmCreateReferenceCoordinates(dm_swarm, &is_points, &X_ref)); 49b6972d74SZach Atkins PetscCall(DMGetDimension(dm_mesh, &dim)); 50b6972d74SZach Atkins CeedElemRestrictionGetNumElements(elem_restr_u_mesh, &num_elem); 51b6972d74SZach Atkins CeedElemRestrictionGetNumComponents(elem_restr_u_mesh, &num_comp); 52b6972d74SZach Atkins 53b6972d74SZach Atkins PetscCall(ISGetIndices(is_points, &cell_points)); 54b6972d74SZach Atkins PetscInt num_points = cell_points[num_elem + 1] - num_elem - 2; 55b6972d74SZach Atkins CeedInt offsets[num_elem + 1 + num_points]; 56b6972d74SZach Atkins 57b6972d74SZach Atkins for (PetscInt i = 0; i < num_elem + 1; i++) offsets[i] = cell_points[i + 1] - 1; 58b6972d74SZach Atkins for (PetscInt i = num_elem + 1; i < num_points + num_elem + 1; i++) offsets[i] = cell_points[i + 1]; 59b6972d74SZach Atkins PetscCall(ISRestoreIndices(is_points, &cell_points)); 60b6972d74SZach Atkins 61b6972d74SZach Atkins // -- Points restrictions 62b6972d74SZach Atkins CeedElemRestrictionCreateAtPoints((*ctx)->ceed, num_elem, num_points, num_comp, num_points * num_comp, CEED_MEM_HOST, CEED_COPY_VALUES, offsets, 63b6972d74SZach Atkins &elem_restr_u_points); 64b6972d74SZach Atkins CeedElemRestrictionCreateAtPoints((*ctx)->ceed, num_elem, num_points, dim, num_points * dim, CEED_MEM_HOST, CEED_COPY_VALUES, offsets, 65b6972d74SZach Atkins &elem_restr_x_points); 66b6972d74SZach Atkins CeedElemRestrictionCreateAtPoints((*ctx)->ceed, num_elem, num_points, 1, num_points, CEED_MEM_HOST, CEED_COPY_VALUES, offsets, 67b6972d74SZach Atkins &elem_restr_q_data_points); 68b6972d74SZach Atkins 69b6972d74SZach Atkins // -- Points vectors 70b6972d74SZach Atkins CeedElemRestrictionCreateVector(elem_restr_u_points, &(*ctx)->u_points, NULL); 71b6972d74SZach Atkins CeedElemRestrictionCreateVector(elem_restr_q_data_points, &q_data_points, NULL); 72b6972d74SZach Atkins 73b6972d74SZach Atkins // -- Ref coordinates 74b6972d74SZach Atkins { 75b6972d74SZach Atkins PetscMemType X_mem_type; 76b6972d74SZach Atkins const PetscScalar *x; 77b6972d74SZach Atkins 78b6972d74SZach Atkins CeedVectorCreate((*ctx)->ceed, num_points * dim, &x_ref_points); 79b6972d74SZach Atkins 80b6972d74SZach Atkins PetscCall(VecGetArrayReadAndMemType(X_ref, (const PetscScalar **)&x, &X_mem_type)); 81b6972d74SZach Atkins CeedVectorSetArray(x_ref_points, MemTypeP2C(X_mem_type), CEED_COPY_VALUES, (CeedScalar *)x); 82b6972d74SZach Atkins PetscCall(VecRestoreArrayReadAndMemType(X_ref, (const PetscScalar **)&x)); 83b6972d74SZach Atkins } 84b6972d74SZach Atkins 85b6972d74SZach Atkins // Create Q data 86b6972d74SZach Atkins { 87b6972d74SZach Atkins CeedQFunction qf_setup; 88b6972d74SZach Atkins CeedOperator op_setup; 89b6972d74SZach Atkins CeedVector x_coord; 90b6972d74SZach Atkins 91b6972d74SZach Atkins { 92b6972d74SZach Atkins Vec X_loc; 93b6972d74SZach Atkins CeedInt len; 94b6972d74SZach Atkins const PetscScalar *x; 95b6972d74SZach Atkins 96b6972d74SZach Atkins PetscCall(DMGetCoordinatesLocal(dm_mesh, &X_loc)); 97b6972d74SZach Atkins PetscCall(VecGetLocalSize(X_loc, &len)); 98b6972d74SZach Atkins CeedVectorCreate((*ctx)->ceed, len, &x_coord); 99b6972d74SZach Atkins 100b6972d74SZach Atkins PetscCall(VecGetArrayRead(X_loc, &x)); 101b6972d74SZach Atkins CeedVectorSetArray(x_coord, CEED_MEM_HOST, CEED_COPY_VALUES, (CeedScalar *)x); 102b6972d74SZach Atkins PetscCall(VecRestoreArrayRead(X_loc, &x)); 103b6972d74SZach Atkins } 104b6972d74SZach Atkins 105b6972d74SZach Atkins // Setup geometric scaling 106b6972d74SZach Atkins CeedQFunctionCreateInterior((*ctx)->ceed, 1, SetupMass, SetupMass_loc, &qf_setup); 107b6972d74SZach Atkins CeedQFunctionAddInput(qf_setup, "x", dim * dim, CEED_EVAL_GRAD); 108b6972d74SZach Atkins CeedQFunctionAddInput(qf_setup, "weight", 1, CEED_EVAL_WEIGHT); 109b6972d74SZach Atkins CeedQFunctionAddOutput(qf_setup, "rho", 1, CEED_EVAL_NONE); 110b6972d74SZach Atkins 111b6972d74SZach Atkins CeedOperatorCreateAtPoints((*ctx)->ceed, qf_setup, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setup); 112b6972d74SZach Atkins CeedOperatorSetField(op_setup, "x", elem_restr_x_mesh, basis_x, CEED_VECTOR_ACTIVE); 113b6972d74SZach Atkins CeedOperatorSetField(op_setup, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE); 114b6972d74SZach Atkins CeedOperatorSetField(op_setup, "rho", elem_restr_q_data_points, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE); 115b6972d74SZach Atkins CeedOperatorAtPointsSetPoints(op_setup, elem_restr_x_points, x_ref_points); 116b6972d74SZach Atkins 117b6972d74SZach Atkins CeedOperatorApply(op_setup, x_coord, q_data_points, CEED_REQUEST_IMMEDIATE); 118b6972d74SZach Atkins 119b6972d74SZach Atkins // Cleanup 120b6972d74SZach Atkins CeedVectorDestroy(&x_coord); 121b6972d74SZach Atkins CeedQFunctionDestroy(&qf_setup); 122b6972d74SZach Atkins CeedOperatorDestroy(&op_setup); 123b6972d74SZach Atkins } 124b6972d74SZach Atkins 125b6972d74SZach Atkins // -- Cleanup 126b6972d74SZach Atkins PetscCall(ISDestroy(&is_points)); 127b6972d74SZach Atkins PetscCall(VecDestroy(&X_ref)); 128b6972d74SZach Atkins } 129b6972d74SZach Atkins 130b6972d74SZach Atkins PetscCall(DMSetApplicationContext(dm_mesh, (void *)(*ctx))); 131b6972d74SZach Atkins 132b6972d74SZach Atkins // Create operators 133b6972d74SZach Atkins // Mesh to points interpolation operator 134b6972d74SZach Atkins { 135b6972d74SZach Atkins CeedQFunction qf_mesh_to_points; 136b6972d74SZach Atkins 137b6972d74SZach Atkins // -- Create operator 138b6972d74SZach Atkins CeedQFunctionCreateIdentity((*ctx)->ceed, num_comp, CEED_EVAL_INTERP, CEED_EVAL_NONE, &qf_mesh_to_points); 139b6972d74SZach Atkins 140b6972d74SZach Atkins CeedOperatorCreateAtPoints((*ctx)->ceed, qf_mesh_to_points, NULL, NULL, &(*ctx)->op_mesh_to_points); 141b6972d74SZach Atkins CeedOperatorSetField((*ctx)->op_mesh_to_points, "input", elem_restr_u_mesh, basis_u, CEED_VECTOR_ACTIVE); 142b6972d74SZach Atkins CeedOperatorSetField((*ctx)->op_mesh_to_points, "output", elem_restr_u_points, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE); 143b6972d74SZach Atkins CeedOperatorAtPointsSetPoints((*ctx)->op_mesh_to_points, elem_restr_x_points, x_ref_points); 144b6972d74SZach Atkins 145b6972d74SZach Atkins // -- Cleanup 146b6972d74SZach Atkins CeedQFunctionDestroy(&qf_mesh_to_points); 147b6972d74SZach Atkins } 148b6972d74SZach Atkins 149b6972d74SZach Atkins // RHS operator 150b6972d74SZach Atkins { 151b6972d74SZach Atkins CeedQFunction qf_pts_to_mesh; 152b6972d74SZach Atkins CeedQFunctionContext qf_ctx; 153b6972d74SZach Atkins 154b6972d74SZach Atkins // -- Mass QFunction 155b6972d74SZach Atkins CeedQFunctionCreateInterior((*ctx)->ceed, 1, Mass, Mass_loc, &qf_pts_to_mesh); 156b6972d74SZach Atkins CeedQFunctionAddInput(qf_pts_to_mesh, "q data", 1, CEED_EVAL_NONE); 157b6972d74SZach Atkins CeedQFunctionAddInput(qf_pts_to_mesh, "u", num_comp, CEED_EVAL_NONE); 158b6972d74SZach Atkins CeedQFunctionAddOutput(qf_pts_to_mesh, "v", num_comp, CEED_EVAL_INTERP); 159b6972d74SZach Atkins 160b6972d74SZach Atkins // -- QFunction context 161b6972d74SZach Atkins CeedQFunctionContextCreate((*ctx)->ceed, &qf_ctx); 162b6972d74SZach Atkins CeedQFunctionContextSetData(qf_ctx, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof(num_comp), &num_comp); 163b6972d74SZach Atkins CeedQFunctionSetContext(qf_pts_to_mesh, qf_ctx); 164b6972d74SZach Atkins 165b6972d74SZach Atkins // -- Mass Operator 166b6972d74SZach Atkins CeedOperatorCreateAtPoints((*ctx)->ceed, qf_pts_to_mesh, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &(*ctx)->op_points_to_mesh); 167b6972d74SZach Atkins CeedOperatorSetField((*ctx)->op_points_to_mesh, "q data", elem_restr_q_data_points, CEED_BASIS_NONE, q_data_points); 168b6972d74SZach Atkins CeedOperatorSetField((*ctx)->op_points_to_mesh, "u", elem_restr_u_points, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE); 169b6972d74SZach Atkins CeedOperatorSetField((*ctx)->op_points_to_mesh, "v", elem_restr_u_mesh, basis_u, CEED_VECTOR_ACTIVE); 170b6972d74SZach Atkins CeedOperatorAtPointsSetPoints((*ctx)->op_points_to_mesh, elem_restr_x_points, x_ref_points); 171b6972d74SZach Atkins 172b6972d74SZach Atkins // -- Cleanup 173b6972d74SZach Atkins CeedQFunctionContextDestroy(&qf_ctx); 174b6972d74SZach Atkins CeedQFunctionDestroy(&qf_pts_to_mesh); 175b6972d74SZach Atkins } 176b6972d74SZach Atkins 177b6972d74SZach Atkins // Mass operator 178b6972d74SZach Atkins { 179b6972d74SZach Atkins CeedQFunction qf_mass; 180b6972d74SZach Atkins CeedQFunctionContext ctx_mass; 181b6972d74SZach Atkins 182b6972d74SZach Atkins // -- Mass QFunction 183b6972d74SZach Atkins CeedQFunctionCreateInterior((*ctx)->ceed, 1, Mass, Mass_loc, &qf_mass); 184b6972d74SZach Atkins CeedQFunctionAddInput(qf_mass, "q data", 1, CEED_EVAL_NONE); 185b6972d74SZach Atkins CeedQFunctionAddInput(qf_mass, "u", num_comp, CEED_EVAL_INTERP); 186b6972d74SZach Atkins CeedQFunctionAddOutput(qf_mass, "v", num_comp, CEED_EVAL_INTERP); 187b6972d74SZach Atkins 188b6972d74SZach Atkins // -- QFunction context 189b6972d74SZach Atkins CeedQFunctionContextCreate((*ctx)->ceed, &ctx_mass); 190b6972d74SZach Atkins CeedQFunctionContextSetData(ctx_mass, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof(num_comp), &num_comp); 191b6972d74SZach Atkins CeedQFunctionSetContext(qf_mass, ctx_mass); 192b6972d74SZach Atkins 193b6972d74SZach Atkins // -- Mass Operator 194b6972d74SZach Atkins CeedOperatorCreateAtPoints((*ctx)->ceed, qf_mass, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &(*ctx)->op_mass); 195b6972d74SZach Atkins CeedOperatorSetField((*ctx)->op_mass, "q data", elem_restr_q_data_points, CEED_BASIS_NONE, q_data_points); 196b6972d74SZach Atkins CeedOperatorSetField((*ctx)->op_mass, "u", elem_restr_u_mesh, basis_u, CEED_VECTOR_ACTIVE); 197b6972d74SZach Atkins CeedOperatorSetField((*ctx)->op_mass, "v", elem_restr_u_mesh, basis_u, CEED_VECTOR_ACTIVE); 198b6972d74SZach Atkins CeedOperatorAtPointsSetPoints((*ctx)->op_mass, elem_restr_x_points, x_ref_points); 199b6972d74SZach Atkins 200b6972d74SZach Atkins // -- Cleanup 201b6972d74SZach Atkins CeedQFunctionContextDestroy(&ctx_mass); 202b6972d74SZach Atkins CeedQFunctionDestroy(&qf_mass); 203b6972d74SZach Atkins } 204b6972d74SZach Atkins 205b6972d74SZach Atkins // Cleanup 206b6972d74SZach Atkins CeedElemRestrictionDestroy(&elem_restr_u_mesh); 207b6972d74SZach Atkins CeedElemRestrictionDestroy(&elem_restr_x_mesh); 208b6972d74SZach Atkins CeedElemRestrictionDestroy(&elem_restr_u_points); 209b6972d74SZach Atkins CeedElemRestrictionDestroy(&elem_restr_x_points); 210b6972d74SZach Atkins CeedElemRestrictionDestroy(&elem_restr_q_data_points); 211b6972d74SZach Atkins CeedBasisDestroy(&basis_u); 212b6972d74SZach Atkins CeedBasisDestroy(&basis_x); 213b6972d74SZach Atkins CeedVectorDestroy(&x_ref_points); 214b6972d74SZach Atkins CeedVectorDestroy(&q_data_points); 215b6972d74SZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 216b6972d74SZach Atkins } 217b6972d74SZach Atkins 218b6972d74SZach Atkins PetscErrorCode DMSwarmCeedContextDestroy(DMSwarmCeedContext *ctx) { 219b6972d74SZach Atkins PetscFunctionBeginUser; 220b6972d74SZach Atkins CeedDestroy(&(*ctx)->ceed); 221b6972d74SZach Atkins CeedVectorDestroy(&(*ctx)->u_mesh); 222b6972d74SZach Atkins CeedVectorDestroy(&(*ctx)->v_mesh); 223b6972d74SZach Atkins CeedVectorDestroy(&(*ctx)->u_points); 224b6972d74SZach Atkins CeedOperatorDestroy(&(*ctx)->op_mesh_to_points); 225b6972d74SZach Atkins CeedOperatorDestroy(&(*ctx)->op_points_to_mesh); 226b6972d74SZach Atkins CeedOperatorDestroy(&(*ctx)->op_mass); 227b6972d74SZach Atkins PetscCall(PetscFree(*ctx)); 228b6972d74SZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 229b6972d74SZach Atkins } 230b6972d74SZach Atkins 231b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 232b6972d74SZach Atkins // PETSc-libCEED memory space utilities 233b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 234b6972d74SZach Atkins PetscErrorCode DMSwarmPICFieldP2C(DM dm_swarm, const char *field, CeedVector x_ceed) { 235b6972d74SZach Atkins PetscScalar *x; 236b6972d74SZach Atkins 237b6972d74SZach Atkins PetscFunctionBeginUser; 238b6972d74SZach Atkins PetscCall(DMSwarmGetField(dm_swarm, field, NULL, NULL, (void **)&x)); 239b6972d74SZach Atkins CeedVectorSetArray(x_ceed, CEED_MEM_HOST, CEED_USE_POINTER, (CeedScalar *)x); 240b6972d74SZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 241b6972d74SZach Atkins } 242b6972d74SZach Atkins 243b6972d74SZach Atkins PetscErrorCode DMSwarmPICFieldC2P(DM dm_swarm, const char *field, CeedVector x_ceed) { 244b6972d74SZach Atkins PetscScalar *x; 245b6972d74SZach Atkins 246b6972d74SZach Atkins PetscFunctionBeginUser; 247b6972d74SZach Atkins CeedVectorTakeArray(x_ceed, CEED_MEM_HOST, (CeedScalar **)&x); 248b6972d74SZach Atkins PetscCall(DMSwarmRestoreField(dm_swarm, field, NULL, NULL, (void **)&x)); 249b6972d74SZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 250b6972d74SZach Atkins } 251b6972d74SZach Atkins 252b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 253b6972d74SZach Atkins // Swarm point location utility 254b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 255b6972d74SZach Atkins PetscErrorCode DMSwarmInitalizePointLocations(DM dm_swarm, PointSwarmType point_swarm_type, PetscInt num_points, PetscInt num_points_per_cell) { 256b6972d74SZach Atkins PetscFunctionBeginUser; 257b6972d74SZach Atkins switch (point_swarm_type) { 258b6972d74SZach Atkins case SWARM_GAUSS: 259b6972d74SZach Atkins case SWARM_UNIFORM: { 260b6972d74SZach Atkins // -- Set gauss or uniform point locations in each cell 261b6972d74SZach Atkins PetscInt num_points_per_cell_1d = round(cbrt(num_points_per_cell * 1.0)), dim = 3; 262b6972d74SZach Atkins PetscScalar point_coords[num_points_per_cell * 3]; 263b6972d74SZach Atkins CeedScalar points_1d[num_points_per_cell_1d], weights_1d[num_points_per_cell_1d]; 264b6972d74SZach Atkins 265b6972d74SZach Atkins if (point_swarm_type == SWARM_GAUSS) { 266b6972d74SZach Atkins PetscCall(CeedGaussQuadrature(num_points_per_cell_1d, points_1d, weights_1d)); 267b6972d74SZach Atkins } else { 268b6972d74SZach Atkins for (PetscInt i = 0; i < num_points_per_cell_1d; i++) points_1d[i] = 2.0 * (PetscReal)(i + 1) / (PetscReal)(num_points_per_cell_1d + 1) - 1; 269b6972d74SZach Atkins } 270b6972d74SZach Atkins for (PetscInt i = 0; i < num_points_per_cell_1d; i++) { 271b6972d74SZach Atkins for (PetscInt j = 0; j < num_points_per_cell_1d; j++) { 272b6972d74SZach Atkins for (PetscInt k = 0; k < num_points_per_cell_1d; k++) { 273b6972d74SZach Atkins PetscInt p = (i * num_points_per_cell_1d + j) * num_points_per_cell_1d + k; 274b6972d74SZach Atkins 275b6972d74SZach Atkins point_coords[p * dim + 0] = points_1d[i]; 276b6972d74SZach Atkins point_coords[p * dim + 1] = points_1d[j]; 277b6972d74SZach Atkins point_coords[p * dim + 2] = points_1d[k]; 278b6972d74SZach Atkins } 279b6972d74SZach Atkins } 280b6972d74SZach Atkins } 281b6972d74SZach Atkins PetscCall(DMSwarmSetPointCoordinatesCellwise(dm_swarm, num_points_per_cell_1d * num_points_per_cell_1d * num_points_per_cell_1d, point_coords)); 282b6972d74SZach Atkins } break; 283b6972d74SZach Atkins case SWARM_CELL_RANDOM: { 284b6972d74SZach Atkins // -- Set points randomly in each cell 285b6972d74SZach Atkins PetscInt dim = 3, num_cells_total = 1, num_cells[] = {1, 1, 1}; 286b6972d74SZach Atkins PetscScalar *point_coords; 287b6972d74SZach Atkins PetscRandom rng; 288b6972d74SZach Atkins 289b6972d74SZach Atkins PetscOptionsBegin(PetscObjectComm((PetscObject)dm_swarm), NULL, "libCEED example using PETSc with DMSwarm", NULL); 290b6972d74SZach Atkins 291b6972d74SZach Atkins PetscCall(PetscOptionsInt("-dm_plex_dim", "Background mesh dimension", NULL, dim, &dim, NULL)); 292b6972d74SZach Atkins PetscCall(PetscOptionsIntArray("-dm_plex_box_faces", "Number of cells", NULL, num_cells, &dim, NULL)); 293b6972d74SZach Atkins 294b6972d74SZach Atkins PetscOptionsEnd(); 295b6972d74SZach Atkins 296b6972d74SZach Atkins PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm_swarm), &rng)); 297b6972d74SZach Atkins 298b6972d74SZach Atkins num_cells_total = num_cells[0] * num_cells[1] * num_cells[2]; 299b6972d74SZach Atkins PetscCall(DMSwarmGetField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&point_coords)); 300b6972d74SZach Atkins for (PetscInt c = 0; c < num_cells_total; c++) { 301b6972d74SZach Atkins PetscInt cell_index[3] = {c % num_cells[0], (c / num_cells[0]) % num_cells[1], (c / num_cells[0] / num_cells[1]) % num_cells[2]}; 302b6972d74SZach Atkins 303b6972d74SZach Atkins for (PetscInt p = 0; p < num_points_per_cell; p++) { 304b6972d74SZach Atkins PetscInt point_index = c * num_points_per_cell + p; 305b6972d74SZach Atkins PetscScalar random_value; 306b6972d74SZach Atkins 307b6972d74SZach Atkins for (PetscInt i = 0; i < dim; i++) { 308b6972d74SZach Atkins PetscCall(PetscRandomGetValue(rng, &random_value)); 309b6972d74SZach Atkins point_coords[point_index * dim + i] = -1.0 + cell_index[i] * 2.0 / (num_cells[i] + 1.0) + random_value; 310b6972d74SZach Atkins } 311b6972d74SZach Atkins } 312b6972d74SZach Atkins } 313b6972d74SZach Atkins PetscCall(DMSwarmRestoreField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&point_coords)); 314b6972d74SZach Atkins PetscCall(PetscRandomDestroy(&rng)); 315b6972d74SZach Atkins } break; 316b6972d74SZach Atkins case SWARM_SINUSOIDAL: { 317b6972d74SZach Atkins // -- Set points distributed per sinusoidal functions 318b6972d74SZach Atkins PetscInt dim = 3; 319b6972d74SZach Atkins PetscScalar *point_coords; 320b6972d74SZach Atkins DM dm_mesh; 321b6972d74SZach Atkins 322b6972d74SZach Atkins PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh)); 323b6972d74SZach Atkins PetscCall(DMGetDimension(dm_mesh, &dim)); 324b6972d74SZach Atkins 325b6972d74SZach Atkins PetscCall(DMSwarmGetField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&point_coords)); 326b6972d74SZach Atkins for (PetscInt p = 0; p < num_points; p++) { 327b6972d74SZach Atkins point_coords[p * dim + 0] = -PetscCosReal((PetscReal)(p + 1) / (PetscReal)(num_points + 1) * PETSC_PI); 328b6972d74SZach Atkins if (dim > 1) point_coords[p * dim + 1] = -PetscSinReal((PetscReal)(p + 1) / (PetscReal)(num_points + 1) * PETSC_PI); 329b6972d74SZach Atkins if (dim > 2) point_coords[p * dim + 2] = PetscSinReal((PetscReal)(p + 1) / (PetscReal)(num_points + 1) * PETSC_PI); 330b6972d74SZach Atkins } 331b6972d74SZach Atkins PetscCall(DMSwarmRestoreField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&point_coords)); 332b6972d74SZach Atkins } break; 333b6972d74SZach Atkins } 334b6972d74SZach Atkins PetscCall(DMSwarmMigrate(dm_swarm, PETSC_TRUE)); 335b6972d74SZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 336b6972d74SZach Atkins } 337b6972d74SZach Atkins 338b6972d74SZach Atkins /*@C 339b6972d74SZach Atkins DMSwarmCreateReferenceCoordinates - Compute the cell reference coordinates for local DMSwarm points. 340b6972d74SZach Atkins 341b6972d74SZach Atkins Collective 342b6972d74SZach Atkins 343b6972d74SZach Atkins Input Parameter: 344b6972d74SZach Atkins . dm_swarm - the `DMSwarm` 345b6972d74SZach Atkins 346b6972d74SZach Atkins Output Parameters: 347b6972d74SZach Atkins + is_points - The IS object for indexing into points per cell 348b6972d74SZach Atkins - X_points_ref - Vec holding the cell reference coordinates for local DMSwarm points 349b6972d74SZach Atkins 350b6972d74SZach Atkins The index set contains ranges of indices for each local cell. This range contains the indices of every point in the cell. 351b6972d74SZach Atkins 352b6972d74SZach Atkins ``` 353b6972d74SZach Atkins total_num_cells 354b6972d74SZach Atkins cell_0_start_index 355b6972d74SZach Atkins cell_1_start_index 356b6972d74SZach Atkins cell_2_start_index 357b6972d74SZach Atkins ... 358b6972d74SZach Atkins cell_n_start_index 359b6972d74SZach Atkins cell_n_stop_index 360b6972d74SZach Atkins cell_0_point_0 361b6972d74SZach Atkins cell_0_point_0 362b6972d74SZach Atkins ... 363b6972d74SZach Atkins cell_n_point_m 364b6972d74SZach Atkins ``` 365b6972d74SZach Atkins 366b6972d74SZach Atkins Level: beginner 367b6972d74SZach Atkins 368b6972d74SZach Atkins .seealso: `DMSwarm` 369b6972d74SZach Atkins @*/ 370b6972d74SZach Atkins PetscErrorCode DMSwarmCreateReferenceCoordinates(DM dm_swarm, IS *is_points, Vec *X_points_ref) { 371b6972d74SZach Atkins PetscInt cell_start, cell_end, num_cells_local, dim, num_points_local, *cell_points, points_offset; 372b6972d74SZach Atkins PetscScalar *coords_points_ref; 373b6972d74SZach Atkins const PetscScalar *coords_points_true; 374b6972d74SZach Atkins DM dm_mesh; 375b6972d74SZach Atkins 376b6972d74SZach Atkins PetscFunctionBeginUser; 377b6972d74SZach Atkins PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh)); 378b6972d74SZach Atkins 379b6972d74SZach Atkins // Create vector to hold reference coordinates 380b6972d74SZach Atkins { 381b6972d74SZach Atkins Vec X_points_true; 382b6972d74SZach Atkins 383b6972d74SZach Atkins PetscCall(DMSwarmCreateLocalVectorFromField(dm_swarm, DMSwarmPICField_coor, &X_points_true)); 384b6972d74SZach Atkins PetscCall(VecDuplicate(X_points_true, X_points_ref)); 385b6972d74SZach Atkins PetscCall(DMSwarmDestroyLocalVectorFromField(dm_swarm, DMSwarmPICField_coor, &X_points_true)); 386b6972d74SZach Atkins } 387b6972d74SZach Atkins 388b6972d74SZach Atkins // Allocate index set array 389b6972d74SZach Atkins PetscCall(DMPlexGetHeightStratum(dm_mesh, 0, &cell_start, &cell_end)); 390b6972d74SZach Atkins num_cells_local = cell_end - cell_start; 391b6972d74SZach Atkins points_offset = num_cells_local + 2; 392b6972d74SZach Atkins PetscCall(VecGetLocalSize(*X_points_ref, &num_points_local)); 393b6972d74SZach Atkins PetscCall(DMGetDimension(dm_mesh, &dim)); 394b6972d74SZach Atkins num_points_local /= dim; 395b6972d74SZach Atkins PetscCall(PetscMalloc1(num_points_local + num_cells_local + 2, &cell_points)); 396b6972d74SZach Atkins cell_points[0] = num_cells_local; 397b6972d74SZach Atkins 398b6972d74SZach Atkins // Get reference coordinates for each swarm point wrt the elements in the background mesh 399b6972d74SZach Atkins PetscCall(DMSwarmSortGetAccess(dm_swarm)); 400b6972d74SZach Atkins PetscCall(DMSwarmGetField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords_points_true)); 401b6972d74SZach Atkins PetscCall(VecGetArray(*X_points_ref, &coords_points_ref)); 402b6972d74SZach Atkins for (PetscInt cell = cell_start, num_points_processed = 0; cell < cell_end; cell++) { 403b6972d74SZach Atkins PetscInt *points_in_cell, num_points_in_cell, local_cell = cell - cell_start; 404b6972d74SZach Atkins PetscReal v[3], J[9], invJ[9], detJ, v0_ref[3] = {-1.0, -1.0, -1.0}; 405b6972d74SZach Atkins 406b6972d74SZach Atkins PetscCall(DMSwarmSortGetPointsPerCell(dm_swarm, cell, &num_points_in_cell, &points_in_cell)); 407b6972d74SZach Atkins // -- Reference coordinates for swarm points in background mesh element 408b6972d74SZach Atkins PetscCall(DMPlexComputeCellGeometryFEM(dm_mesh, cell, NULL, v, J, invJ, &detJ)); 409b6972d74SZach Atkins cell_points[local_cell + 1] = num_points_processed + points_offset; 410b6972d74SZach Atkins for (PetscInt p = 0; p < num_points_in_cell; p++) { 411b6972d74SZach Atkins const CeedInt point = points_in_cell[p]; 412b6972d74SZach Atkins 413b6972d74SZach Atkins cell_points[points_offset + (num_points_processed++)] = point; 414b6972d74SZach Atkins CoordinatesRealToRef(dim, dim, v0_ref, v, invJ, &coords_points_true[point * dim], &coords_points_ref[point * dim]); 415b6972d74SZach Atkins } 416b6972d74SZach Atkins 417b6972d74SZach Atkins // -- Cleanup 418b6972d74SZach Atkins PetscCall(PetscFree(points_in_cell)); 419b6972d74SZach Atkins } 420b6972d74SZach Atkins cell_points[points_offset - 1] = num_points_local + points_offset; 421b6972d74SZach Atkins 422b6972d74SZach Atkins // Cleanup 423b6972d74SZach Atkins PetscCall(DMSwarmRestoreField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords_points_true)); 424b6972d74SZach Atkins PetscCall(VecRestoreArray(*X_points_ref, &coords_points_ref)); 425b6972d74SZach Atkins PetscCall(DMSwarmSortRestoreAccess(dm_swarm)); 426b6972d74SZach Atkins 427b6972d74SZach Atkins // Create index set 428b6972d74SZach Atkins PetscCall(ISCreateGeneral(PETSC_COMM_SELF, num_points_local + points_offset, cell_points, PETSC_OWN_POINTER, is_points)); 429b6972d74SZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 430b6972d74SZach Atkins } 431b6972d74SZach Atkins 432b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 433b6972d74SZach Atkins // RHS for Swarm to Mesh projection 434b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 435b6972d74SZach Atkins PetscErrorCode DMSwarmCreateProjectionRHS(DM dm_swarm, const char *field, Vec B_mesh) { 436b6972d74SZach Atkins PetscMemType B_mem_type; 437b6972d74SZach Atkins DM dm_mesh; 438b6972d74SZach Atkins Vec B_mesh_loc; 439b6972d74SZach Atkins DMSwarmCeedContext swarm_ceed_context; 440b6972d74SZach Atkins 441b6972d74SZach Atkins PetscFunctionBeginUser; 442b6972d74SZach Atkins // Get mesh DM 443b6972d74SZach Atkins PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh)); 444b6972d74SZach Atkins PetscCall(DMGetApplicationContext(dm_mesh, (void *)&swarm_ceed_context)); 445b6972d74SZach Atkins 446b6972d74SZach Atkins // Get mesh values 447b6972d74SZach Atkins PetscCall(DMGetLocalVector(dm_mesh, &B_mesh_loc)); 448b6972d74SZach Atkins PetscCall(VecZeroEntries(B_mesh_loc)); 449b6972d74SZach Atkins PetscCall(VecP2C(B_mesh_loc, &B_mem_type, swarm_ceed_context->v_mesh)); 450b6972d74SZach Atkins 451b6972d74SZach Atkins // Get swarm access 452b6972d74SZach Atkins PetscCall(DMSwarmSortGetAccess(dm_swarm)); 453b6972d74SZach Atkins PetscCall(DMSwarmPICFieldP2C(dm_swarm, field, swarm_ceed_context->u_points)); 454b6972d74SZach Atkins 455b6972d74SZach Atkins // Interpolate field from swarm points to mesh 456b6972d74SZach Atkins CeedOperatorApply(swarm_ceed_context->op_points_to_mesh, swarm_ceed_context->u_points, swarm_ceed_context->v_mesh, CEED_REQUEST_IMMEDIATE); 457b6972d74SZach Atkins 458b6972d74SZach Atkins // Restore PETSc Vecs and Local to Global 459b6972d74SZach Atkins PetscCall(VecC2P(swarm_ceed_context->v_mesh, B_mem_type, B_mesh_loc)); 460b6972d74SZach Atkins PetscCall(VecZeroEntries(B_mesh)); 461b6972d74SZach Atkins PetscCall(DMLocalToGlobal(dm_mesh, B_mesh_loc, ADD_VALUES, B_mesh)); 462b6972d74SZach Atkins 463b6972d74SZach Atkins // Cleanup 464b6972d74SZach Atkins PetscCall(DMSwarmPICFieldC2P(dm_swarm, field, swarm_ceed_context->u_points)); 465b6972d74SZach Atkins PetscCall(DMSwarmSortRestoreAccess(dm_swarm)); 466b6972d74SZach Atkins PetscCall(DMRestoreLocalVector(dm_mesh, &B_mesh_loc)); 467b6972d74SZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 468b6972d74SZach Atkins } 469b6972d74SZach Atkins 470b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 471b6972d74SZach Atkins // Swarm "mass matrix" 472b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 473b6972d74SZach Atkins PetscErrorCode MatMult_SwarmMass(Mat A, Vec U_mesh, Vec V_mesh) { 474b6972d74SZach Atkins PetscMemType U_mem_type, V_mem_type; 475b6972d74SZach Atkins DM dm_mesh; 476b6972d74SZach Atkins Vec U_mesh_loc, V_mesh_loc; 477b6972d74SZach Atkins DMSwarmCeedContext swarm_ceed_context; 478b6972d74SZach Atkins 479b6972d74SZach Atkins PetscFunctionBeginUser; 480b6972d74SZach Atkins // Get mesh DM 481b6972d74SZach Atkins PetscCall(MatGetDM(A, &dm_mesh)); 482b6972d74SZach Atkins PetscCall(DMGetApplicationContext(dm_mesh, (void *)&swarm_ceed_context)); 483b6972d74SZach Atkins 484b6972d74SZach Atkins // Global to Local and get PETSc Vec access 485b6972d74SZach Atkins PetscCall(DMGetLocalVector(dm_mesh, &U_mesh_loc)); 486b6972d74SZach Atkins PetscCall(VecZeroEntries(U_mesh_loc)); 487b6972d74SZach Atkins PetscCall(DMGlobalToLocal(dm_mesh, U_mesh, INSERT_VALUES, U_mesh_loc)); 488b6972d74SZach Atkins PetscCall(VecReadP2C(U_mesh_loc, &U_mem_type, swarm_ceed_context->u_mesh)); 489b6972d74SZach Atkins PetscCall(DMGetLocalVector(dm_mesh, &V_mesh_loc)); 490b6972d74SZach Atkins PetscCall(VecZeroEntries(V_mesh_loc)); 491b6972d74SZach Atkins PetscCall(VecP2C(V_mesh_loc, &V_mem_type, swarm_ceed_context->v_mesh)); 492b6972d74SZach Atkins 493b6972d74SZach Atkins // Apply swarm mass operator 494b6972d74SZach Atkins CeedOperatorApply(swarm_ceed_context->op_mass, swarm_ceed_context->u_mesh, swarm_ceed_context->v_mesh, CEED_REQUEST_IMMEDIATE); 495b6972d74SZach Atkins 496b6972d74SZach Atkins // Restore PETSc Vecs and Local to Global 497b6972d74SZach Atkins PetscCall(VecReadC2P(swarm_ceed_context->u_mesh, U_mem_type, U_mesh_loc)); 498b6972d74SZach Atkins PetscCall(VecC2P(swarm_ceed_context->v_mesh, V_mem_type, V_mesh_loc)); 499b6972d74SZach Atkins PetscCall(VecZeroEntries(V_mesh)); 500b6972d74SZach Atkins PetscCall(DMLocalToGlobal(dm_mesh, V_mesh_loc, ADD_VALUES, V_mesh)); 501b6972d74SZach Atkins 502b6972d74SZach Atkins // Cleanup 503b6972d74SZach Atkins PetscCall(DMRestoreLocalVector(dm_mesh, &U_mesh_loc)); 504b6972d74SZach Atkins PetscCall(DMRestoreLocalVector(dm_mesh, &V_mesh_loc)); 505b6972d74SZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 506b6972d74SZach Atkins } 507b6972d74SZach Atkins 508b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 509b6972d74SZach Atkins // Swarm to mesh projection 510b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 511b6972d74SZach Atkins PetscErrorCode DMSwarmProjectFromSwarmToCells(DM dm_swarm, const char *field, Vec U_mesh) { 512b6972d74SZach Atkins PetscBool test_mode; 513b6972d74SZach Atkins Vec B_mesh; 514b6972d74SZach Atkins Mat M; 515b6972d74SZach Atkins KSP ksp; 516b6972d74SZach Atkins DM dm_mesh; 517b6972d74SZach Atkins DMSwarmCeedContext swarm_ceed_context; 518b6972d74SZach Atkins MPI_Comm comm; 519b6972d74SZach Atkins 520b6972d74SZach Atkins PetscFunctionBeginUser; 521b6972d74SZach Atkins PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swarm-to-Mesh Projection Options", NULL); 522b6972d74SZach Atkins PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, test_mode, &test_mode, NULL)); 523b6972d74SZach Atkins PetscOptionsEnd(); 524b6972d74SZach Atkins 525b6972d74SZach Atkins comm = PetscObjectComm((PetscObject)dm_swarm); 526b6972d74SZach Atkins PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh)); 527b6972d74SZach Atkins PetscCall(DMGetApplicationContext(dm_mesh, (void *)&swarm_ceed_context)); 528b6972d74SZach Atkins PetscCall(VecDuplicate(U_mesh, &B_mesh)); 529b6972d74SZach Atkins 530b6972d74SZach Atkins // Setup "mass matrix" 531b6972d74SZach Atkins { 532b6972d74SZach Atkins PetscInt l_size, g_size; 533b6972d74SZach Atkins 534b6972d74SZach Atkins PetscCall(VecGetLocalSize(U_mesh, &l_size)); 535b6972d74SZach Atkins PetscCall(VecGetSize(U_mesh, &g_size)); 536b6972d74SZach Atkins PetscCall(MatCreateShell(comm, l_size, l_size, g_size, g_size, swarm_ceed_context, &M)); 537b6972d74SZach Atkins PetscCall(MatSetDM(M, dm_mesh)); 538b6972d74SZach Atkins PetscCall(MatShellSetOperation(M, MATOP_MULT, (void (*)(void))MatMult_SwarmMass)); 539b6972d74SZach Atkins } 540b6972d74SZach Atkins 541b6972d74SZach Atkins // Setup KSP 542b6972d74SZach Atkins { 543b6972d74SZach Atkins PC pc; 544b6972d74SZach Atkins 545b6972d74SZach Atkins PetscCall(KSPCreate(comm, &ksp)); 546b6972d74SZach Atkins PetscCall(KSPGetPC(ksp, &pc)); 547b6972d74SZach Atkins PetscCall(PCSetType(pc, PCJACOBI)); 548b6972d74SZach Atkins PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 549b6972d74SZach Atkins PetscCall(KSPSetType(ksp, KSPCG)); 550b6972d74SZach Atkins PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL)); 551b6972d74SZach Atkins PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 552b6972d74SZach Atkins PetscCall(KSPSetOperators(ksp, M, M)); 553b6972d74SZach Atkins PetscCall(KSPSetFromOptions(ksp)); 554b6972d74SZach Atkins PetscCall(PetscObjectSetName((PetscObject)ksp, "Swarm-to-Mesh Projection")); 555b6972d74SZach Atkins PetscCall(KSPViewFromOptions(ksp, NULL, "-ksp_projection_view")); 556b6972d74SZach Atkins } 557b6972d74SZach Atkins 558b6972d74SZach Atkins // Setup RHS 559b6972d74SZach Atkins PetscCall(DMSwarmCreateProjectionRHS(dm_swarm, field, B_mesh)); 560b6972d74SZach Atkins 561b6972d74SZach Atkins // Solve 562b6972d74SZach Atkins PetscCall(VecZeroEntries(U_mesh)); 563b6972d74SZach Atkins PetscCall(KSPSolve(ksp, B_mesh, U_mesh)); 564b6972d74SZach Atkins 565b6972d74SZach Atkins // KSP summary 566b6972d74SZach Atkins { 567b6972d74SZach Atkins KSPType ksp_type; 568b6972d74SZach Atkins KSPConvergedReason reason; 569b6972d74SZach Atkins PetscReal rnorm; 570b6972d74SZach Atkins PetscInt its; 571b6972d74SZach Atkins PetscCall(KSPGetType(ksp, &ksp_type)); 572b6972d74SZach Atkins PetscCall(KSPGetConvergedReason(ksp, &reason)); 573b6972d74SZach Atkins PetscCall(KSPGetIterationNumber(ksp, &its)); 574b6972d74SZach Atkins PetscCall(KSPGetResidualNorm(ksp, &rnorm)); 575b6972d74SZach Atkins 576b6972d74SZach Atkins if (!test_mode || reason < 0 || rnorm > 1e-8) { 577b6972d74SZach Atkins PetscCall(PetscPrintf(comm, 578b6972d74SZach Atkins "Swarm-to-Mesh Projection KSP Solve:\n" 579b6972d74SZach Atkins " KSP type: %s\n" 580b6972d74SZach Atkins " KSP convergence: %s\n" 581b6972d74SZach Atkins " Total KSP iterations: %" PetscInt_FMT "\n" 582b6972d74SZach Atkins " Final rnorm: %e\n", 583b6972d74SZach Atkins ksp_type, KSPConvergedReasons[reason], its, (double)rnorm)); 584b6972d74SZach Atkins } 585b6972d74SZach Atkins } 586b6972d74SZach Atkins 587b6972d74SZach Atkins // Optional viewing 588b6972d74SZach Atkins PetscCall(KSPViewFromOptions(ksp, NULL, "-ksp_view")); 589b6972d74SZach Atkins 590b6972d74SZach Atkins // Cleanup 591b6972d74SZach Atkins PetscCall(VecDestroy(&B_mesh)); 592b6972d74SZach Atkins PetscCall(MatDestroy(&M)); 593b6972d74SZach Atkins PetscCall(KSPDestroy(&ksp)); 594b6972d74SZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 595b6972d74SZach Atkins } 596