18c85b835SJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 28c85b835SJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 38c85b835SJames Wright /// @file 48c85b835SJames Wright /// Functions for setting up and projecting the divergence of the diffusive flux 58c85b835SJames Wright 68c85b835SJames Wright #include "../qfunctions/diff_flux_projection.h" 78c85b835SJames Wright 88c85b835SJames Wright #include <petscdmplex.h> 98c85b835SJames Wright 108c85b835SJames Wright #include <navierstokes.h> 118c85b835SJames Wright 12e5a8cae0SJames Wright const char *const DivDiffFluxProjectionMethods[] = {"NONE", "DIRECT", "INDIRECT", "DivDiffFluxProjectionMethod", "DIV_DIFF_FLUX_PROJ_", NULL}; 13e5a8cae0SJames Wright 148c85b835SJames Wright /** 150c373b74SJames Wright @brief Create `DivDiffFluxProjectionData` for solution DM in `honee` 168c85b835SJames Wright 170c373b74SJames Wright @param[in] honee `Honee` context 18d9e69621SJames Wright @param[in] divFdiffproj_method Method used to perform the divergence of diffusive flux method (can be `DIV_DIFF_FLUX_PROJ_NONE`) 1936038bbcSJames Wright @param[in] num_diff_flux_comps Number of components that makes up the diffusive flux (e.g. 1 for scalar advection-diffusion) 20d9e69621SJames Wright @param[out] diff_flux_proj The `DivDiffFluxProjectionData` object created, or set to `NULL` if `divFdiffproj_method = DIV_DIFF_FLUX_PROJ_NONE` 218c85b835SJames Wright **/ 22d9e69621SJames Wright PetscErrorCode DivDiffFluxProjectionCreate(Honee honee, DivDiffFluxProjectionMethod divFdiffproj_method, PetscInt num_diff_flux_comps, 23d9e69621SJames Wright DivDiffFluxProjectionData *diff_flux_proj) { 240c373b74SJames Wright PetscInt label_value = 0, height = 0, dm_field = 0, dim, degree = honee->app_ctx->degree, q_extra = honee->app_ctx->q_extra; 258c85b835SJames Wright DMLabel domain_label = NULL; 26d9e69621SJames Wright DivDiffFluxProjectionData diff_flux_proj_; 278c85b835SJames Wright NodalProjectionData projection; 288c85b835SJames Wright 298c85b835SJames Wright PetscFunctionBeginUser; 30d9e69621SJames Wright if (divFdiffproj_method == DIV_DIFF_FLUX_PROJ_NONE) { 31d9e69621SJames Wright *diff_flux_proj = NULL; 3236038bbcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3336038bbcSJames Wright } 34f5dc303cSJames Wright PetscCall(PetscNew(&projection)); 35d9e69621SJames Wright PetscCall(PetscNew(&diff_flux_proj_)); 36f5dc303cSJames Wright *diff_flux_proj_ = (struct DivDiffFluxProjectionData_){ 37f5dc303cSJames Wright .method = divFdiffproj_method, 38f5dc303cSJames Wright .num_diff_flux_comps = num_diff_flux_comps, 39f5dc303cSJames Wright .projection = projection, 40f5dc303cSJames Wright }; 418c85b835SJames Wright 420c373b74SJames Wright PetscCall(DMClone(honee->dm, &projection->dm)); 438c85b835SJames Wright PetscCall(DMSetMatrixPreallocateSkip(projection->dm, PETSC_TRUE)); 448c85b835SJames Wright PetscCall(DMGetDimension(projection->dm, &dim)); 45d9e69621SJames Wright switch (diff_flux_proj_->method) { 468c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_DIRECT: { 47d9e69621SJames Wright projection->num_comp = diff_flux_proj_->num_diff_flux_comps; 488c85b835SJames Wright PetscCall(PetscObjectSetName((PetscObject)projection->dm, "DivDiffFluxProj")); 4936038bbcSJames Wright PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &projection->num_comp, projection->dm)); 508c85b835SJames Wright 510c373b74SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(honee->ceed, projection->dm, domain_label, label_value, height, dm_field, 52d9e69621SJames Wright &diff_flux_proj_->elem_restr_div_diff_flux)); 530c373b74SJames Wright PetscCallCeed(honee->ceed, 54d9e69621SJames Wright CeedElemRestrictionCreateVector(diff_flux_proj_->elem_restr_div_diff_flux, &diff_flux_proj_->div_diff_flux_ceed, NULL)); 55d9e69621SJames Wright PetscCall(CreateBasisFromPlex(honee->ceed, projection->dm, domain_label, label_value, height, dm_field, &diff_flux_proj_->basis_div_diff_flux)); 56d9e69621SJames Wright diff_flux_proj_->eval_mode_div_diff_flux = CEED_EVAL_INTERP; 571af555e8SJames Wright 581af555e8SJames Wright { // Create face labels on projection->dm for boundary integrals 591af555e8SJames Wright DMLabel face_sets_label; 601af555e8SJames Wright PetscInt num_face_set_values, *face_set_values; 611af555e8SJames Wright 620c373b74SJames Wright PetscCall(DMGetLabel(honee->dm, "Face Sets", &face_sets_label)); 630c373b74SJames Wright PetscCall(DMLabelCreateGlobalValueArray(honee->dm, face_sets_label, &num_face_set_values, &face_set_values)); 641af555e8SJames Wright for (PetscInt f = 0; f < num_face_set_values; f++) { 651af555e8SJames Wright DMLabel face_orientation_label; 661af555e8SJames Wright char *face_orientation_label_name; 671af555e8SJames Wright 680c373b74SJames Wright PetscCall(DMPlexCreateFaceLabel(honee->dm, face_set_values[f], &face_orientation_label_name)); 690c373b74SJames Wright PetscCall(DMGetLabel(honee->dm, face_orientation_label_name, &face_orientation_label)); 701af555e8SJames Wright PetscCall(DMAddLabel(projection->dm, face_orientation_label)); 711af555e8SJames Wright PetscCall(PetscFree(face_orientation_label_name)); 721af555e8SJames Wright } 7336038bbcSJames Wright PetscCall(PetscFree(face_set_values)); 741af555e8SJames Wright } 758c85b835SJames Wright } break; 768c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_INDIRECT: { 77d9e69621SJames Wright projection->num_comp = diff_flux_proj_->num_diff_flux_comps * dim; 788c85b835SJames Wright PetscCall(PetscObjectSetName((PetscObject)projection->dm, "DiffFluxProj")); 7936038bbcSJames Wright PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &projection->num_comp, projection->dm)); 808c85b835SJames Wright 810c373b74SJames Wright PetscCall(DMPlexCeedElemRestrictionQDataCreate(honee->ceed, projection->dm, domain_label, label_value, height, 82d9e69621SJames Wright diff_flux_proj_->num_diff_flux_comps, &diff_flux_proj_->elem_restr_div_diff_flux)); 830c373b74SJames Wright PetscCallCeed(honee->ceed, 84d9e69621SJames Wright CeedElemRestrictionCreateVector(diff_flux_proj_->elem_restr_div_diff_flux, &diff_flux_proj_->div_diff_flux_ceed, NULL)); 85d9e69621SJames Wright diff_flux_proj_->basis_div_diff_flux = CEED_BASIS_NONE; 86d9e69621SJames Wright diff_flux_proj_->eval_mode_div_diff_flux = CEED_EVAL_NONE; 878c85b835SJames Wright } break; 888c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_NONE: 89d9e69621SJames Wright SETERRQ(honee->comm, PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s", 90d9e69621SJames Wright DivDiffFluxProjectionMethods[divFdiffproj_method]); 918c85b835SJames Wright break; 928c85b835SJames Wright } 93d9e69621SJames Wright *diff_flux_proj = diff_flux_proj_; 948c85b835SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 958c85b835SJames Wright }; 968c85b835SJames Wright 978c85b835SJames Wright /** 980880fbb6SJames Wright @brief Return the objects required for the Divergence of Diffusive flux to be read by a `CeedOperator` 990880fbb6SJames Wright 1000880fbb6SJames Wright @param[in] diff_flux_proj Projection object 1010880fbb6SJames Wright @param[out] elem_restr Element restriction for the divergence of diffusive flux, or `NULL` 1020880fbb6SJames Wright @param[out] basis Basis for the divergence of diffusive flux, or `NULL` 1030880fbb6SJames Wright @param[out] vector Vector for the divergence of diffusive flux, or `NULL` 1040880fbb6SJames Wright @param[out] eval_mode Eval mode for the divergence of diffusive flux, or `NULL` 1050880fbb6SJames Wright **/ 1060880fbb6SJames Wright PetscErrorCode DivDiffFluxProjectionGetOperatorFieldData(DivDiffFluxProjectionData diff_flux_proj, CeedElemRestriction *elem_restr, CeedBasis *basis, 1070880fbb6SJames Wright CeedVector *vector, CeedEvalMode *eval_mode) { 1080880fbb6SJames Wright Ceed ceed = CeedVectorReturnCeed(diff_flux_proj->div_diff_flux_ceed); 1090880fbb6SJames Wright 1100880fbb6SJames Wright PetscFunctionBeginUser; 1110880fbb6SJames Wright if (elem_restr) PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(diff_flux_proj->elem_restr_div_diff_flux, elem_restr)); 1120880fbb6SJames Wright if (basis) PetscCallCeed(ceed, CeedBasisReferenceCopy(diff_flux_proj->basis_div_diff_flux, basis)); 1130880fbb6SJames Wright if (vector) PetscCallCeed(ceed, CeedVectorReferenceCopy(diff_flux_proj->div_diff_flux_ceed, vector)); 1140880fbb6SJames Wright if (eval_mode) *eval_mode = diff_flux_proj->eval_mode_div_diff_flux; 1150880fbb6SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1160880fbb6SJames Wright } 1170880fbb6SJames Wright 1180880fbb6SJames Wright /** 1198c85b835SJames Wright @brief Setup direct projection of divergence of diffusive flux 1208c85b835SJames Wright 1210c373b74SJames Wright @param[in] honee `Honee` context 1228561fee2SJames Wright @param[in,out] diff_flux_proj Flux projection object to setup 1238c85b835SJames Wright **/ 124e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionSetup_Direct(Honee honee, DivDiffFluxProjectionData diff_flux_proj) { 1250c373b74SJames Wright Ceed ceed = honee->ceed; 1268c85b835SJames Wright NodalProjectionData projection = diff_flux_proj->projection; 12736038bbcSJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)projection->dm); 1288c85b835SJames Wright 1298c85b835SJames Wright PetscFunctionBeginUser; 13036038bbcSJames Wright { // Create Projection RHS OperatorApplyContext 13136038bbcSJames Wright CeedOperator op_rhs; 1328c85b835SJames Wright 13336038bbcSJames Wright PetscCheck(diff_flux_proj->CreateRHSOperator_Direct, comm, PETSC_ERR_ARG_WRONGSTATE, 13436038bbcSJames Wright "Must define CreateRHSOperator_Direct to use indirect div_diff_flux projection"); 135e3663b90SJames Wright PetscCall(diff_flux_proj->CreateRHSOperator_Direct(honee, diff_flux_proj, &op_rhs)); 1368c85b835SJames Wright PetscCall(DMCreateLocalVector(projection->dm, &diff_flux_proj->DivDiffFlux_loc)); 1378c85b835SJames Wright diff_flux_proj->ceed_vec_has_array = PETSC_FALSE; 1380c373b74SJames Wright PetscCall(OperatorApplyContextCreate(honee->dm, projection->dm, ceed, op_rhs, NULL, NULL, NULL, diff_flux_proj->DivDiffFlux_loc, 13936038bbcSJames Wright &projection->l2_rhs_ctx)); 14036038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs)); 14136038bbcSJames Wright } 1428c85b835SJames Wright 1438c85b835SJames Wright { // -- Build Mass operator 1448c85b835SJames Wright CeedQFunction qf_mass; 1458c85b835SJames Wright CeedOperator op_mass; 146590d9cddSJames Wright CeedBasis basis_div_diff_flux = NULL; 147590d9cddSJames Wright CeedElemRestriction elem_restr_div_diff_flux_volume = NULL, elem_restr_qd; 1480880fbb6SJames Wright CeedVector q_data; 1490880fbb6SJames Wright CeedInt q_data_size; 1500880fbb6SJames Wright PetscInt label_value = 0; 1510880fbb6SJames Wright DMLabel domain_label = NULL; 1520880fbb6SJames Wright 153590d9cddSJames Wright PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_div_diff_flux_volume, &basis_div_diff_flux, NULL, NULL)); 154e3663b90SJames Wright PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, &q_data, 155e3663b90SJames Wright &q_data_size)); 1568c85b835SJames Wright 15764dd23feSJames Wright PetscCall(HoneeMassQFunctionCreate(ceed, projection->num_comp, q_data_size, &qf_mass)); 1588c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass)); 159590d9cddSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_div_diff_flux_volume, basis_div_diff_flux, CEED_VECTOR_ACTIVE)); 1608c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 161590d9cddSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_div_diff_flux_volume, basis_div_diff_flux, CEED_VECTOR_ACTIVE)); 1628c85b835SJames Wright 1638c85b835SJames Wright { // -- Setup KSP for L^2 projection 1648c85b835SJames Wright Mat mat_mass; 1658c85b835SJames Wright 1668c85b835SJames Wright PetscCall(MatCreateCeed(projection->dm, projection->dm, op_mass, NULL, &mat_mass)); 1678c85b835SJames Wright 1688c85b835SJames Wright PetscCall(KSPCreate(comm, &projection->ksp)); 1698c85b835SJames Wright PetscCall(KSPSetOptionsPrefix(projection->ksp, "div_diff_flux_projection_")); 1708c85b835SJames Wright { // lumped by default 1718c85b835SJames Wright PC pc; 1728c85b835SJames Wright PetscCall(KSPGetPC(projection->ksp, &pc)); 1738c85b835SJames Wright PetscCall(PCSetType(pc, PCJACOBI)); 1748c85b835SJames Wright PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 1758c85b835SJames Wright PetscCall(KSPSetType(projection->ksp, KSPPREONLY)); 1768c85b835SJames Wright } 1778c85b835SJames Wright PetscCall(KSPSetFromOptions_WithMatCeed(projection->ksp, mat_mass)); 1788c85b835SJames Wright PetscCall(MatDestroy(&mat_mass)); 1798c85b835SJames Wright } 1800880fbb6SJames Wright 181590d9cddSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_div_diff_flux_volume)); 182590d9cddSJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_div_diff_flux)); 1830880fbb6SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 1840880fbb6SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 1858c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 1868c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 1878c85b835SJames Wright } 1888c85b835SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1898c85b835SJames Wright } 1908c85b835SJames Wright 1918c85b835SJames Wright /** 1928c85b835SJames Wright @brief Setup indirect projection of divergence of diffusive flux 1938c85b835SJames Wright 1940c373b74SJames Wright @param[in] honee `Honee` context 1958561fee2SJames Wright @param[in,out] diff_flux_proj Flux projection object to setup 1968c85b835SJames Wright **/ 197e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionSetup_Indirect(Honee honee, DivDiffFluxProjectionData diff_flux_proj) { 1980c373b74SJames Wright Ceed ceed = honee->ceed; 1998c85b835SJames Wright NodalProjectionData projection = diff_flux_proj->projection; 2008c85b835SJames Wright CeedBasis basis_diff_flux; 2018c85b835SJames Wright CeedElemRestriction elem_restr_diff_flux, elem_restr_qd; 2028c85b835SJames Wright CeedVector q_data; 2030880fbb6SJames Wright CeedInt q_data_size; 20436038bbcSJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)projection->dm); 2058c85b835SJames Wright 2068c85b835SJames Wright PetscFunctionBeginUser; 2070880fbb6SJames Wright { 2080880fbb6SJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 2090880fbb6SJames Wright DMLabel domain_label = NULL; 2108c85b835SJames Wright 2118c85b835SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, domain_label, label_value, height, dm_field, &elem_restr_diff_flux)); 2128c85b835SJames Wright PetscCall(CreateBasisFromPlex(ceed, projection->dm, domain_label, label_value, height, dm_field, &basis_diff_flux)); 213e3663b90SJames Wright PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, &q_data, 214e3663b90SJames Wright &q_data_size)); 2150880fbb6SJames Wright } 2168c85b835SJames Wright 21736038bbcSJames Wright { 2188c85b835SJames Wright CeedOperator op_rhs; 2198c85b835SJames Wright 22036038bbcSJames Wright PetscCheck(diff_flux_proj->CreateRHSOperator_Indirect, comm, PETSC_ERR_ARG_WRONGSTATE, 22136038bbcSJames Wright "Must define CreateRHSOperator_Indirect to use indirect div_diff_flux projection"); 222e3663b90SJames Wright PetscCall(diff_flux_proj->CreateRHSOperator_Indirect(honee, diff_flux_proj, &op_rhs)); 2230c373b74SJames Wright PetscCall(OperatorApplyContextCreate(honee->dm, projection->dm, ceed, op_rhs, NULL, NULL, NULL, NULL, &projection->l2_rhs_ctx)); 2248c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs)); 2258c85b835SJames Wright } 2268c85b835SJames Wright 2278c85b835SJames Wright { // -- Build Mass operator 2288c85b835SJames Wright CeedQFunction qf_mass; 2298c85b835SJames Wright CeedOperator op_mass; 2308c85b835SJames Wright 23164dd23feSJames Wright PetscCall(HoneeMassQFunctionCreate(ceed, projection->num_comp, q_data_size, &qf_mass)); 2328c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass)); 2338c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE)); 2348c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 2358c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE)); 2368c85b835SJames Wright 2378c85b835SJames Wright { // -- Setup KSP for L^2 projection 2388c85b835SJames Wright Mat mat_mass; 2398c85b835SJames Wright 2408c85b835SJames Wright PetscCall(MatCreateCeed(projection->dm, projection->dm, op_mass, NULL, &mat_mass)); 2418c85b835SJames Wright 2428c85b835SJames Wright PetscCall(KSPCreate(comm, &projection->ksp)); 2438c85b835SJames Wright PetscCall(KSPSetOptionsPrefix(projection->ksp, "div_diff_flux_projection_")); 2448c85b835SJames Wright { // lumped by default 2458c85b835SJames Wright PC pc; 2468c85b835SJames Wright PetscCall(KSPGetPC(projection->ksp, &pc)); 2478c85b835SJames Wright PetscCall(PCSetType(pc, PCJACOBI)); 2488c85b835SJames Wright PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 2498c85b835SJames Wright PetscCall(KSPSetType(projection->ksp, KSPPREONLY)); 2508c85b835SJames Wright } 2518c85b835SJames Wright PetscCall(KSPSetFromOptions_WithMatCeed(projection->ksp, mat_mass)); 2528c85b835SJames Wright PetscCall(MatDestroy(&mat_mass)); 2538c85b835SJames Wright } 2548c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 2558c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 2568c85b835SJames Wright } 2578c85b835SJames Wright 2588c85b835SJames Wright { // Create OperatorApplyContext to calculate divergence at quadrature points 25940b78511SJames Wright CeedQFunction qf_calc_divergence = NULL; 2608c85b835SJames Wright CeedOperator op_calc_divergence; 2610880fbb6SJames Wright CeedElemRestriction elem_restr_div_diff_flux = NULL; 2620880fbb6SJames Wright PetscInt dim; 2638c85b835SJames Wright 2640880fbb6SJames Wright PetscCall(DMGetDimension(projection->dm, &dim)); 2650880fbb6SJames Wright PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_div_diff_flux, NULL, NULL, NULL)); 2668c85b835SJames Wright 26740b78511SJames Wright switch (dim) { 26840b78511SJames Wright case 2: 26940b78511SJames Wright switch (diff_flux_proj->num_diff_flux_comps) { 27040b78511SJames Wright case 1: 27140b78511SJames Wright PetscCallCeed(ceed, 27240b78511SJames Wright CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux2D_1, ComputeDivDiffusiveFlux2D_1_loc, &qf_calc_divergence)); 27340b78511SJames Wright break; 27440b78511SJames Wright } 27540b78511SJames Wright break; 27640b78511SJames Wright case 3: 27740b78511SJames Wright switch (diff_flux_proj->num_diff_flux_comps) { 27840b78511SJames Wright case 1: 27940b78511SJames Wright PetscCallCeed(ceed, 28040b78511SJames Wright CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux3D_1, ComputeDivDiffusiveFlux3D_1_loc, &qf_calc_divergence)); 28140b78511SJames Wright break; 28240b78511SJames Wright case 4: 28340b78511SJames Wright PetscCallCeed(ceed, 28440b78511SJames Wright CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux3D_4, ComputeDivDiffusiveFlux3D_4_loc, &qf_calc_divergence)); 28540b78511SJames Wright break; 28640b78511SJames Wright } 28740b78511SJames Wright break; 28840b78511SJames Wright } 28940b78511SJames Wright PetscCheck(qf_calc_divergence, comm, PETSC_ERR_SUP, 29040b78511SJames Wright "QFunction for calculating divergence of diffusive flux does not exist for" 29140b78511SJames Wright " %" PetscInt_FMT " dimensional grid and %" PetscInt_FMT 29240b78511SJames Wright " number of components.\nA new qfunction can be easily added; see source code for pattern.", 29340b78511SJames Wright dim, diff_flux_proj->num_diff_flux_comps); 2948c85b835SJames Wright 2958c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_calc_divergence, "Grad F_diff", projection->num_comp * dim, CEED_EVAL_GRAD)); 2968c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_calc_divergence, "qdata", q_data_size, CEED_EVAL_NONE)); 29740b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_calc_divergence, "Div F_diff", diff_flux_proj->num_diff_flux_comps, CEED_EVAL_NONE)); 2988c85b835SJames Wright 2998c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_calc_divergence, NULL, NULL, &op_calc_divergence)); 3008c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "Grad F_diff", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE)); 3018c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 3028c85b835SJames Wright PetscCallCeed( 3038c85b835SJames Wright ceed, CeedOperatorSetField(op_calc_divergence, "Div F_diff", elem_restr_div_diff_flux, CEED_BASIS_NONE, diff_flux_proj->div_diff_flux_ceed)); 3048c85b835SJames Wright 305*65fa3c0aSJames Wright PetscCall(OperatorApplyContextCreate(projection->dm, NULL, ceed, op_calc_divergence, NULL, CEED_VECTOR_NONE, NULL, NULL, 306*65fa3c0aSJames Wright &diff_flux_proj->calc_div_diff_flux)); 3078c85b835SJames Wright 3080880fbb6SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_div_diff_flux)); 3098c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_calc_divergence)); 3108c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_calc_divergence)); 3118c85b835SJames Wright } 3128c85b835SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 3138c85b835SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 3148c85b835SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 3158c85b835SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux)); 3168c85b835SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3178c85b835SJames Wright } 3188c85b835SJames Wright 3198c85b835SJames Wright /** 3208c85b835SJames Wright @brief Setup projection of divergence of diffusive flux 3218c85b835SJames Wright 3220c373b74SJames Wright @param[in] honee `Honee` context 3238561fee2SJames Wright @param[in,out] diff_flux_proj Flux projection object to setup 3248c85b835SJames Wright **/ 325e3663b90SJames Wright PetscErrorCode DivDiffFluxProjectionSetup(Honee honee, DivDiffFluxProjectionData diff_flux_proj) { 3268c85b835SJames Wright PetscFunctionBeginUser; 3270c373b74SJames Wright switch (honee->app_ctx->divFdiffproj_method) { 3288c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_DIRECT: 329e3663b90SJames Wright PetscCall(DivDiffFluxProjectionSetup_Direct(honee, diff_flux_proj)); 3308c85b835SJames Wright break; 3318c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_INDIRECT: 332e3663b90SJames Wright PetscCall(DivDiffFluxProjectionSetup_Indirect(honee, diff_flux_proj)); 3338c85b835SJames Wright break; 3348c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_NONE: 3350c373b74SJames Wright SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s", 3360c373b74SJames Wright DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]); 3378c85b835SJames Wright break; 3388c85b835SJames Wright } 3398c85b835SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3408c85b835SJames Wright } 3418c85b835SJames Wright 3428c85b835SJames Wright /** 3438c85b835SJames Wright @brief Project the divergence of diffusive flux 3448c85b835SJames Wright 3458c85b835SJames Wright This implicitly sets the `CeedVector` input (`div_diff_flux_ceed`) to the divergence of diffusive flux. 3468c85b835SJames Wright 3478c85b835SJames Wright @param[in] diff_flux_proj `NodalProjectionData` for the projection 3488c85b835SJames Wright @param[in] Q_loc Localized solution vector 3498c85b835SJames Wright **/ 35036038bbcSJames Wright PetscErrorCode DivDiffFluxProjectionApply(DivDiffFluxProjectionData diff_flux_proj, Vec Q_loc) { 3518c85b835SJames Wright NodalProjectionData projection = diff_flux_proj->projection; 3528c85b835SJames Wright 3538c85b835SJames Wright PetscFunctionBeginUser; 354ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_DivDiffFluxProjection, Q_loc, 0, 0, 0)); 3558c85b835SJames Wright switch (diff_flux_proj->method) { 3568c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_DIRECT: { 357c8cb0c9bSJames Wright Vec DivDiffFlux, RHS; 3588c85b835SJames Wright 3598c85b835SJames Wright PetscCall(DMGetGlobalVector(projection->dm, &DivDiffFlux)); 360c8cb0c9bSJames Wright PetscCall(DMGetGlobalVector(projection->dm, &RHS)); 3618c85b835SJames Wright if (diff_flux_proj->ceed_vec_has_array) { 3628c85b835SJames Wright PetscCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->DivDiffFlux_loc)); 3638c85b835SJames Wright diff_flux_proj->ceed_vec_has_array = PETSC_FALSE; 3648c85b835SJames Wright } 365c8cb0c9bSJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, RHS, projection->l2_rhs_ctx)); 3668c85b835SJames Wright PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_rhs_view")); 3678c85b835SJames Wright 368c8cb0c9bSJames Wright { 369c8cb0c9bSJames Wright // Run PCApply manually if using ksp_type preonly -pc_type jacobi 370c8cb0c9bSJames Wright // This is to avoid an AllReduce call in KSPSolve_Preonly, which causes significant slowdowns for lumped mass matrix solves. 371c8cb0c9bSJames Wright // See https://gitlab.com/petsc/petsc/-/merge_requests/8048 for more details and a possible fix 372c8cb0c9bSJames Wright PC pc; 373c8cb0c9bSJames Wright PetscBool ispreonly, isjacobi; 374c8cb0c9bSJames Wright PetscCall(KSPGetPC(projection->ksp, &pc)); 375c8cb0c9bSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)projection->ksp, KSPPREONLY, &ispreonly)); 376c8cb0c9bSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &isjacobi)); 377c8cb0c9bSJames Wright if (ispreonly && isjacobi) PetscCall(PCApply(pc, RHS, DivDiffFlux)); 378c8cb0c9bSJames Wright else PetscCall(KSPSolve(projection->ksp, RHS, DivDiffFlux)); 379c8cb0c9bSJames Wright } 3808c85b835SJames Wright PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_view")); 3818c85b835SJames Wright 3828c85b835SJames Wright PetscCall(DMGlobalToLocal(projection->dm, DivDiffFlux, INSERT_VALUES, diff_flux_proj->DivDiffFlux_loc)); 3838c85b835SJames Wright PetscCall(VecReadPetscToCeed(diff_flux_proj->DivDiffFlux_loc, &diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->div_diff_flux_ceed)); 3848c85b835SJames Wright diff_flux_proj->ceed_vec_has_array = PETSC_TRUE; 3858c85b835SJames Wright 386c8cb0c9bSJames Wright PetscCall(DMRestoreGlobalVector(projection->dm, &RHS)); 3878c85b835SJames Wright PetscCall(DMRestoreGlobalVector(projection->dm, &DivDiffFlux)); 3888c85b835SJames Wright break; 3898c85b835SJames Wright } 3908c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_INDIRECT: { 391c8cb0c9bSJames Wright Vec DiffFlux, RHS; 3928c85b835SJames Wright 3938c85b835SJames Wright PetscCall(DMGetGlobalVector(projection->dm, &DiffFlux)); 394c8cb0c9bSJames Wright PetscCall(DMGetGlobalVector(projection->dm, &RHS)); 395c8cb0c9bSJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, RHS, projection->l2_rhs_ctx)); 3968c85b835SJames Wright PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_rhs_view")); 3978c85b835SJames Wright 398c8cb0c9bSJames Wright { 399c8cb0c9bSJames Wright // Run PCApply manually if using -ksp_type preonly -pc_type jacobi 400c8cb0c9bSJames Wright // This is to avoid an AllReduce call in KSPSolve_Preonly, which causes significant slowdowns for lumped mass matrix solves. 401c8cb0c9bSJames Wright // See https://gitlab.com/petsc/petsc/-/merge_requests/8048 for more details and a possible fix 402c8cb0c9bSJames Wright PC pc; 403c8cb0c9bSJames Wright PetscBool ispreonly, isjacobi; 404c8cb0c9bSJames Wright PetscCall(KSPGetPC(projection->ksp, &pc)); 405c8cb0c9bSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)projection->ksp, KSPPREONLY, &ispreonly)); 406c8cb0c9bSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &isjacobi)); 407c8cb0c9bSJames Wright if (ispreonly && isjacobi) PetscCall(PCApply(pc, RHS, DiffFlux)); 408c8cb0c9bSJames Wright else PetscCall(KSPSolve(projection->ksp, RHS, DiffFlux)); 409c8cb0c9bSJames Wright } 4108c85b835SJames Wright PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_view")); 4118c85b835SJames Wright 4128c85b835SJames Wright PetscCall(ApplyCeedOperatorGlobalToLocal(DiffFlux, NULL, diff_flux_proj->calc_div_diff_flux)); 413c8cb0c9bSJames Wright PetscCall(DMRestoreGlobalVector(projection->dm, &RHS)); 4148c85b835SJames Wright PetscCall(DMRestoreGlobalVector(projection->dm, &DiffFlux)); 4158c85b835SJames Wright } break; 4168c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_NONE: 4178c85b835SJames Wright SETERRQ(PetscObjectComm((PetscObject)projection->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s", 4188c85b835SJames Wright DivDiffFluxProjectionMethods[diff_flux_proj->method]); 4198c85b835SJames Wright break; 4208c85b835SJames Wright } 421ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_DivDiffFluxProjection, Q_loc, 0, 0, 0)); 4228c85b835SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 4238c85b835SJames Wright } 4248c85b835SJames Wright 4258c85b835SJames Wright /** 4268c85b835SJames Wright @brief Destroy `DivDiffFluxProjectionData` object 4278c85b835SJames Wright 4288c85b835SJames Wright @param[in,out] diff_flux_proj Object to destroy 4298c85b835SJames Wright **/ 4308c85b835SJames Wright PetscErrorCode DivDiffFluxProjectionDataDestroy(DivDiffFluxProjectionData diff_flux_proj) { 4318c85b835SJames Wright PetscFunctionBeginUser; 4328c85b835SJames Wright if (diff_flux_proj == NULL) PetscFunctionReturn(PETSC_SUCCESS); 4330880fbb6SJames Wright Ceed ceed = CeedVectorReturnCeed(diff_flux_proj->div_diff_flux_ceed); 4340880fbb6SJames Wright 4354ea616f4SJames Wright PetscCall(NodalProjectionDataDestroy(&diff_flux_proj->projection)); 4368c85b835SJames Wright PetscCall(OperatorApplyContextDestroy(diff_flux_proj->calc_div_diff_flux)); 4378c85b835SJames Wright if (diff_flux_proj->ceed_vec_has_array) { 4388c85b835SJames Wright PetscCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->DivDiffFlux_loc)); 4398c85b835SJames Wright diff_flux_proj->ceed_vec_has_array = PETSC_FALSE; 4408c85b835SJames Wright } 4410880fbb6SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&diff_flux_proj->div_diff_flux_ceed)); 4420880fbb6SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&diff_flux_proj->elem_restr_div_diff_flux)); 4430880fbb6SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&diff_flux_proj->basis_div_diff_flux)); 4448c85b835SJames Wright PetscCall(VecDestroy(&diff_flux_proj->DivDiffFlux_loc)); 4458c85b835SJames Wright PetscCall(PetscFree(diff_flux_proj)); 4468c85b835SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 4478c85b835SJames Wright } 448