198285ab4SZach Atkins // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 298285ab4SZach Atkins // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 398285ab4SZach Atkins // 498285ab4SZach Atkins // SPDX-License-Identifier: BSD-2-Clause 598285ab4SZach Atkins // 698285ab4SZach Atkins // This file is part of CEED: http://github.com/ceed 798285ab4SZach Atkins 8b6972d74SZach Atkins #include "../include/swarmutils.h" 9*78f7fce3SZach Atkins #include "../include/matops.h" 10b6972d74SZach Atkins #include "../qfunctions/swarm/swarmmass.h" 11b6972d74SZach Atkins 12b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 13b6972d74SZach Atkins // Context utilities 14b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 15b6972d74SZach Atkins PetscErrorCode DMSwarmCeedContextCreate(DM dm_swarm, const char *ceed_resource, DMSwarmCeedContext *ctx) { 16b6972d74SZach Atkins DM dm_mesh, dm_coord; 17b6972d74SZach Atkins CeedElemRestriction elem_restr_u_mesh, elem_restr_x_mesh, elem_restr_x_points, elem_restr_u_points, elem_restr_q_data_points; 18b6972d74SZach Atkins CeedBasis basis_u, basis_x; 19b6972d74SZach Atkins CeedVector x_ref_points, q_data_points; 20b6972d74SZach Atkins CeedInt num_comp; 21b6972d74SZach Atkins 22b6972d74SZach Atkins PetscFunctionBeginUser; 23b6972d74SZach Atkins PetscCall(PetscNew(ctx)); 24b6972d74SZach Atkins PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh)); 25b6972d74SZach Atkins PetscCall(DMGetCoordinateDM(dm_mesh, &dm_coord)); 26b6972d74SZach Atkins 27b6972d74SZach Atkins CeedInit(ceed_resource, &(*ctx)->ceed); 28b6972d74SZach Atkins // Background mesh objects 29b6972d74SZach Atkins { 30b6972d74SZach Atkins BPData bp_data = {.q_mode = CEED_GAUSS}; 31b6972d74SZach Atkins 32b6972d74SZach Atkins PetscCall(CreateBasisFromPlex((*ctx)->ceed, dm_mesh, NULL, 0, 0, 0, bp_data, &basis_u)); 33b6972d74SZach Atkins PetscCall(CreateBasisFromPlex((*ctx)->ceed, dm_coord, NULL, 0, 0, 0, bp_data, &basis_x)); 34b6972d74SZach Atkins PetscCall(CreateRestrictionFromPlex((*ctx)->ceed, dm_mesh, 0, NULL, 0, &elem_restr_u_mesh)); 35b6972d74SZach Atkins PetscCall(CreateRestrictionFromPlex((*ctx)->ceed, dm_coord, 0, NULL, 0, &elem_restr_x_mesh)); 36b6972d74SZach Atkins 37b6972d74SZach Atkins // -- Mesh vectors 38b6972d74SZach Atkins CeedElemRestrictionCreateVector(elem_restr_u_mesh, &(*ctx)->u_mesh, NULL); 39b6972d74SZach Atkins CeedElemRestrictionCreateVector(elem_restr_u_mesh, &(*ctx)->v_mesh, NULL); 40b6972d74SZach Atkins } 41b6972d74SZach Atkins // Swarm objects 42b6972d74SZach Atkins { 43b6972d74SZach Atkins PetscInt dim; 44b6972d74SZach Atkins const PetscInt *cell_points; 45b6972d74SZach Atkins IS is_points; 46b6972d74SZach Atkins Vec X_ref; 47b6972d74SZach Atkins CeedInt num_elem; 48b6972d74SZach Atkins 49b6972d74SZach Atkins PetscCall(DMSwarmCreateReferenceCoordinates(dm_swarm, &is_points, &X_ref)); 50b6972d74SZach Atkins PetscCall(DMGetDimension(dm_mesh, &dim)); 51b6972d74SZach Atkins CeedElemRestrictionGetNumElements(elem_restr_u_mesh, &num_elem); 52b6972d74SZach Atkins CeedElemRestrictionGetNumComponents(elem_restr_u_mesh, &num_comp); 53b6972d74SZach Atkins 54b6972d74SZach Atkins PetscCall(ISGetIndices(is_points, &cell_points)); 55b6972d74SZach Atkins PetscInt num_points = cell_points[num_elem + 1] - num_elem - 2; 56b6972d74SZach Atkins CeedInt offsets[num_elem + 1 + num_points]; 57b6972d74SZach Atkins 58b6972d74SZach Atkins for (PetscInt i = 0; i < num_elem + 1; i++) offsets[i] = cell_points[i + 1] - 1; 59b6972d74SZach Atkins for (PetscInt i = num_elem + 1; i < num_points + num_elem + 1; i++) offsets[i] = cell_points[i + 1]; 60b6972d74SZach Atkins PetscCall(ISRestoreIndices(is_points, &cell_points)); 61b6972d74SZach Atkins 62b6972d74SZach Atkins // -- Points restrictions 63b6972d74SZach Atkins CeedElemRestrictionCreateAtPoints((*ctx)->ceed, num_elem, num_points, num_comp, num_points * num_comp, CEED_MEM_HOST, CEED_COPY_VALUES, offsets, 64b6972d74SZach Atkins &elem_restr_u_points); 65b6972d74SZach Atkins CeedElemRestrictionCreateAtPoints((*ctx)->ceed, num_elem, num_points, dim, num_points * dim, CEED_MEM_HOST, CEED_COPY_VALUES, offsets, 66b6972d74SZach Atkins &elem_restr_x_points); 67b6972d74SZach Atkins CeedElemRestrictionCreateAtPoints((*ctx)->ceed, num_elem, num_points, 1, num_points, CEED_MEM_HOST, CEED_COPY_VALUES, offsets, 68b6972d74SZach Atkins &elem_restr_q_data_points); 69b6972d74SZach Atkins 70b6972d74SZach Atkins // -- Points vectors 71b6972d74SZach Atkins CeedElemRestrictionCreateVector(elem_restr_u_points, &(*ctx)->u_points, NULL); 72b6972d74SZach Atkins CeedElemRestrictionCreateVector(elem_restr_q_data_points, &q_data_points, NULL); 73b6972d74SZach Atkins 74b6972d74SZach Atkins // -- Ref coordinates 75b6972d74SZach Atkins { 76b6972d74SZach Atkins PetscMemType X_mem_type; 77b6972d74SZach Atkins const PetscScalar *x; 78b6972d74SZach Atkins 79b6972d74SZach Atkins CeedVectorCreate((*ctx)->ceed, num_points * dim, &x_ref_points); 80b6972d74SZach Atkins 81b6972d74SZach Atkins PetscCall(VecGetArrayReadAndMemType(X_ref, (const PetscScalar **)&x, &X_mem_type)); 82b6972d74SZach Atkins CeedVectorSetArray(x_ref_points, MemTypeP2C(X_mem_type), CEED_COPY_VALUES, (CeedScalar *)x); 83b6972d74SZach Atkins PetscCall(VecRestoreArrayReadAndMemType(X_ref, (const PetscScalar **)&x)); 84b6972d74SZach Atkins } 85b6972d74SZach Atkins 86b6972d74SZach Atkins // Create Q data 87b6972d74SZach Atkins { 88b6972d74SZach Atkins CeedQFunction qf_setup; 89b6972d74SZach Atkins CeedOperator op_setup; 90b6972d74SZach Atkins CeedVector x_coord; 91b6972d74SZach Atkins 92b6972d74SZach Atkins { 93b6972d74SZach Atkins Vec X_loc; 94b6972d74SZach Atkins CeedInt len; 95b6972d74SZach Atkins const PetscScalar *x; 96b6972d74SZach Atkins 97b6972d74SZach Atkins PetscCall(DMGetCoordinatesLocal(dm_mesh, &X_loc)); 98b6972d74SZach Atkins PetscCall(VecGetLocalSize(X_loc, &len)); 99b6972d74SZach Atkins CeedVectorCreate((*ctx)->ceed, len, &x_coord); 100b6972d74SZach Atkins 101b6972d74SZach Atkins PetscCall(VecGetArrayRead(X_loc, &x)); 102b6972d74SZach Atkins CeedVectorSetArray(x_coord, CEED_MEM_HOST, CEED_COPY_VALUES, (CeedScalar *)x); 103b6972d74SZach Atkins PetscCall(VecRestoreArrayRead(X_loc, &x)); 104b6972d74SZach Atkins } 105b6972d74SZach Atkins 106b6972d74SZach Atkins // Setup geometric scaling 107b6972d74SZach Atkins CeedQFunctionCreateInterior((*ctx)->ceed, 1, SetupMass, SetupMass_loc, &qf_setup); 108b6972d74SZach Atkins CeedQFunctionAddInput(qf_setup, "x", dim * dim, CEED_EVAL_GRAD); 109b6972d74SZach Atkins CeedQFunctionAddInput(qf_setup, "weight", 1, CEED_EVAL_WEIGHT); 110b6972d74SZach Atkins CeedQFunctionAddOutput(qf_setup, "rho", 1, CEED_EVAL_NONE); 111b6972d74SZach Atkins 112b6972d74SZach Atkins CeedOperatorCreateAtPoints((*ctx)->ceed, qf_setup, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setup); 113b6972d74SZach Atkins CeedOperatorSetField(op_setup, "x", elem_restr_x_mesh, basis_x, CEED_VECTOR_ACTIVE); 114b6972d74SZach Atkins CeedOperatorSetField(op_setup, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE); 115b6972d74SZach Atkins CeedOperatorSetField(op_setup, "rho", elem_restr_q_data_points, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE); 116b6972d74SZach Atkins CeedOperatorAtPointsSetPoints(op_setup, elem_restr_x_points, x_ref_points); 117b6972d74SZach Atkins 118b6972d74SZach Atkins CeedOperatorApply(op_setup, x_coord, q_data_points, CEED_REQUEST_IMMEDIATE); 119b6972d74SZach Atkins 120b6972d74SZach Atkins // Cleanup 121b6972d74SZach Atkins CeedVectorDestroy(&x_coord); 122b6972d74SZach Atkins CeedQFunctionDestroy(&qf_setup); 123b6972d74SZach Atkins CeedOperatorDestroy(&op_setup); 124b6972d74SZach Atkins } 125b6972d74SZach Atkins 126b6972d74SZach Atkins // -- Cleanup 127b6972d74SZach Atkins PetscCall(ISDestroy(&is_points)); 128b6972d74SZach Atkins PetscCall(VecDestroy(&X_ref)); 129b6972d74SZach Atkins } 130b6972d74SZach Atkins 131b6972d74SZach Atkins PetscCall(DMSetApplicationContext(dm_mesh, (void *)(*ctx))); 132b6972d74SZach Atkins 133b6972d74SZach Atkins // Create operators 134b6972d74SZach Atkins // Mesh to points interpolation operator 135b6972d74SZach Atkins { 136b6972d74SZach Atkins CeedQFunction qf_mesh_to_points; 137b6972d74SZach Atkins 138b6972d74SZach Atkins // -- Create operator 139b6972d74SZach Atkins CeedQFunctionCreateIdentity((*ctx)->ceed, num_comp, CEED_EVAL_INTERP, CEED_EVAL_NONE, &qf_mesh_to_points); 140b6972d74SZach Atkins 141b6972d74SZach Atkins CeedOperatorCreateAtPoints((*ctx)->ceed, qf_mesh_to_points, NULL, NULL, &(*ctx)->op_mesh_to_points); 142b6972d74SZach Atkins CeedOperatorSetField((*ctx)->op_mesh_to_points, "input", elem_restr_u_mesh, basis_u, CEED_VECTOR_ACTIVE); 143b6972d74SZach Atkins CeedOperatorSetField((*ctx)->op_mesh_to_points, "output", elem_restr_u_points, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE); 144b6972d74SZach Atkins CeedOperatorAtPointsSetPoints((*ctx)->op_mesh_to_points, elem_restr_x_points, x_ref_points); 145b6972d74SZach Atkins 146b6972d74SZach Atkins // -- Cleanup 147b6972d74SZach Atkins CeedQFunctionDestroy(&qf_mesh_to_points); 148b6972d74SZach Atkins } 149b6972d74SZach Atkins 150b6972d74SZach Atkins // RHS operator 151b6972d74SZach Atkins { 152b6972d74SZach Atkins CeedQFunction qf_pts_to_mesh; 153b6972d74SZach Atkins CeedQFunctionContext qf_ctx; 154b6972d74SZach Atkins 155b6972d74SZach Atkins // -- Mass QFunction 156b6972d74SZach Atkins CeedQFunctionCreateInterior((*ctx)->ceed, 1, Mass, Mass_loc, &qf_pts_to_mesh); 157b6972d74SZach Atkins CeedQFunctionAddInput(qf_pts_to_mesh, "q data", 1, CEED_EVAL_NONE); 158b6972d74SZach Atkins CeedQFunctionAddInput(qf_pts_to_mesh, "u", num_comp, CEED_EVAL_NONE); 159b6972d74SZach Atkins CeedQFunctionAddOutput(qf_pts_to_mesh, "v", num_comp, CEED_EVAL_INTERP); 160b6972d74SZach Atkins 161b6972d74SZach Atkins // -- QFunction context 162b6972d74SZach Atkins CeedQFunctionContextCreate((*ctx)->ceed, &qf_ctx); 163b6972d74SZach Atkins CeedQFunctionContextSetData(qf_ctx, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof(num_comp), &num_comp); 164b6972d74SZach Atkins CeedQFunctionSetContext(qf_pts_to_mesh, qf_ctx); 165b6972d74SZach Atkins 166b6972d74SZach Atkins // -- Mass Operator 167b6972d74SZach Atkins CeedOperatorCreateAtPoints((*ctx)->ceed, qf_pts_to_mesh, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &(*ctx)->op_points_to_mesh); 168b6972d74SZach Atkins CeedOperatorSetField((*ctx)->op_points_to_mesh, "q data", elem_restr_q_data_points, CEED_BASIS_NONE, q_data_points); 169b6972d74SZach Atkins CeedOperatorSetField((*ctx)->op_points_to_mesh, "u", elem_restr_u_points, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE); 170b6972d74SZach Atkins CeedOperatorSetField((*ctx)->op_points_to_mesh, "v", elem_restr_u_mesh, basis_u, CEED_VECTOR_ACTIVE); 171b6972d74SZach Atkins CeedOperatorAtPointsSetPoints((*ctx)->op_points_to_mesh, elem_restr_x_points, x_ref_points); 172b6972d74SZach Atkins 173b6972d74SZach Atkins // -- Cleanup 174b6972d74SZach Atkins CeedQFunctionContextDestroy(&qf_ctx); 175b6972d74SZach Atkins CeedQFunctionDestroy(&qf_pts_to_mesh); 176b6972d74SZach Atkins } 177b6972d74SZach Atkins 178b6972d74SZach Atkins // Mass operator 179b6972d74SZach Atkins { 180b6972d74SZach Atkins CeedQFunction qf_mass; 181b6972d74SZach Atkins CeedQFunctionContext ctx_mass; 182b6972d74SZach Atkins 183b6972d74SZach Atkins // -- Mass QFunction 184b6972d74SZach Atkins CeedQFunctionCreateInterior((*ctx)->ceed, 1, Mass, Mass_loc, &qf_mass); 185b6972d74SZach Atkins CeedQFunctionAddInput(qf_mass, "q data", 1, CEED_EVAL_NONE); 186b6972d74SZach Atkins CeedQFunctionAddInput(qf_mass, "u", num_comp, CEED_EVAL_INTERP); 187b6972d74SZach Atkins CeedQFunctionAddOutput(qf_mass, "v", num_comp, CEED_EVAL_INTERP); 188b6972d74SZach Atkins 189b6972d74SZach Atkins // -- QFunction context 190b6972d74SZach Atkins CeedQFunctionContextCreate((*ctx)->ceed, &ctx_mass); 191b6972d74SZach Atkins CeedQFunctionContextSetData(ctx_mass, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof(num_comp), &num_comp); 192b6972d74SZach Atkins CeedQFunctionSetContext(qf_mass, ctx_mass); 193b6972d74SZach Atkins 194b6972d74SZach Atkins // -- Mass Operator 195b6972d74SZach Atkins CeedOperatorCreateAtPoints((*ctx)->ceed, qf_mass, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &(*ctx)->op_mass); 196b6972d74SZach Atkins CeedOperatorSetField((*ctx)->op_mass, "q data", elem_restr_q_data_points, CEED_BASIS_NONE, q_data_points); 197b6972d74SZach Atkins CeedOperatorSetField((*ctx)->op_mass, "u", elem_restr_u_mesh, basis_u, CEED_VECTOR_ACTIVE); 198b6972d74SZach Atkins CeedOperatorSetField((*ctx)->op_mass, "v", elem_restr_u_mesh, basis_u, CEED_VECTOR_ACTIVE); 199b6972d74SZach Atkins CeedOperatorAtPointsSetPoints((*ctx)->op_mass, elem_restr_x_points, x_ref_points); 200b6972d74SZach Atkins 201b6972d74SZach Atkins // -- Cleanup 202b6972d74SZach Atkins CeedQFunctionContextDestroy(&ctx_mass); 203b6972d74SZach Atkins CeedQFunctionDestroy(&qf_mass); 204b6972d74SZach Atkins } 205b6972d74SZach Atkins 206b6972d74SZach Atkins // Cleanup 207b6972d74SZach Atkins CeedElemRestrictionDestroy(&elem_restr_u_mesh); 208b6972d74SZach Atkins CeedElemRestrictionDestroy(&elem_restr_x_mesh); 209b6972d74SZach Atkins CeedElemRestrictionDestroy(&elem_restr_u_points); 210b6972d74SZach Atkins CeedElemRestrictionDestroy(&elem_restr_x_points); 211b6972d74SZach Atkins CeedElemRestrictionDestroy(&elem_restr_q_data_points); 212b6972d74SZach Atkins CeedBasisDestroy(&basis_u); 213b6972d74SZach Atkins CeedBasisDestroy(&basis_x); 214b6972d74SZach Atkins CeedVectorDestroy(&x_ref_points); 215b6972d74SZach Atkins CeedVectorDestroy(&q_data_points); 216b6972d74SZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 217b6972d74SZach Atkins } 218b6972d74SZach Atkins 219b6972d74SZach Atkins PetscErrorCode DMSwarmCeedContextDestroy(DMSwarmCeedContext *ctx) { 220b6972d74SZach Atkins PetscFunctionBeginUser; 221b6972d74SZach Atkins CeedDestroy(&(*ctx)->ceed); 222b6972d74SZach Atkins CeedVectorDestroy(&(*ctx)->u_mesh); 223b6972d74SZach Atkins CeedVectorDestroy(&(*ctx)->v_mesh); 224b6972d74SZach Atkins CeedVectorDestroy(&(*ctx)->u_points); 225b6972d74SZach Atkins CeedOperatorDestroy(&(*ctx)->op_mesh_to_points); 226b6972d74SZach Atkins CeedOperatorDestroy(&(*ctx)->op_points_to_mesh); 227b6972d74SZach Atkins CeedOperatorDestroy(&(*ctx)->op_mass); 228b6972d74SZach Atkins PetscCall(PetscFree(*ctx)); 229b6972d74SZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 230b6972d74SZach Atkins } 231b6972d74SZach Atkins 232b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 233b6972d74SZach Atkins // PETSc-libCEED memory space utilities 234b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 235b6972d74SZach Atkins PetscErrorCode DMSwarmPICFieldP2C(DM dm_swarm, const char *field, CeedVector x_ceed) { 236b6972d74SZach Atkins PetscScalar *x; 237b6972d74SZach Atkins 238b6972d74SZach Atkins PetscFunctionBeginUser; 239b6972d74SZach Atkins PetscCall(DMSwarmGetField(dm_swarm, field, NULL, NULL, (void **)&x)); 240b6972d74SZach Atkins CeedVectorSetArray(x_ceed, CEED_MEM_HOST, CEED_USE_POINTER, (CeedScalar *)x); 241b6972d74SZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 242b6972d74SZach Atkins } 243b6972d74SZach Atkins 244b6972d74SZach Atkins PetscErrorCode DMSwarmPICFieldC2P(DM dm_swarm, const char *field, CeedVector x_ceed) { 245b6972d74SZach Atkins PetscScalar *x; 246b6972d74SZach Atkins 247b6972d74SZach Atkins PetscFunctionBeginUser; 248b6972d74SZach Atkins CeedVectorTakeArray(x_ceed, CEED_MEM_HOST, (CeedScalar **)&x); 249b6972d74SZach Atkins PetscCall(DMSwarmRestoreField(dm_swarm, field, NULL, NULL, (void **)&x)); 250b6972d74SZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 251b6972d74SZach Atkins } 252b6972d74SZach Atkins 253b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 254b6972d74SZach Atkins // Swarm point location utility 255b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 256b6972d74SZach Atkins PetscErrorCode DMSwarmInitalizePointLocations(DM dm_swarm, PointSwarmType point_swarm_type, PetscInt num_points, PetscInt num_points_per_cell) { 257b6972d74SZach Atkins PetscFunctionBeginUser; 258b6972d74SZach Atkins switch (point_swarm_type) { 259b6972d74SZach Atkins case SWARM_GAUSS: 260b6972d74SZach Atkins case SWARM_UNIFORM: { 261b6972d74SZach Atkins // -- Set gauss or uniform point locations in each cell 262b6972d74SZach Atkins PetscInt num_points_per_cell_1d = round(cbrt(num_points_per_cell * 1.0)), dim = 3; 263b6972d74SZach Atkins PetscScalar point_coords[num_points_per_cell * 3]; 264b6972d74SZach Atkins CeedScalar points_1d[num_points_per_cell_1d], weights_1d[num_points_per_cell_1d]; 265b6972d74SZach Atkins 266b6972d74SZach Atkins if (point_swarm_type == SWARM_GAUSS) { 267b6972d74SZach Atkins PetscCall(CeedGaussQuadrature(num_points_per_cell_1d, points_1d, weights_1d)); 268b6972d74SZach Atkins } else { 269b6972d74SZach 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; 270b6972d74SZach Atkins } 271b6972d74SZach Atkins for (PetscInt i = 0; i < num_points_per_cell_1d; i++) { 272b6972d74SZach Atkins for (PetscInt j = 0; j < num_points_per_cell_1d; j++) { 273b6972d74SZach Atkins for (PetscInt k = 0; k < num_points_per_cell_1d; k++) { 274b6972d74SZach Atkins PetscInt p = (i * num_points_per_cell_1d + j) * num_points_per_cell_1d + k; 275b6972d74SZach Atkins 276b6972d74SZach Atkins point_coords[p * dim + 0] = points_1d[i]; 277b6972d74SZach Atkins point_coords[p * dim + 1] = points_1d[j]; 278b6972d74SZach Atkins point_coords[p * dim + 2] = points_1d[k]; 279b6972d74SZach Atkins } 280b6972d74SZach Atkins } 281b6972d74SZach Atkins } 282b6972d74SZach Atkins PetscCall(DMSwarmSetPointCoordinatesCellwise(dm_swarm, num_points_per_cell_1d * num_points_per_cell_1d * num_points_per_cell_1d, point_coords)); 283b6972d74SZach Atkins } break; 284b6972d74SZach Atkins case SWARM_CELL_RANDOM: { 285b6972d74SZach Atkins // -- Set points randomly in each cell 286b6972d74SZach Atkins PetscInt dim = 3, num_cells_total = 1, num_cells[] = {1, 1, 1}; 287b6972d74SZach Atkins PetscScalar *point_coords; 288b6972d74SZach Atkins PetscRandom rng; 289b6972d74SZach Atkins 290b6972d74SZach Atkins PetscOptionsBegin(PetscObjectComm((PetscObject)dm_swarm), NULL, "libCEED example using PETSc with DMSwarm", NULL); 291b6972d74SZach Atkins 292b6972d74SZach Atkins PetscCall(PetscOptionsInt("-dm_plex_dim", "Background mesh dimension", NULL, dim, &dim, NULL)); 293b6972d74SZach Atkins PetscCall(PetscOptionsIntArray("-dm_plex_box_faces", "Number of cells", NULL, num_cells, &dim, NULL)); 294b6972d74SZach Atkins 295b6972d74SZach Atkins PetscOptionsEnd(); 296b6972d74SZach Atkins 297b6972d74SZach Atkins PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm_swarm), &rng)); 298b6972d74SZach Atkins 299b6972d74SZach Atkins num_cells_total = num_cells[0] * num_cells[1] * num_cells[2]; 300b6972d74SZach Atkins PetscCall(DMSwarmGetField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&point_coords)); 301b6972d74SZach Atkins for (PetscInt c = 0; c < num_cells_total; c++) { 302b6972d74SZach 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]}; 303b6972d74SZach Atkins 304b6972d74SZach Atkins for (PetscInt p = 0; p < num_points_per_cell; p++) { 305b6972d74SZach Atkins PetscInt point_index = c * num_points_per_cell + p; 306b6972d74SZach Atkins PetscScalar random_value; 307b6972d74SZach Atkins 308b6972d74SZach Atkins for (PetscInt i = 0; i < dim; i++) { 309b6972d74SZach Atkins PetscCall(PetscRandomGetValue(rng, &random_value)); 310b6972d74SZach Atkins point_coords[point_index * dim + i] = -1.0 + cell_index[i] * 2.0 / (num_cells[i] + 1.0) + random_value; 311b6972d74SZach Atkins } 312b6972d74SZach Atkins } 313b6972d74SZach Atkins } 314b6972d74SZach Atkins PetscCall(DMSwarmRestoreField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&point_coords)); 315b6972d74SZach Atkins PetscCall(PetscRandomDestroy(&rng)); 316b6972d74SZach Atkins } break; 317b6972d74SZach Atkins case SWARM_SINUSOIDAL: { 318b6972d74SZach Atkins // -- Set points distributed per sinusoidal functions 319b6972d74SZach Atkins PetscInt dim = 3; 320b6972d74SZach Atkins PetscScalar *point_coords; 321b6972d74SZach Atkins DM dm_mesh; 322b6972d74SZach Atkins 323b6972d74SZach Atkins PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh)); 324b6972d74SZach Atkins PetscCall(DMGetDimension(dm_mesh, &dim)); 325b6972d74SZach Atkins 326b6972d74SZach Atkins PetscCall(DMSwarmGetField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&point_coords)); 327b6972d74SZach Atkins for (PetscInt p = 0; p < num_points; p++) { 328b6972d74SZach Atkins point_coords[p * dim + 0] = -PetscCosReal((PetscReal)(p + 1) / (PetscReal)(num_points + 1) * PETSC_PI); 329b6972d74SZach Atkins if (dim > 1) point_coords[p * dim + 1] = -PetscSinReal((PetscReal)(p + 1) / (PetscReal)(num_points + 1) * PETSC_PI); 330b6972d74SZach Atkins if (dim > 2) point_coords[p * dim + 2] = PetscSinReal((PetscReal)(p + 1) / (PetscReal)(num_points + 1) * PETSC_PI); 331b6972d74SZach Atkins } 332b6972d74SZach Atkins PetscCall(DMSwarmRestoreField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&point_coords)); 333b6972d74SZach Atkins } break; 334b6972d74SZach Atkins } 335b6972d74SZach Atkins PetscCall(DMSwarmMigrate(dm_swarm, PETSC_TRUE)); 336b6972d74SZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 337b6972d74SZach Atkins } 338b6972d74SZach Atkins 339b6972d74SZach Atkins /*@C 340b6972d74SZach Atkins DMSwarmCreateReferenceCoordinates - Compute the cell reference coordinates for local DMSwarm points. 341b6972d74SZach Atkins 342b6972d74SZach Atkins Collective 343b6972d74SZach Atkins 344b6972d74SZach Atkins Input Parameter: 345b6972d74SZach Atkins . dm_swarm - the `DMSwarm` 346b6972d74SZach Atkins 347b6972d74SZach Atkins Output Parameters: 348b6972d74SZach Atkins + is_points - The IS object for indexing into points per cell 349b6972d74SZach Atkins - X_points_ref - Vec holding the cell reference coordinates for local DMSwarm points 350b6972d74SZach Atkins 351b6972d74SZach Atkins The index set contains ranges of indices for each local cell. This range contains the indices of every point in the cell. 352b6972d74SZach Atkins 353b6972d74SZach Atkins ``` 354b6972d74SZach Atkins total_num_cells 355b6972d74SZach Atkins cell_0_start_index 356b6972d74SZach Atkins cell_1_start_index 357b6972d74SZach Atkins cell_2_start_index 358b6972d74SZach Atkins ... 359b6972d74SZach Atkins cell_n_start_index 360b6972d74SZach Atkins cell_n_stop_index 361b6972d74SZach Atkins cell_0_point_0 362b6972d74SZach Atkins cell_0_point_0 363b6972d74SZach Atkins ... 364b6972d74SZach Atkins cell_n_point_m 365b6972d74SZach Atkins ``` 366b6972d74SZach Atkins 367b6972d74SZach Atkins Level: beginner 368b6972d74SZach Atkins 369b6972d74SZach Atkins .seealso: `DMSwarm` 370b6972d74SZach Atkins @*/ 371b6972d74SZach Atkins PetscErrorCode DMSwarmCreateReferenceCoordinates(DM dm_swarm, IS *is_points, Vec *X_points_ref) { 372b6972d74SZach Atkins PetscInt cell_start, cell_end, num_cells_local, dim, num_points_local, *cell_points, points_offset; 373b6972d74SZach Atkins PetscScalar *coords_points_ref; 374b6972d74SZach Atkins const PetscScalar *coords_points_true; 375b6972d74SZach Atkins DM dm_mesh; 376b6972d74SZach Atkins 377b6972d74SZach Atkins PetscFunctionBeginUser; 378b6972d74SZach Atkins PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh)); 379b6972d74SZach Atkins 380b6972d74SZach Atkins // Create vector to hold reference coordinates 381b6972d74SZach Atkins { 382b6972d74SZach Atkins Vec X_points_true; 383b6972d74SZach Atkins 384b6972d74SZach Atkins PetscCall(DMSwarmCreateLocalVectorFromField(dm_swarm, DMSwarmPICField_coor, &X_points_true)); 385b6972d74SZach Atkins PetscCall(VecDuplicate(X_points_true, X_points_ref)); 386b6972d74SZach Atkins PetscCall(DMSwarmDestroyLocalVectorFromField(dm_swarm, DMSwarmPICField_coor, &X_points_true)); 387b6972d74SZach Atkins } 388b6972d74SZach Atkins 389b6972d74SZach Atkins // Allocate index set array 390b6972d74SZach Atkins PetscCall(DMPlexGetHeightStratum(dm_mesh, 0, &cell_start, &cell_end)); 391b6972d74SZach Atkins num_cells_local = cell_end - cell_start; 392b6972d74SZach Atkins points_offset = num_cells_local + 2; 393b6972d74SZach Atkins PetscCall(VecGetLocalSize(*X_points_ref, &num_points_local)); 394b6972d74SZach Atkins PetscCall(DMGetDimension(dm_mesh, &dim)); 395b6972d74SZach Atkins num_points_local /= dim; 396b6972d74SZach Atkins PetscCall(PetscMalloc1(num_points_local + num_cells_local + 2, &cell_points)); 397b6972d74SZach Atkins cell_points[0] = num_cells_local; 398b6972d74SZach Atkins 399b6972d74SZach Atkins // Get reference coordinates for each swarm point wrt the elements in the background mesh 400b6972d74SZach Atkins PetscCall(DMSwarmSortGetAccess(dm_swarm)); 401b6972d74SZach Atkins PetscCall(DMSwarmGetField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords_points_true)); 402b6972d74SZach Atkins PetscCall(VecGetArray(*X_points_ref, &coords_points_ref)); 403b6972d74SZach Atkins for (PetscInt cell = cell_start, num_points_processed = 0; cell < cell_end; cell++) { 404b6972d74SZach Atkins PetscInt *points_in_cell, num_points_in_cell, local_cell = cell - cell_start; 405b6972d74SZach Atkins PetscReal v[3], J[9], invJ[9], detJ, v0_ref[3] = {-1.0, -1.0, -1.0}; 406b6972d74SZach Atkins 407b6972d74SZach Atkins PetscCall(DMSwarmSortGetPointsPerCell(dm_swarm, cell, &num_points_in_cell, &points_in_cell)); 408b6972d74SZach Atkins // -- Reference coordinates for swarm points in background mesh element 409b6972d74SZach Atkins PetscCall(DMPlexComputeCellGeometryFEM(dm_mesh, cell, NULL, v, J, invJ, &detJ)); 410b6972d74SZach Atkins cell_points[local_cell + 1] = num_points_processed + points_offset; 411b6972d74SZach Atkins for (PetscInt p = 0; p < num_points_in_cell; p++) { 412b6972d74SZach Atkins const CeedInt point = points_in_cell[p]; 413b6972d74SZach Atkins 414b6972d74SZach Atkins cell_points[points_offset + (num_points_processed++)] = point; 415b6972d74SZach Atkins CoordinatesRealToRef(dim, dim, v0_ref, v, invJ, &coords_points_true[point * dim], &coords_points_ref[point * dim]); 416b6972d74SZach Atkins } 417b6972d74SZach Atkins 418b6972d74SZach Atkins // -- Cleanup 419b6972d74SZach Atkins PetscCall(PetscFree(points_in_cell)); 420b6972d74SZach Atkins } 421b6972d74SZach Atkins cell_points[points_offset - 1] = num_points_local + points_offset; 422b6972d74SZach Atkins 423b6972d74SZach Atkins // Cleanup 424b6972d74SZach Atkins PetscCall(DMSwarmRestoreField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords_points_true)); 425b6972d74SZach Atkins PetscCall(VecRestoreArray(*X_points_ref, &coords_points_ref)); 426b6972d74SZach Atkins PetscCall(DMSwarmSortRestoreAccess(dm_swarm)); 427b6972d74SZach Atkins 428b6972d74SZach Atkins // Create index set 429b6972d74SZach Atkins PetscCall(ISCreateGeneral(PETSC_COMM_SELF, num_points_local + points_offset, cell_points, PETSC_OWN_POINTER, is_points)); 430b6972d74SZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 431b6972d74SZach Atkins } 432b6972d74SZach Atkins 433b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 434b6972d74SZach Atkins // RHS for Swarm to Mesh projection 435b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 436*78f7fce3SZach Atkins PetscErrorCode DMSwarmCreateProjectionRHS(DM dm_swarm, const char *field, Vec U_points, Vec B_mesh) { 437*78f7fce3SZach Atkins PetscMemType B_mem_type, U_mem_type; 438b6972d74SZach Atkins DM dm_mesh; 439b6972d74SZach Atkins Vec B_mesh_loc; 440*78f7fce3SZach Atkins PetscBool has_u_points; 441b6972d74SZach Atkins DMSwarmCeedContext swarm_ceed_context; 442b6972d74SZach Atkins 443b6972d74SZach Atkins PetscFunctionBeginUser; 444b6972d74SZach Atkins // Get mesh DM 445b6972d74SZach Atkins PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh)); 446b6972d74SZach Atkins PetscCall(DMGetApplicationContext(dm_mesh, (void *)&swarm_ceed_context)); 447b6972d74SZach Atkins 448*78f7fce3SZach Atkins // Get swarm access 449*78f7fce3SZach Atkins has_u_points = !!U_points; 450*78f7fce3SZach Atkins if (!has_u_points) { 451*78f7fce3SZach Atkins PetscCall(DMSwarmSortGetAccess(dm_swarm)); 452*78f7fce3SZach Atkins PetscCall(DMSwarmCreateLocalVectorFromField(dm_swarm, field, &U_points)); 453*78f7fce3SZach Atkins } 454*78f7fce3SZach Atkins 455b6972d74SZach Atkins // Get mesh values 456*78f7fce3SZach Atkins PetscCall(VecReadP2C(U_points, &U_mem_type, swarm_ceed_context->u_points)); 457b6972d74SZach Atkins PetscCall(DMGetLocalVector(dm_mesh, &B_mesh_loc)); 458b6972d74SZach Atkins PetscCall(VecZeroEntries(B_mesh_loc)); 459b6972d74SZach Atkins PetscCall(VecP2C(B_mesh_loc, &B_mem_type, swarm_ceed_context->v_mesh)); 460b6972d74SZach Atkins 461b6972d74SZach Atkins // Interpolate field from swarm points to mesh 462b6972d74SZach Atkins CeedOperatorApply(swarm_ceed_context->op_points_to_mesh, swarm_ceed_context->u_points, swarm_ceed_context->v_mesh, CEED_REQUEST_IMMEDIATE); 463b6972d74SZach Atkins 464b6972d74SZach Atkins // Restore PETSc Vecs and Local to Global 465*78f7fce3SZach Atkins PetscCall(VecReadC2P(swarm_ceed_context->u_points, U_mem_type, U_points)); 466b6972d74SZach Atkins PetscCall(VecC2P(swarm_ceed_context->v_mesh, B_mem_type, B_mesh_loc)); 467b6972d74SZach Atkins PetscCall(VecZeroEntries(B_mesh)); 468b6972d74SZach Atkins PetscCall(DMLocalToGlobal(dm_mesh, B_mesh_loc, ADD_VALUES, B_mesh)); 469b6972d74SZach Atkins 470*78f7fce3SZach Atkins // Restore swarm access 471*78f7fce3SZach Atkins if (!has_u_points) { 472*78f7fce3SZach Atkins PetscCall(DMSwarmDestroyLocalVectorFromField(dm_swarm, field, &U_points)); 473b6972d74SZach Atkins PetscCall(DMSwarmSortRestoreAccess(dm_swarm)); 474*78f7fce3SZach Atkins } 475*78f7fce3SZach Atkins 476*78f7fce3SZach Atkins // Cleanup 477b6972d74SZach Atkins PetscCall(DMRestoreLocalVector(dm_mesh, &B_mesh_loc)); 478b6972d74SZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 479b6972d74SZach Atkins } 480b6972d74SZach Atkins 481b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 482b6972d74SZach Atkins // Swarm "mass matrix" 483b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 484b6972d74SZach Atkins PetscErrorCode MatMult_SwarmMass(Mat A, Vec U_mesh, Vec V_mesh) { 485b6972d74SZach Atkins PetscMemType U_mem_type, V_mem_type; 486b6972d74SZach Atkins DM dm_mesh; 487b6972d74SZach Atkins Vec U_mesh_loc, V_mesh_loc; 488b6972d74SZach Atkins DMSwarmCeedContext swarm_ceed_context; 489b6972d74SZach Atkins 490b6972d74SZach Atkins PetscFunctionBeginUser; 491b6972d74SZach Atkins // Get mesh DM 492b6972d74SZach Atkins PetscCall(MatGetDM(A, &dm_mesh)); 493b6972d74SZach Atkins PetscCall(DMGetApplicationContext(dm_mesh, (void *)&swarm_ceed_context)); 494b6972d74SZach Atkins 495b6972d74SZach Atkins // Global to Local and get PETSc Vec access 496b6972d74SZach Atkins PetscCall(DMGetLocalVector(dm_mesh, &U_mesh_loc)); 497b6972d74SZach Atkins PetscCall(VecZeroEntries(U_mesh_loc)); 498b6972d74SZach Atkins PetscCall(DMGlobalToLocal(dm_mesh, U_mesh, INSERT_VALUES, U_mesh_loc)); 499b6972d74SZach Atkins PetscCall(VecReadP2C(U_mesh_loc, &U_mem_type, swarm_ceed_context->u_mesh)); 500b6972d74SZach Atkins PetscCall(DMGetLocalVector(dm_mesh, &V_mesh_loc)); 501b6972d74SZach Atkins PetscCall(VecZeroEntries(V_mesh_loc)); 502b6972d74SZach Atkins PetscCall(VecP2C(V_mesh_loc, &V_mem_type, swarm_ceed_context->v_mesh)); 503b6972d74SZach Atkins 504b6972d74SZach Atkins // Apply swarm mass operator 505b6972d74SZach Atkins CeedOperatorApply(swarm_ceed_context->op_mass, swarm_ceed_context->u_mesh, swarm_ceed_context->v_mesh, CEED_REQUEST_IMMEDIATE); 506b6972d74SZach Atkins 507b6972d74SZach Atkins // Restore PETSc Vecs and Local to Global 508b6972d74SZach Atkins PetscCall(VecReadC2P(swarm_ceed_context->u_mesh, U_mem_type, U_mesh_loc)); 509b6972d74SZach Atkins PetscCall(VecC2P(swarm_ceed_context->v_mesh, V_mem_type, V_mesh_loc)); 510b6972d74SZach Atkins PetscCall(VecZeroEntries(V_mesh)); 511b6972d74SZach Atkins PetscCall(DMLocalToGlobal(dm_mesh, V_mesh_loc, ADD_VALUES, V_mesh)); 512b6972d74SZach Atkins 513b6972d74SZach Atkins // Cleanup 514b6972d74SZach Atkins PetscCall(DMRestoreLocalVector(dm_mesh, &U_mesh_loc)); 515b6972d74SZach Atkins PetscCall(DMRestoreLocalVector(dm_mesh, &V_mesh_loc)); 516b6972d74SZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 517b6972d74SZach Atkins } 518b6972d74SZach Atkins 519b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 520b6972d74SZach Atkins // Swarm to mesh projection 521b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------ 522*78f7fce3SZach Atkins PetscErrorCode DMSwarmProjectFromSwarmToCells(DM dm_swarm, const char *field, Vec U_points, Vec U_mesh) { 523b6972d74SZach Atkins PetscBool test_mode; 524b6972d74SZach Atkins Vec B_mesh; 525b6972d74SZach Atkins Mat M; 526b6972d74SZach Atkins KSP ksp; 527b6972d74SZach Atkins DM dm_mesh; 528b6972d74SZach Atkins DMSwarmCeedContext swarm_ceed_context; 529b6972d74SZach Atkins MPI_Comm comm; 530b6972d74SZach Atkins 531b6972d74SZach Atkins PetscFunctionBeginUser; 532b6972d74SZach Atkins PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swarm-to-Mesh Projection Options", NULL); 533b6972d74SZach Atkins PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, test_mode, &test_mode, NULL)); 534b6972d74SZach Atkins PetscOptionsEnd(); 535b6972d74SZach Atkins 536b6972d74SZach Atkins comm = PetscObjectComm((PetscObject)dm_swarm); 537b6972d74SZach Atkins PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh)); 538b6972d74SZach Atkins PetscCall(DMGetApplicationContext(dm_mesh, (void *)&swarm_ceed_context)); 539b6972d74SZach Atkins PetscCall(VecDuplicate(U_mesh, &B_mesh)); 540b6972d74SZach Atkins 541b6972d74SZach Atkins // Setup "mass matrix" 542b6972d74SZach Atkins { 543b6972d74SZach Atkins PetscInt l_size, g_size; 544b6972d74SZach Atkins 545b6972d74SZach Atkins PetscCall(VecGetLocalSize(U_mesh, &l_size)); 546b6972d74SZach Atkins PetscCall(VecGetSize(U_mesh, &g_size)); 547b6972d74SZach Atkins PetscCall(MatCreateShell(comm, l_size, l_size, g_size, g_size, swarm_ceed_context, &M)); 548b6972d74SZach Atkins PetscCall(MatSetDM(M, dm_mesh)); 549b6972d74SZach Atkins PetscCall(MatShellSetOperation(M, MATOP_MULT, (void (*)(void))MatMult_SwarmMass)); 550b6972d74SZach Atkins } 551b6972d74SZach Atkins 552b6972d74SZach Atkins // Setup KSP 553b6972d74SZach Atkins { 554b6972d74SZach Atkins PC pc; 555b6972d74SZach Atkins 556b6972d74SZach Atkins PetscCall(KSPCreate(comm, &ksp)); 557b6972d74SZach Atkins PetscCall(KSPGetPC(ksp, &pc)); 558b6972d74SZach Atkins PetscCall(PCSetType(pc, PCJACOBI)); 559b6972d74SZach Atkins PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 560b6972d74SZach Atkins PetscCall(KSPSetType(ksp, KSPCG)); 561b6972d74SZach Atkins PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL)); 562b6972d74SZach Atkins PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 563b6972d74SZach Atkins PetscCall(KSPSetOperators(ksp, M, M)); 564b6972d74SZach Atkins PetscCall(KSPSetFromOptions(ksp)); 565b6972d74SZach Atkins PetscCall(PetscObjectSetName((PetscObject)ksp, "Swarm-to-Mesh Projection")); 566b6972d74SZach Atkins PetscCall(KSPViewFromOptions(ksp, NULL, "-ksp_projection_view")); 567b6972d74SZach Atkins } 568b6972d74SZach Atkins 569b6972d74SZach Atkins // Setup RHS 570*78f7fce3SZach Atkins PetscCall(DMSwarmCreateProjectionRHS(dm_swarm, field, U_points, B_mesh)); 571b6972d74SZach Atkins 572b6972d74SZach Atkins // Solve 573b6972d74SZach Atkins PetscCall(VecZeroEntries(U_mesh)); 574b6972d74SZach Atkins PetscCall(KSPSolve(ksp, B_mesh, U_mesh)); 575b6972d74SZach Atkins 576b6972d74SZach Atkins // KSP summary 577b6972d74SZach Atkins { 578b6972d74SZach Atkins KSPType ksp_type; 579b6972d74SZach Atkins KSPConvergedReason reason; 580b6972d74SZach Atkins PetscReal rnorm; 581b6972d74SZach Atkins PetscInt its; 582b6972d74SZach Atkins PetscCall(KSPGetType(ksp, &ksp_type)); 583b6972d74SZach Atkins PetscCall(KSPGetConvergedReason(ksp, &reason)); 584b6972d74SZach Atkins PetscCall(KSPGetIterationNumber(ksp, &its)); 585b6972d74SZach Atkins PetscCall(KSPGetResidualNorm(ksp, &rnorm)); 586b6972d74SZach Atkins 587b6972d74SZach Atkins if (!test_mode || reason < 0 || rnorm > 1e-8) { 588b6972d74SZach Atkins PetscCall(PetscPrintf(comm, 589b6972d74SZach Atkins "Swarm-to-Mesh Projection KSP Solve:\n" 590b6972d74SZach Atkins " KSP type: %s\n" 591b6972d74SZach Atkins " KSP convergence: %s\n" 592b6972d74SZach Atkins " Total KSP iterations: %" PetscInt_FMT "\n" 593b6972d74SZach Atkins " Final rnorm: %e\n", 594b6972d74SZach Atkins ksp_type, KSPConvergedReasons[reason], its, (double)rnorm)); 595b6972d74SZach Atkins } 596b6972d74SZach Atkins } 597b6972d74SZach Atkins 598b6972d74SZach Atkins // Optional viewing 599b6972d74SZach Atkins PetscCall(KSPViewFromOptions(ksp, NULL, "-ksp_view")); 600b6972d74SZach Atkins 601b6972d74SZach Atkins // Cleanup 602b6972d74SZach Atkins PetscCall(VecDestroy(&B_mesh)); 603b6972d74SZach Atkins PetscCall(MatDestroy(&M)); 604b6972d74SZach Atkins PetscCall(KSPDestroy(&ksp)); 605b6972d74SZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 606b6972d74SZach Atkins } 607*78f7fce3SZach Atkins 608*78f7fce3SZach Atkins // ------------------------------------------------------------------------------------------------ 609*78f7fce3SZach Atkins // BP setup 610*78f7fce3SZach Atkins // ------------------------------------------------------------------------------------------------ 611*78f7fce3SZach Atkins PetscErrorCode SetupProblemSwarm(DM dm_swarm, Ceed ceed, BPData bp_data, CeedData data, PetscBool setup_rhs, Vec rhs, Vec target) { 612*78f7fce3SZach Atkins DM dm_mesh, dm_coord; 613*78f7fce3SZach Atkins CeedElemRestriction elem_restr_u_mesh, elem_restr_x_mesh, elem_restr_x_points, elem_restr_u_points, elem_restr_q_data_points; 614*78f7fce3SZach Atkins CeedBasis basis_u, basis_x; 615*78f7fce3SZach Atkins CeedVector x_coord, x_ref_points, q_data_points; 616*78f7fce3SZach Atkins CeedInt num_comp, q_data_size = bp_data.q_data_size, dim, X_loc_len; 617*78f7fce3SZach Atkins CeedScalar R = 1; // radius of the sphere 618*78f7fce3SZach Atkins CeedScalar l = 1.0 / PetscSqrtReal(3.0); // half edge of the inscribed cube 619*78f7fce3SZach Atkins Vec X_loc; 620*78f7fce3SZach Atkins PetscMemType X_mem_type; 621*78f7fce3SZach Atkins 622*78f7fce3SZach Atkins PetscFunctionBeginUser; 623*78f7fce3SZach Atkins PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh)); 624*78f7fce3SZach Atkins PetscCall(DMGetCoordinateDM(dm_mesh, &dm_coord)); 625*78f7fce3SZach Atkins 626*78f7fce3SZach Atkins // Get coordinates 627*78f7fce3SZach Atkins PetscCall(DMGetCoordinatesLocal(dm_mesh, &X_loc)); 628*78f7fce3SZach Atkins PetscCall(VecGetLocalSize(X_loc, &X_loc_len)); 629*78f7fce3SZach Atkins CeedVectorCreate(ceed, X_loc_len, &x_coord); 630*78f7fce3SZach Atkins PetscCall(VecReadP2C(X_loc, &X_mem_type, x_coord)); 631*78f7fce3SZach Atkins 632*78f7fce3SZach Atkins // Background mesh objects 633*78f7fce3SZach Atkins PetscCall(CreateBasisFromPlex(ceed, dm_mesh, NULL, 0, 0, 0, bp_data, &basis_u)); 634*78f7fce3SZach Atkins PetscCall(CreateBasisFromPlex(ceed, dm_coord, NULL, 0, 0, 0, bp_data, &basis_x)); 635*78f7fce3SZach Atkins PetscCall(CreateRestrictionFromPlex(ceed, dm_mesh, 0, NULL, 0, &elem_restr_u_mesh)); 636*78f7fce3SZach Atkins PetscCall(CreateRestrictionFromPlex(ceed, dm_coord, 0, NULL, 0, &elem_restr_x_mesh)); 637*78f7fce3SZach Atkins 638*78f7fce3SZach Atkins CeedElemRestrictionCreateVector(elem_restr_u_mesh, &data->x_ceed, NULL); 639*78f7fce3SZach Atkins CeedElemRestrictionCreateVector(elem_restr_u_mesh, &data->y_ceed, NULL); 640*78f7fce3SZach Atkins 641*78f7fce3SZach Atkins // Swarm objects 642*78f7fce3SZach Atkins { 643*78f7fce3SZach Atkins const PetscInt *cell_points; 644*78f7fce3SZach Atkins IS is_points; 645*78f7fce3SZach Atkins Vec X_ref; 646*78f7fce3SZach Atkins CeedInt num_elem; 647*78f7fce3SZach Atkins 648*78f7fce3SZach Atkins PetscCall(DMSwarmCreateReferenceCoordinates(dm_swarm, &is_points, &X_ref)); 649*78f7fce3SZach Atkins PetscCall(DMGetDimension(dm_mesh, &dim)); 650*78f7fce3SZach Atkins CeedElemRestrictionGetNumElements(elem_restr_u_mesh, &num_elem); 651*78f7fce3SZach Atkins CeedElemRestrictionGetNumComponents(elem_restr_u_mesh, &num_comp); 652*78f7fce3SZach Atkins 653*78f7fce3SZach Atkins PetscCall(ISGetIndices(is_points, &cell_points)); 654*78f7fce3SZach Atkins PetscInt num_points = cell_points[num_elem + 1] - num_elem - 2; 655*78f7fce3SZach Atkins CeedInt offsets[num_elem + 1 + num_points]; 656*78f7fce3SZach Atkins 657*78f7fce3SZach Atkins for (PetscInt i = 0; i < num_elem + 1; i++) offsets[i] = cell_points[i + 1] - 1; 658*78f7fce3SZach Atkins for (PetscInt i = num_elem + 1; i < num_points + num_elem + 1; i++) offsets[i] = cell_points[i + 1]; 659*78f7fce3SZach Atkins PetscCall(ISRestoreIndices(is_points, &cell_points)); 660*78f7fce3SZach Atkins 661*78f7fce3SZach Atkins // -- Points restrictions 662*78f7fce3SZach Atkins CeedElemRestrictionCreateAtPoints(ceed, num_elem, num_points, num_comp, num_points * num_comp, CEED_MEM_HOST, CEED_COPY_VALUES, offsets, 663*78f7fce3SZach Atkins &elem_restr_u_points); 664*78f7fce3SZach Atkins CeedElemRestrictionCreateAtPoints(ceed, num_elem, num_points, dim, num_points * dim, CEED_MEM_HOST, CEED_COPY_VALUES, offsets, 665*78f7fce3SZach Atkins &elem_restr_x_points); 666*78f7fce3SZach Atkins CeedElemRestrictionCreateAtPoints(ceed, num_elem, num_points, q_data_size, num_points * q_data_size, CEED_MEM_HOST, CEED_COPY_VALUES, offsets, 667*78f7fce3SZach Atkins &elem_restr_q_data_points); 668*78f7fce3SZach Atkins 669*78f7fce3SZach Atkins // -- Points vectors 670*78f7fce3SZach Atkins CeedElemRestrictionCreateVector(elem_restr_q_data_points, &q_data_points, NULL); 671*78f7fce3SZach Atkins 672*78f7fce3SZach Atkins // -- Ref coordinates 673*78f7fce3SZach Atkins { 674*78f7fce3SZach Atkins PetscMemType X_mem_type; 675*78f7fce3SZach Atkins const PetscScalar *x; 676*78f7fce3SZach Atkins 677*78f7fce3SZach Atkins CeedVectorCreate(ceed, num_points * dim, &x_ref_points); 678*78f7fce3SZach Atkins 679*78f7fce3SZach Atkins PetscCall(VecGetArrayReadAndMemType(X_ref, (const PetscScalar **)&x, &X_mem_type)); 680*78f7fce3SZach Atkins CeedVectorSetArray(x_ref_points, MemTypeP2C(X_mem_type), CEED_COPY_VALUES, (CeedScalar *)x); 681*78f7fce3SZach Atkins PetscCall(VecRestoreArrayReadAndMemType(X_ref, (const PetscScalar **)&x)); 682*78f7fce3SZach Atkins } 683*78f7fce3SZach Atkins 684*78f7fce3SZach Atkins // Create Q data 685*78f7fce3SZach Atkins { 686*78f7fce3SZach Atkins CeedQFunction qf_setup; 687*78f7fce3SZach Atkins CeedOperator op_setup; 688*78f7fce3SZach Atkins 689*78f7fce3SZach Atkins // Setup geometric scaling 690*78f7fce3SZach Atkins CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_geo, bp_data.setup_geo_loc, &qf_setup); 691*78f7fce3SZach Atkins CeedQFunctionAddInput(qf_setup, "x", dim, CEED_EVAL_INTERP); 692*78f7fce3SZach Atkins CeedQFunctionAddInput(qf_setup, "dx", dim * dim, CEED_EVAL_GRAD); 693*78f7fce3SZach Atkins CeedQFunctionAddInput(qf_setup, "weight", 1, CEED_EVAL_WEIGHT); 694*78f7fce3SZach Atkins CeedQFunctionAddOutput(qf_setup, "qdata", q_data_size, CEED_EVAL_NONE); 695*78f7fce3SZach Atkins 696*78f7fce3SZach Atkins CeedOperatorCreateAtPoints(ceed, qf_setup, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setup); 697*78f7fce3SZach Atkins CeedOperatorSetField(op_setup, "x", elem_restr_x_mesh, basis_x, CEED_VECTOR_ACTIVE); 698*78f7fce3SZach Atkins CeedOperatorSetField(op_setup, "dx", elem_restr_x_mesh, basis_x, CEED_VECTOR_ACTIVE); 699*78f7fce3SZach Atkins CeedOperatorSetField(op_setup, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE); 700*78f7fce3SZach Atkins CeedOperatorSetField(op_setup, "qdata", elem_restr_q_data_points, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE); 701*78f7fce3SZach Atkins CeedOperatorAtPointsSetPoints(op_setup, elem_restr_x_points, x_ref_points); 702*78f7fce3SZach Atkins 703*78f7fce3SZach Atkins CeedOperatorApply(op_setup, x_coord, q_data_points, CEED_REQUEST_IMMEDIATE); 704*78f7fce3SZach Atkins 705*78f7fce3SZach Atkins // Cleanup 706*78f7fce3SZach Atkins CeedQFunctionDestroy(&qf_setup); 707*78f7fce3SZach Atkins CeedOperatorDestroy(&op_setup); 708*78f7fce3SZach Atkins } 709*78f7fce3SZach Atkins 710*78f7fce3SZach Atkins // Cleanup 711*78f7fce3SZach Atkins PetscCall(ISDestroy(&is_points)); 712*78f7fce3SZach Atkins PetscCall(VecDestroy(&X_ref)); 713*78f7fce3SZach Atkins } 714*78f7fce3SZach Atkins 715*78f7fce3SZach Atkins // Set up PDE operator 716*78f7fce3SZach Atkins 717*78f7fce3SZach Atkins CeedQFunction qf_apply; 718*78f7fce3SZach Atkins CeedOperator op_apply; 719*78f7fce3SZach Atkins CeedInt in_scale = bp_data.in_mode == CEED_EVAL_GRAD ? dim : 1; 720*78f7fce3SZach Atkins CeedInt out_scale = bp_data.out_mode == CEED_EVAL_GRAD ? dim : 1; 721*78f7fce3SZach Atkins 722*78f7fce3SZach Atkins CeedQFunctionCreateInterior(ceed, 1, bp_data.apply, bp_data.apply_loc, &qf_apply); 723*78f7fce3SZach Atkins CeedQFunctionAddInput(qf_apply, "u", num_comp * in_scale, bp_data.in_mode); 724*78f7fce3SZach Atkins CeedQFunctionAddInput(qf_apply, "qdata", q_data_size, CEED_EVAL_NONE); 725*78f7fce3SZach Atkins CeedQFunctionAddOutput(qf_apply, "v", num_comp * out_scale, bp_data.out_mode); 726*78f7fce3SZach Atkins 727*78f7fce3SZach Atkins // Create the mass or diff operator 728*78f7fce3SZach Atkins CeedOperatorCreateAtPoints(ceed, qf_apply, NULL, NULL, &op_apply); 729*78f7fce3SZach Atkins CeedOperatorSetField(op_apply, "u", elem_restr_u_mesh, basis_u, CEED_VECTOR_ACTIVE); 730*78f7fce3SZach Atkins CeedOperatorSetField(op_apply, "qdata", elem_restr_q_data_points, CEED_BASIS_NONE, q_data_points); 731*78f7fce3SZach Atkins CeedOperatorSetField(op_apply, "v", elem_restr_u_mesh, basis_u, CEED_VECTOR_ACTIVE); 732*78f7fce3SZach Atkins CeedOperatorAtPointsSetPoints(op_apply, elem_restr_x_points, x_ref_points); 733*78f7fce3SZach Atkins 734*78f7fce3SZach Atkins // Set up RHS if needed 735*78f7fce3SZach Atkins if (setup_rhs) { 736*78f7fce3SZach Atkins CeedQFunction qf_setup_rhs; 737*78f7fce3SZach Atkins CeedOperator op_setup_rhs; 738*78f7fce3SZach Atkins Vec rhs_loc; 739*78f7fce3SZach Atkins CeedVector rhs_ceed, target_ceed; 740*78f7fce3SZach Atkins PetscMemType rhs_mem_type, target_mem_type; 741*78f7fce3SZach Atkins 742*78f7fce3SZach Atkins // Create RHS vector 743*78f7fce3SZach Atkins PetscCall(DMCreateLocalVector(dm_mesh, &rhs_loc)); 744*78f7fce3SZach Atkins 745*78f7fce3SZach Atkins CeedElemRestrictionCreateVector(elem_restr_u_points, &target_ceed, NULL); 746*78f7fce3SZach Atkins CeedElemRestrictionCreateVector(elem_restr_u_mesh, &rhs_ceed, NULL); 747*78f7fce3SZach Atkins 748*78f7fce3SZach Atkins // Create the q-function that sets up the RHS and true solution 749*78f7fce3SZach Atkins CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_rhs, bp_data.setup_rhs_loc, &qf_setup_rhs); 750*78f7fce3SZach Atkins CeedQFunctionAddInput(qf_setup_rhs, "x", dim, CEED_EVAL_INTERP); 751*78f7fce3SZach Atkins CeedQFunctionAddInput(qf_setup_rhs, "qdata", q_data_size, CEED_EVAL_NONE); 752*78f7fce3SZach Atkins CeedQFunctionAddOutput(qf_setup_rhs, "true solution", num_comp, CEED_EVAL_NONE); 753*78f7fce3SZach Atkins CeedQFunctionAddOutput(qf_setup_rhs, "rhs", num_comp, CEED_EVAL_INTERP); 754*78f7fce3SZach Atkins 755*78f7fce3SZach Atkins // Create the operator that builds the RHS and true solution 756*78f7fce3SZach Atkins CeedOperatorCreateAtPoints(ceed, qf_setup_rhs, NULL, NULL, &op_setup_rhs); 757*78f7fce3SZach Atkins CeedOperatorSetField(op_setup_rhs, "x", elem_restr_x_mesh, basis_x, CEED_VECTOR_ACTIVE); 758*78f7fce3SZach Atkins CeedOperatorSetField(op_setup_rhs, "qdata", elem_restr_q_data_points, CEED_BASIS_NONE, q_data_points); 759*78f7fce3SZach Atkins CeedOperatorSetField(op_setup_rhs, "true solution", elem_restr_u_points, CEED_BASIS_NONE, target_ceed); 760*78f7fce3SZach Atkins CeedOperatorSetField(op_setup_rhs, "rhs", elem_restr_u_mesh, basis_u, CEED_VECTOR_ACTIVE); 761*78f7fce3SZach Atkins CeedOperatorAtPointsSetPoints(op_setup_rhs, elem_restr_x_points, x_ref_points); 762*78f7fce3SZach Atkins 763*78f7fce3SZach Atkins // Set up the libCEED context 764*78f7fce3SZach Atkins CeedQFunctionContext ctx_rhs_setup; 765*78f7fce3SZach Atkins CeedQFunctionContextCreate(ceed, &ctx_rhs_setup); 766*78f7fce3SZach Atkins CeedScalar rhs_setup_data[2] = {R, l}; 767*78f7fce3SZach Atkins CeedQFunctionContextSetData(ctx_rhs_setup, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof rhs_setup_data, &rhs_setup_data); 768*78f7fce3SZach Atkins CeedQFunctionSetContext(qf_setup_rhs, ctx_rhs_setup); 769*78f7fce3SZach Atkins CeedQFunctionContextDestroy(&ctx_rhs_setup); 770*78f7fce3SZach Atkins 771*78f7fce3SZach Atkins // Setup RHS and target 772*78f7fce3SZach Atkins PetscCall(VecP2C(rhs_loc, &rhs_mem_type, rhs_ceed)); 773*78f7fce3SZach Atkins PetscCall(VecP2C(target, &target_mem_type, target_ceed)); 774*78f7fce3SZach Atkins CeedOperatorApply(op_setup_rhs, x_coord, rhs_ceed, CEED_REQUEST_IMMEDIATE); 775*78f7fce3SZach Atkins PetscCall(VecC2P(rhs_ceed, rhs_mem_type, rhs_loc)); 776*78f7fce3SZach Atkins PetscCall(VecC2P(target_ceed, target_mem_type, target)); 777*78f7fce3SZach Atkins 778*78f7fce3SZach Atkins // Local-to-global 779*78f7fce3SZach Atkins PetscCall(VecZeroEntries(rhs)); 780*78f7fce3SZach Atkins PetscCall(DMLocalToGlobal(dm_mesh, rhs_loc, ADD_VALUES, rhs)); 781*78f7fce3SZach Atkins 782*78f7fce3SZach Atkins PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view")); 783*78f7fce3SZach Atkins 784*78f7fce3SZach Atkins // Cleanup 785*78f7fce3SZach Atkins PetscCall(DMRestoreLocalVector(dm_mesh, &rhs_loc)); 786*78f7fce3SZach Atkins CeedVectorDestroy(&rhs_ceed); 787*78f7fce3SZach Atkins CeedVectorDestroy(&target_ceed); 788*78f7fce3SZach Atkins CeedQFunctionDestroy(&qf_setup_rhs); 789*78f7fce3SZach Atkins CeedOperatorDestroy(&op_setup_rhs); 790*78f7fce3SZach Atkins } 791*78f7fce3SZach Atkins 792*78f7fce3SZach Atkins // Save libCEED data 793*78f7fce3SZach Atkins data->basis_x = basis_x; 794*78f7fce3SZach Atkins data->basis_u = basis_u; 795*78f7fce3SZach Atkins data->elem_restr_x = elem_restr_x_mesh; 796*78f7fce3SZach Atkins data->elem_restr_u = elem_restr_u_mesh; 797*78f7fce3SZach Atkins data->elem_restr_u_i = elem_restr_u_points; 798*78f7fce3SZach Atkins data->elem_restr_qd_i = elem_restr_q_data_points; 799*78f7fce3SZach Atkins data->qf_apply = qf_apply; 800*78f7fce3SZach Atkins data->op_apply = op_apply; 801*78f7fce3SZach Atkins data->q_data = q_data_points; 802*78f7fce3SZach Atkins data->q_data_size = q_data_size; 803*78f7fce3SZach Atkins 804*78f7fce3SZach Atkins // Cleanup 805*78f7fce3SZach Atkins PetscCall(VecReadC2P(x_coord, X_mem_type, X_loc)); 806*78f7fce3SZach Atkins CeedVectorDestroy(&x_coord); 807*78f7fce3SZach Atkins CeedVectorDestroy(&x_ref_points); 808*78f7fce3SZach Atkins CeedElemRestrictionDestroy(&elem_restr_x_points); 809*78f7fce3SZach Atkins 810*78f7fce3SZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 811*78f7fce3SZach Atkins } 812