1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3 /// @file 4 /// Functions for setting up and projecting the divergence of the diffusive flux 5 6 #include "../qfunctions/diff_flux_projection.h" 7 8 #include <petscdmplex.h> 9 10 #include <navierstokes.h> 11 12 const char *const DivDiffFluxProjectionMethods[] = {"NONE", "DIRECT", "INDIRECT", "DivDiffFluxProjectionMethod", "DIV_DIFF_FLUX_PROJ_", NULL}; 13 14 /** 15 @brief Create `DivDiffFluxProjectionData` for solution DM in `honee` 16 17 @param[in] honee `Honee` context 18 @param[in] divFdiffproj_method Method used to perform the divergence of diffusive flux method (can be `DIV_DIFF_FLUX_PROJ_NONE`) 19 @param[in] num_diff_flux_comps Number of components that makes up the diffusive flux (e.g. 1 for scalar advection-diffusion) 20 @param[out] diff_flux_proj The `DivDiffFluxProjectionData` object created, or set to `NULL` if `divFdiffproj_method = DIV_DIFF_FLUX_PROJ_NONE` 21 **/ 22 PetscErrorCode DivDiffFluxProjectionCreate(Honee honee, DivDiffFluxProjectionMethod divFdiffproj_method, PetscInt num_diff_flux_comps, 23 DivDiffFluxProjectionData *diff_flux_proj) { 24 PetscInt height = 0, dm_field = 0, dim, degree = honee->app_ctx->degree, q_extra = honee->app_ctx->q_extra; 25 DivDiffFluxProjectionData diff_flux_proj_; 26 NodalProjectionData projection; 27 28 PetscFunctionBeginUser; 29 if (divFdiffproj_method == DIV_DIFF_FLUX_PROJ_NONE) { 30 *diff_flux_proj = NULL; 31 PetscFunctionReturn(PETSC_SUCCESS); 32 } 33 PetscCall(PetscNew(&projection)); 34 PetscCall(PetscNew(&diff_flux_proj_)); 35 *diff_flux_proj_ = (struct DivDiffFluxProjectionData_){ 36 .method = divFdiffproj_method, 37 .num_diff_flux_comps = num_diff_flux_comps, 38 .projection = projection, 39 }; 40 41 PetscCall(DMClone(honee->dm, &projection->dm)); 42 PetscCall(DMSetMatrixPreallocateSkip(projection->dm, PETSC_TRUE)); 43 PetscCall(DMGetDimension(projection->dm, &dim)); 44 switch (diff_flux_proj_->method) { 45 case DIV_DIFF_FLUX_PROJ_DIRECT: { 46 projection->num_comp = diff_flux_proj_->num_diff_flux_comps; 47 PetscCall(PetscObjectSetName((PetscObject)projection->dm, "DivDiffFluxProj")); 48 PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &projection->num_comp, projection->dm)); 49 50 PetscCall(DMPlexCeedElemRestrictionCreate(honee->ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, 51 &diff_flux_proj_->elem_restr_div_diff_flux)); 52 PetscCallCeed(honee->ceed, 53 CeedElemRestrictionCreateVector(diff_flux_proj_->elem_restr_div_diff_flux, &diff_flux_proj_->div_diff_flux_ceed, NULL)); 54 PetscCall(DMPlexCeedBasisCreate(honee->ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, 55 &diff_flux_proj_->basis_div_diff_flux)); 56 diff_flux_proj_->eval_mode_div_diff_flux = CEED_EVAL_INTERP; 57 58 { // Create face labels on projection->dm for boundary integrals 59 DMLabel face_sets_label; 60 PetscInt num_face_set_values, *face_set_values; 61 62 PetscCall(DMGetLabel(honee->dm, "Face Sets", &face_sets_label)); 63 PetscCall(DMLabelCreateGlobalValueArray(honee->dm, face_sets_label, &num_face_set_values, &face_set_values)); 64 for (PetscInt f = 0; f < num_face_set_values; f++) { 65 DMLabel face_orientation_label; 66 char *face_orientation_label_name; 67 68 PetscCall(DMPlexCreateFaceLabel(honee->dm, face_set_values[f], &face_orientation_label_name)); 69 PetscCall(DMGetLabel(honee->dm, face_orientation_label_name, &face_orientation_label)); 70 PetscCall(DMAddLabel(projection->dm, face_orientation_label)); 71 PetscCall(PetscFree(face_orientation_label_name)); 72 } 73 PetscCall(PetscFree(face_set_values)); 74 } 75 } break; 76 case DIV_DIFF_FLUX_PROJ_INDIRECT: { 77 projection->num_comp = diff_flux_proj_->num_diff_flux_comps * dim; 78 PetscCall(PetscObjectSetName((PetscObject)projection->dm, "DiffFluxProj")); 79 PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &projection->num_comp, projection->dm)); 80 81 PetscCall(DMPlexCeedElemRestrictionQDataCreate(honee->ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, 82 diff_flux_proj_->num_diff_flux_comps, &diff_flux_proj_->elem_restr_div_diff_flux)); 83 PetscCallCeed(honee->ceed, 84 CeedElemRestrictionCreateVector(diff_flux_proj_->elem_restr_div_diff_flux, &diff_flux_proj_->div_diff_flux_ceed, NULL)); 85 diff_flux_proj_->basis_div_diff_flux = CEED_BASIS_NONE; 86 diff_flux_proj_->eval_mode_div_diff_flux = CEED_EVAL_NONE; 87 } break; 88 case DIV_DIFF_FLUX_PROJ_NONE: 89 SETERRQ(honee->comm, PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s", 90 DivDiffFluxProjectionMethods[divFdiffproj_method]); 91 break; 92 } 93 *diff_flux_proj = diff_flux_proj_; 94 PetscFunctionReturn(PETSC_SUCCESS); 95 }; 96 97 /** 98 @brief Return the objects required for the Divergence of Diffusive flux to be read by a `CeedOperator` 99 100 @param[in] diff_flux_proj Projection object 101 @param[out] elem_restr Element restriction for the divergence of diffusive flux, or `NULL` 102 @param[out] basis Basis for the divergence of diffusive flux, or `NULL` 103 @param[out] vector Vector for the divergence of diffusive flux, or `NULL` 104 @param[out] eval_mode Eval mode for the divergence of diffusive flux, or `NULL` 105 **/ 106 PetscErrorCode DivDiffFluxProjectionGetOperatorFieldData(DivDiffFluxProjectionData diff_flux_proj, CeedElemRestriction *elem_restr, CeedBasis *basis, 107 CeedVector *vector, CeedEvalMode *eval_mode) { 108 Ceed ceed = CeedVectorReturnCeed(diff_flux_proj->div_diff_flux_ceed); 109 110 PetscFunctionBeginUser; 111 if (elem_restr) PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(diff_flux_proj->elem_restr_div_diff_flux, elem_restr)); 112 if (basis) PetscCallCeed(ceed, CeedBasisReferenceCopy(diff_flux_proj->basis_div_diff_flux, basis)); 113 if (vector) PetscCallCeed(ceed, CeedVectorReferenceCopy(diff_flux_proj->div_diff_flux_ceed, vector)); 114 if (eval_mode) *eval_mode = diff_flux_proj->eval_mode_div_diff_flux; 115 PetscFunctionReturn(PETSC_SUCCESS); 116 } 117 118 /** 119 @brief Setup direct projection of divergence of diffusive flux 120 121 @param[in] honee `Honee` context 122 @param[in,out] diff_flux_proj Flux projection object to setup 123 **/ 124 static PetscErrorCode DivDiffFluxProjectionSetup_Direct(Honee honee, DivDiffFluxProjectionData diff_flux_proj) { 125 Ceed ceed = honee->ceed; 126 NodalProjectionData projection = diff_flux_proj->projection; 127 MPI_Comm comm = PetscObjectComm((PetscObject)projection->dm); 128 129 PetscFunctionBeginUser; 130 { // Create Projection RHS OperatorApplyContext 131 CeedOperator op_rhs; 132 133 PetscCheck(diff_flux_proj->CreateRHSOperator_Direct, comm, PETSC_ERR_ARG_WRONGSTATE, 134 "Must define CreateRHSOperator_Direct to use indirect div_diff_flux projection"); 135 PetscCall(diff_flux_proj->CreateRHSOperator_Direct(honee, diff_flux_proj, &op_rhs)); 136 PetscCall(DMCreateLocalVector(projection->dm, &diff_flux_proj->DivDiffFlux_loc)); 137 diff_flux_proj->ceed_vec_has_array = PETSC_FALSE; 138 PetscCall(OperatorApplyContextCreate(honee->dm, projection->dm, ceed, op_rhs, NULL, NULL, NULL, diff_flux_proj->DivDiffFlux_loc, 139 &projection->l2_rhs_ctx)); 140 PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs)); 141 } 142 143 { // -- Build Mass operator 144 CeedQFunction qf_mass; 145 CeedOperator op_mass; 146 CeedBasis basis_div_diff_flux = NULL; 147 CeedElemRestriction elem_restr_div_diff_flux_volume = NULL, elem_restr_qd; 148 CeedVector q_data; 149 CeedInt q_data_size; 150 151 PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_div_diff_flux_volume, &basis_div_diff_flux, NULL, NULL)); 152 PetscCall(QDataGet(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, honee->elem_restr_x, honee->basis_x, honee->x_coord, 153 &elem_restr_qd, &q_data, &q_data_size)); 154 155 PetscCall(HoneeMassQFunctionCreate(ceed, projection->num_comp, q_data_size, &qf_mass)); 156 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass)); 157 PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_div_diff_flux_volume, basis_div_diff_flux, CEED_VECTOR_ACTIVE)); 158 PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 159 PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_div_diff_flux_volume, basis_div_diff_flux, CEED_VECTOR_ACTIVE)); 160 161 { // -- Setup KSP for L^2 projection 162 Mat mat_mass; 163 164 PetscCall(MatCreateCeed(projection->dm, projection->dm, op_mass, NULL, &mat_mass)); 165 166 PetscCall(KSPCreate(comm, &projection->ksp)); 167 PetscCall(KSPSetOptionsPrefix(projection->ksp, "div_diff_flux_projection_")); 168 { // lumped by default 169 PC pc; 170 PetscCall(KSPGetPC(projection->ksp, &pc)); 171 PetscCall(PCSetType(pc, PCJACOBI)); 172 PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 173 PetscCall(KSPSetType(projection->ksp, KSPPREONLY)); 174 } 175 PetscCall(KSPSetFromOptions_WithMatCeed(projection->ksp, mat_mass)); 176 PetscCall(MatDestroy(&mat_mass)); 177 } 178 179 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_div_diff_flux_volume)); 180 PetscCallCeed(ceed, CeedBasisDestroy(&basis_div_diff_flux)); 181 PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 182 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 183 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 184 PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 185 } 186 PetscFunctionReturn(PETSC_SUCCESS); 187 } 188 189 /** 190 @brief Setup indirect projection of divergence of diffusive flux 191 192 @param[in] honee `Honee` context 193 @param[in,out] diff_flux_proj Flux projection object to setup 194 **/ 195 static PetscErrorCode DivDiffFluxProjectionSetup_Indirect(Honee honee, DivDiffFluxProjectionData diff_flux_proj) { 196 Ceed ceed = honee->ceed; 197 NodalProjectionData projection = diff_flux_proj->projection; 198 CeedBasis basis_diff_flux; 199 CeedElemRestriction elem_restr_diff_flux, elem_restr_qd; 200 CeedVector q_data; 201 CeedInt q_data_size; 202 MPI_Comm comm = PetscObjectComm((PetscObject)projection->dm); 203 204 PetscFunctionBeginUser; 205 { 206 PetscInt height = 0, dm_field = 0; 207 208 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &elem_restr_diff_flux)); 209 PetscCall(DMPlexCeedBasisCreate(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &basis_diff_flux)); 210 PetscCall(QDataGet(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, honee->elem_restr_x, honee->basis_x, honee->x_coord, 211 &elem_restr_qd, &q_data, &q_data_size)); 212 } 213 214 { 215 CeedOperator op_rhs; 216 217 PetscCheck(diff_flux_proj->CreateRHSOperator_Indirect, comm, PETSC_ERR_ARG_WRONGSTATE, 218 "Must define CreateRHSOperator_Indirect to use indirect div_diff_flux projection"); 219 PetscCall(diff_flux_proj->CreateRHSOperator_Indirect(honee, diff_flux_proj, &op_rhs)); 220 PetscCall(OperatorApplyContextCreate(honee->dm, projection->dm, ceed, op_rhs, NULL, NULL, NULL, NULL, &projection->l2_rhs_ctx)); 221 PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs)); 222 } 223 224 { // -- Build Mass operator 225 CeedQFunction qf_mass; 226 CeedOperator op_mass; 227 228 PetscCall(HoneeMassQFunctionCreate(ceed, projection->num_comp, q_data_size, &qf_mass)); 229 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass)); 230 PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE)); 231 PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 232 PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE)); 233 234 { // -- Setup KSP for L^2 projection 235 Mat mat_mass; 236 237 PetscCall(MatCreateCeed(projection->dm, projection->dm, op_mass, NULL, &mat_mass)); 238 239 PetscCall(KSPCreate(comm, &projection->ksp)); 240 PetscCall(KSPSetOptionsPrefix(projection->ksp, "div_diff_flux_projection_")); 241 { // lumped by default 242 PC pc; 243 PetscCall(KSPGetPC(projection->ksp, &pc)); 244 PetscCall(PCSetType(pc, PCJACOBI)); 245 PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 246 PetscCall(KSPSetType(projection->ksp, KSPPREONLY)); 247 } 248 PetscCall(KSPSetFromOptions_WithMatCeed(projection->ksp, mat_mass)); 249 PetscCall(MatDestroy(&mat_mass)); 250 } 251 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 252 PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 253 } 254 255 { // Create OperatorApplyContext to calculate divergence at quadrature points 256 CeedQFunction qf_calc_divergence = NULL; 257 CeedOperator op_calc_divergence; 258 CeedElemRestriction elem_restr_div_diff_flux = NULL; 259 PetscInt dim; 260 261 PetscCall(DMGetDimension(projection->dm, &dim)); 262 PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_div_diff_flux, NULL, NULL, NULL)); 263 264 switch (dim) { 265 case 2: 266 switch (diff_flux_proj->num_diff_flux_comps) { 267 case 1: 268 PetscCallCeed(ceed, 269 CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux2D_1, ComputeDivDiffusiveFlux2D_1_loc, &qf_calc_divergence)); 270 break; 271 } 272 break; 273 case 3: 274 switch (diff_flux_proj->num_diff_flux_comps) { 275 case 1: 276 PetscCallCeed(ceed, 277 CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux3D_1, ComputeDivDiffusiveFlux3D_1_loc, &qf_calc_divergence)); 278 break; 279 case 4: 280 PetscCallCeed(ceed, 281 CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux3D_4, ComputeDivDiffusiveFlux3D_4_loc, &qf_calc_divergence)); 282 break; 283 } 284 break; 285 } 286 PetscCheck(qf_calc_divergence, comm, PETSC_ERR_SUP, 287 "QFunction for calculating divergence of diffusive flux does not exist for" 288 " %" PetscInt_FMT " dimensional grid and %" PetscInt_FMT 289 " number of components.\nA new qfunction can be easily added; see source code for pattern.", 290 dim, diff_flux_proj->num_diff_flux_comps); 291 292 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_calc_divergence, "Grad F_diff", projection->num_comp * dim, CEED_EVAL_GRAD)); 293 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_calc_divergence, "qdata", q_data_size, CEED_EVAL_NONE)); 294 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_calc_divergence, "Div F_diff", diff_flux_proj->num_diff_flux_comps, CEED_EVAL_NONE)); 295 296 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_calc_divergence, NULL, NULL, &op_calc_divergence)); 297 PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "Grad F_diff", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE)); 298 PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 299 PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "Div F_diff", elem_restr_div_diff_flux, CEED_BASIS_NONE, 300 diff_flux_proj->div_diff_flux_ceed)); 301 302 PetscCall(OperatorApplyContextCreate(projection->dm, NULL, ceed, op_calc_divergence, NULL, CEED_VECTOR_NONE, NULL, NULL, 303 &diff_flux_proj->calc_div_diff_flux)); 304 305 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_div_diff_flux)); 306 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_calc_divergence)); 307 PetscCallCeed(ceed, CeedOperatorDestroy(&op_calc_divergence)); 308 } 309 PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 310 PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 311 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 312 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux)); 313 PetscFunctionReturn(PETSC_SUCCESS); 314 } 315 316 /** 317 @brief Setup projection of divergence of diffusive flux 318 319 @param[in] honee `Honee` context 320 @param[in,out] diff_flux_proj Flux projection object to setup 321 **/ 322 PetscErrorCode DivDiffFluxProjectionSetup(Honee honee, DivDiffFluxProjectionData diff_flux_proj) { 323 PetscFunctionBeginUser; 324 switch (honee->app_ctx->divFdiffproj_method) { 325 case DIV_DIFF_FLUX_PROJ_DIRECT: 326 PetscCall(DivDiffFluxProjectionSetup_Direct(honee, diff_flux_proj)); 327 break; 328 case DIV_DIFF_FLUX_PROJ_INDIRECT: 329 PetscCall(DivDiffFluxProjectionSetup_Indirect(honee, diff_flux_proj)); 330 break; 331 case DIV_DIFF_FLUX_PROJ_NONE: 332 SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s", 333 DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]); 334 break; 335 } 336 PetscFunctionReturn(PETSC_SUCCESS); 337 } 338 339 /** 340 @brief Project the divergence of diffusive flux 341 342 This implicitly sets the `CeedVector` input (`div_diff_flux_ceed`) to the divergence of diffusive flux. 343 344 @param[in] diff_flux_proj `NodalProjectionData` for the projection 345 @param[in] Q_loc Localized solution vector 346 **/ 347 PetscErrorCode DivDiffFluxProjectionApply(DivDiffFluxProjectionData diff_flux_proj, Vec Q_loc) { 348 NodalProjectionData projection = diff_flux_proj->projection; 349 350 PetscFunctionBeginUser; 351 PetscCall(PetscLogEventBegin(HONEE_DivDiffFluxProjection, Q_loc, 0, 0, 0)); 352 switch (diff_flux_proj->method) { 353 case DIV_DIFF_FLUX_PROJ_DIRECT: { 354 Vec DivDiffFlux, RHS; 355 356 PetscCall(DMGetGlobalVector(projection->dm, &DivDiffFlux)); 357 PetscCall(DMGetGlobalVector(projection->dm, &RHS)); 358 if (diff_flux_proj->ceed_vec_has_array) { 359 PetscCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->DivDiffFlux_loc)); 360 diff_flux_proj->ceed_vec_has_array = PETSC_FALSE; 361 } 362 PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, RHS, projection->l2_rhs_ctx)); 363 PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_rhs_view")); 364 365 { 366 // Run PCApply manually if using ksp_type preonly -pc_type jacobi 367 // This is to avoid an AllReduce call in KSPSolve_Preonly, which causes significant slowdowns for lumped mass matrix solves. 368 // See https://gitlab.com/petsc/petsc/-/merge_requests/8048 for more details and a possible fix 369 PC pc; 370 PetscBool ispreonly, isjacobi; 371 PetscCall(KSPGetPC(projection->ksp, &pc)); 372 PetscCall(PetscObjectTypeCompare((PetscObject)projection->ksp, KSPPREONLY, &ispreonly)); 373 PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &isjacobi)); 374 if (ispreonly && isjacobi) PetscCall(PCApply(pc, RHS, DivDiffFlux)); 375 else PetscCall(KSPSolve(projection->ksp, RHS, DivDiffFlux)); 376 } 377 PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_view")); 378 379 PetscCall(DMGlobalToLocal(projection->dm, DivDiffFlux, INSERT_VALUES, diff_flux_proj->DivDiffFlux_loc)); 380 PetscCall(VecReadPetscToCeed(diff_flux_proj->DivDiffFlux_loc, &diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->div_diff_flux_ceed)); 381 diff_flux_proj->ceed_vec_has_array = PETSC_TRUE; 382 383 PetscCall(DMRestoreGlobalVector(projection->dm, &RHS)); 384 PetscCall(DMRestoreGlobalVector(projection->dm, &DivDiffFlux)); 385 break; 386 } 387 case DIV_DIFF_FLUX_PROJ_INDIRECT: { 388 Vec DiffFlux, RHS; 389 390 PetscCall(DMGetGlobalVector(projection->dm, &DiffFlux)); 391 PetscCall(DMGetGlobalVector(projection->dm, &RHS)); 392 PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, RHS, projection->l2_rhs_ctx)); 393 PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_rhs_view")); 394 395 { 396 // Run PCApply manually if using -ksp_type preonly -pc_type jacobi 397 // This is to avoid an AllReduce call in KSPSolve_Preonly, which causes significant slowdowns for lumped mass matrix solves. 398 // See https://gitlab.com/petsc/petsc/-/merge_requests/8048 for more details and a possible fix 399 PC pc; 400 PetscBool ispreonly, isjacobi; 401 PetscCall(KSPGetPC(projection->ksp, &pc)); 402 PetscCall(PetscObjectTypeCompare((PetscObject)projection->ksp, KSPPREONLY, &ispreonly)); 403 PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &isjacobi)); 404 if (ispreonly && isjacobi) PetscCall(PCApply(pc, RHS, DiffFlux)); 405 else PetscCall(KSPSolve(projection->ksp, RHS, DiffFlux)); 406 } 407 PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_view")); 408 409 PetscCall(ApplyCeedOperatorGlobalToLocal(DiffFlux, NULL, diff_flux_proj->calc_div_diff_flux)); 410 PetscCall(DMRestoreGlobalVector(projection->dm, &RHS)); 411 PetscCall(DMRestoreGlobalVector(projection->dm, &DiffFlux)); 412 } break; 413 case DIV_DIFF_FLUX_PROJ_NONE: 414 SETERRQ(PetscObjectComm((PetscObject)projection->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s", 415 DivDiffFluxProjectionMethods[diff_flux_proj->method]); 416 break; 417 } 418 PetscCall(PetscLogEventEnd(HONEE_DivDiffFluxProjection, Q_loc, 0, 0, 0)); 419 PetscFunctionReturn(PETSC_SUCCESS); 420 } 421 422 /** 423 @brief Destroy `DivDiffFluxProjectionData` object 424 425 @param[in,out] diff_flux_proj Object to destroy 426 **/ 427 PetscErrorCode DivDiffFluxProjectionDataDestroy(DivDiffFluxProjectionData diff_flux_proj) { 428 PetscFunctionBeginUser; 429 if (diff_flux_proj == NULL) PetscFunctionReturn(PETSC_SUCCESS); 430 Ceed ceed = CeedVectorReturnCeed(diff_flux_proj->div_diff_flux_ceed); 431 432 PetscCall(NodalProjectionDataDestroy(&diff_flux_proj->projection)); 433 PetscCall(OperatorApplyContextDestroy(diff_flux_proj->calc_div_diff_flux)); 434 if (diff_flux_proj->ceed_vec_has_array) { 435 PetscCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->DivDiffFlux_loc)); 436 diff_flux_proj->ceed_vec_has_array = PETSC_FALSE; 437 } 438 PetscCallCeed(ceed, CeedVectorDestroy(&diff_flux_proj->div_diff_flux_ceed)); 439 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&diff_flux_proj->elem_restr_div_diff_flux)); 440 PetscCallCeed(ceed, CeedBasisDestroy(&diff_flux_proj->basis_div_diff_flux)); 441 PetscCall(VecDestroy(&diff_flux_proj->DivDiffFlux_loc)); 442 PetscCall(PetscFree(diff_flux_proj)); 443 PetscFunctionReturn(PETSC_SUCCESS); 444 } 445