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 1836038bbcSJames Wright @param[in] num_diff_flux_comps Number of components that makes up the diffusive flux (e.g. 1 for scalar advection-diffusion) 1936038bbcSJames Wright @param[out] pdiff_flux_proj The `DivDiffFluxProjectionData` object created, or `NULL` if not created 208c85b835SJames Wright **/ 210c373b74SJames Wright PetscErrorCode DivDiffFluxProjectionCreate(Honee honee, PetscInt num_diff_flux_comps, DivDiffFluxProjectionData *pdiff_flux_proj) { 220c373b74SJames Wright PetscInt label_value = 0, height = 0, dm_field = 0, dim, degree = honee->app_ctx->degree, q_extra = honee->app_ctx->q_extra; 238c85b835SJames Wright DMLabel domain_label = NULL; 248c85b835SJames Wright DivDiffFluxProjectionData diff_flux_proj; 258c85b835SJames Wright NodalProjectionData projection; 268c85b835SJames Wright 278c85b835SJames Wright PetscFunctionBeginUser; 280c373b74SJames Wright if (honee->app_ctx->divFdiffproj_method == DIV_DIFF_FLUX_PROJ_NONE) { 2936038bbcSJames Wright *pdiff_flux_proj = NULL; 3036038bbcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3136038bbcSJames Wright } 3236038bbcSJames Wright PetscCall(PetscNew(pdiff_flux_proj)); 3336038bbcSJames Wright diff_flux_proj = *pdiff_flux_proj; 340c373b74SJames Wright PetscCall(PetscNew(&honee->diff_flux_proj->projection)); 350c373b74SJames Wright projection = honee->diff_flux_proj->projection; 360c373b74SJames Wright diff_flux_proj->method = honee->app_ctx->divFdiffproj_method; 3736038bbcSJames Wright diff_flux_proj->num_diff_flux_comps = num_diff_flux_comps; 388c85b835SJames Wright 390c373b74SJames Wright PetscCall(DMClone(honee->dm, &projection->dm)); 408c85b835SJames Wright PetscCall(DMSetMatrixPreallocateSkip(projection->dm, PETSC_TRUE)); 418c85b835SJames Wright PetscCall(DMGetDimension(projection->dm, &dim)); 428c85b835SJames Wright switch (diff_flux_proj->method) { 438c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_DIRECT: { 448c85b835SJames Wright projection->num_comp = diff_flux_proj->num_diff_flux_comps; 458c85b835SJames Wright PetscCall(PetscObjectSetName((PetscObject)projection->dm, "DivDiffFluxProj")); 4636038bbcSJames Wright PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &projection->num_comp, projection->dm)); 478c85b835SJames Wright 480c373b74SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(honee->ceed, projection->dm, domain_label, label_value, height, dm_field, 4936038bbcSJames Wright &diff_flux_proj->elem_restr_div_diff_flux)); 500c373b74SJames Wright PetscCallCeed(honee->ceed, 510c373b74SJames Wright CeedElemRestrictionCreateVector(diff_flux_proj->elem_restr_div_diff_flux, &diff_flux_proj->div_diff_flux_ceed, NULL)); 520c373b74SJames Wright PetscCall(CreateBasisFromPlex(honee->ceed, projection->dm, domain_label, label_value, height, dm_field, &diff_flux_proj->basis_div_diff_flux)); 5336038bbcSJames Wright diff_flux_proj->eval_mode_div_diff_flux = CEED_EVAL_INTERP; 541af555e8SJames Wright 551af555e8SJames Wright { // Create face labels on projection->dm for boundary integrals 561af555e8SJames Wright DMLabel face_sets_label; 571af555e8SJames Wright PetscInt num_face_set_values, *face_set_values; 581af555e8SJames Wright 590c373b74SJames Wright PetscCall(DMGetLabel(honee->dm, "Face Sets", &face_sets_label)); 600c373b74SJames Wright PetscCall(DMLabelCreateGlobalValueArray(honee->dm, face_sets_label, &num_face_set_values, &face_set_values)); 611af555e8SJames Wright for (PetscInt f = 0; f < num_face_set_values; f++) { 621af555e8SJames Wright DMLabel face_orientation_label; 631af555e8SJames Wright char *face_orientation_label_name; 641af555e8SJames Wright 650c373b74SJames Wright PetscCall(DMPlexCreateFaceLabel(honee->dm, face_set_values[f], &face_orientation_label_name)); 660c373b74SJames Wright PetscCall(DMGetLabel(honee->dm, face_orientation_label_name, &face_orientation_label)); 671af555e8SJames Wright PetscCall(DMAddLabel(projection->dm, face_orientation_label)); 681af555e8SJames Wright PetscCall(PetscFree(face_orientation_label_name)); 691af555e8SJames Wright } 7036038bbcSJames Wright PetscCall(PetscFree(face_set_values)); 711af555e8SJames Wright } 728c85b835SJames Wright } break; 738c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_INDIRECT: { 748c85b835SJames Wright projection->num_comp = diff_flux_proj->num_diff_flux_comps * dim; 758c85b835SJames Wright PetscCall(PetscObjectSetName((PetscObject)projection->dm, "DiffFluxProj")); 7636038bbcSJames Wright PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &projection->num_comp, projection->dm)); 778c85b835SJames Wright 780c373b74SJames Wright PetscCall(DMPlexCeedElemRestrictionQDataCreate(honee->ceed, projection->dm, domain_label, label_value, height, 7936038bbcSJames Wright diff_flux_proj->num_diff_flux_comps, &diff_flux_proj->elem_restr_div_diff_flux)); 800c373b74SJames Wright PetscCallCeed(honee->ceed, 810c373b74SJames Wright CeedElemRestrictionCreateVector(diff_flux_proj->elem_restr_div_diff_flux, &diff_flux_proj->div_diff_flux_ceed, NULL)); 8236038bbcSJames Wright diff_flux_proj->basis_div_diff_flux = CEED_BASIS_NONE; 8336038bbcSJames Wright diff_flux_proj->eval_mode_div_diff_flux = CEED_EVAL_NONE; 848c85b835SJames Wright } break; 858c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_NONE: 860c373b74SJames Wright SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s", 870c373b74SJames Wright DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]); 888c85b835SJames Wright break; 898c85b835SJames Wright } 908c85b835SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 918c85b835SJames Wright }; 928c85b835SJames Wright 938c85b835SJames Wright /** 940880fbb6SJames Wright @brief Return the objects required for the Divergence of Diffusive flux to be read by a `CeedOperator` 950880fbb6SJames Wright 960880fbb6SJames Wright @param[in] diff_flux_proj Projection object 970880fbb6SJames Wright @param[out] elem_restr Element restriction for the divergence of diffusive flux, or `NULL` 980880fbb6SJames Wright @param[out] basis Basis for the divergence of diffusive flux, or `NULL` 990880fbb6SJames Wright @param[out] vector Vector for the divergence of diffusive flux, or `NULL` 1000880fbb6SJames Wright @param[out] eval_mode Eval mode for the divergence of diffusive flux, or `NULL` 1010880fbb6SJames Wright **/ 1020880fbb6SJames Wright PetscErrorCode DivDiffFluxProjectionGetOperatorFieldData(DivDiffFluxProjectionData diff_flux_proj, CeedElemRestriction *elem_restr, CeedBasis *basis, 1030880fbb6SJames Wright CeedVector *vector, CeedEvalMode *eval_mode) { 1040880fbb6SJames Wright Ceed ceed = CeedVectorReturnCeed(diff_flux_proj->div_diff_flux_ceed); 1050880fbb6SJames Wright 1060880fbb6SJames Wright PetscFunctionBeginUser; 1070880fbb6SJames Wright if (elem_restr) PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(diff_flux_proj->elem_restr_div_diff_flux, elem_restr)); 1080880fbb6SJames Wright if (basis) PetscCallCeed(ceed, CeedBasisReferenceCopy(diff_flux_proj->basis_div_diff_flux, basis)); 1090880fbb6SJames Wright if (vector) PetscCallCeed(ceed, CeedVectorReferenceCopy(diff_flux_proj->div_diff_flux_ceed, vector)); 1100880fbb6SJames Wright if (eval_mode) *eval_mode = diff_flux_proj->eval_mode_div_diff_flux; 1110880fbb6SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1120880fbb6SJames Wright } 1130880fbb6SJames Wright 1140880fbb6SJames Wright /** 1158c85b835SJames Wright @brief Setup direct projection of divergence of diffusive flux 1168c85b835SJames Wright 1170c373b74SJames Wright @param[in] honee `Honee` context 1188561fee2SJames Wright @param[in,out] diff_flux_proj Flux projection object to setup 1198c85b835SJames Wright **/ 120e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionSetup_Direct(Honee honee, DivDiffFluxProjectionData diff_flux_proj) { 1210c373b74SJames Wright Ceed ceed = honee->ceed; 1228c85b835SJames Wright NodalProjectionData projection = diff_flux_proj->projection; 12336038bbcSJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)projection->dm); 1248c85b835SJames Wright 1258c85b835SJames Wright PetscFunctionBeginUser; 12636038bbcSJames Wright { // Create Projection RHS OperatorApplyContext 12736038bbcSJames Wright CeedOperator op_rhs; 1288c85b835SJames Wright 12936038bbcSJames Wright PetscCheck(diff_flux_proj->CreateRHSOperator_Direct, comm, PETSC_ERR_ARG_WRONGSTATE, 13036038bbcSJames Wright "Must define CreateRHSOperator_Direct to use indirect div_diff_flux projection"); 131e3663b90SJames Wright PetscCall(diff_flux_proj->CreateRHSOperator_Direct(honee, diff_flux_proj, &op_rhs)); 1328c85b835SJames Wright PetscCall(DMCreateLocalVector(projection->dm, &diff_flux_proj->DivDiffFlux_loc)); 1338c85b835SJames Wright diff_flux_proj->ceed_vec_has_array = PETSC_FALSE; 1340c373b74SJames Wright PetscCall(OperatorApplyContextCreate(honee->dm, projection->dm, ceed, op_rhs, NULL, NULL, NULL, diff_flux_proj->DivDiffFlux_loc, 13536038bbcSJames Wright &projection->l2_rhs_ctx)); 13636038bbcSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs)); 13736038bbcSJames Wright } 1388c85b835SJames Wright 1398c85b835SJames Wright { // -- Build Mass operator 1408c85b835SJames Wright CeedQFunction qf_mass; 1418c85b835SJames Wright CeedOperator op_mass; 142590d9cddSJames Wright CeedBasis basis_div_diff_flux = NULL; 143590d9cddSJames Wright CeedElemRestriction elem_restr_div_diff_flux_volume = NULL, elem_restr_qd; 1440880fbb6SJames Wright CeedVector q_data; 1450880fbb6SJames Wright CeedInt q_data_size; 1460880fbb6SJames Wright PetscInt label_value = 0; 1470880fbb6SJames Wright DMLabel domain_label = NULL; 1480880fbb6SJames Wright 149590d9cddSJames Wright PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_div_diff_flux_volume, &basis_div_diff_flux, NULL, NULL)); 150e3663b90SJames 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, 151e3663b90SJames Wright &q_data_size)); 1528c85b835SJames Wright 15364dd23feSJames Wright PetscCall(HoneeMassQFunctionCreate(ceed, projection->num_comp, q_data_size, &qf_mass)); 1548c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass)); 155590d9cddSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_div_diff_flux_volume, basis_div_diff_flux, CEED_VECTOR_ACTIVE)); 1568c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 157590d9cddSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_div_diff_flux_volume, basis_div_diff_flux, CEED_VECTOR_ACTIVE)); 1588c85b835SJames Wright 1598c85b835SJames Wright { // -- Setup KSP for L^2 projection 1608c85b835SJames Wright Mat mat_mass; 1618c85b835SJames Wright 1628c85b835SJames Wright PetscCall(MatCreateCeed(projection->dm, projection->dm, op_mass, NULL, &mat_mass)); 1638c85b835SJames Wright 1648c85b835SJames Wright PetscCall(KSPCreate(comm, &projection->ksp)); 1658c85b835SJames Wright PetscCall(KSPSetOptionsPrefix(projection->ksp, "div_diff_flux_projection_")); 1668c85b835SJames Wright { // lumped by default 1678c85b835SJames Wright PC pc; 1688c85b835SJames Wright PetscCall(KSPGetPC(projection->ksp, &pc)); 1698c85b835SJames Wright PetscCall(PCSetType(pc, PCJACOBI)); 1708c85b835SJames Wright PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 1718c85b835SJames Wright PetscCall(KSPSetType(projection->ksp, KSPPREONLY)); 1728c85b835SJames Wright } 1738c85b835SJames Wright PetscCall(KSPSetFromOptions_WithMatCeed(projection->ksp, mat_mass)); 1748c85b835SJames Wright PetscCall(MatDestroy(&mat_mass)); 1758c85b835SJames Wright } 1760880fbb6SJames Wright 177590d9cddSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_div_diff_flux_volume)); 178590d9cddSJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_div_diff_flux)); 1790880fbb6SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 1800880fbb6SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 1818c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 1828c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 1838c85b835SJames Wright } 1848c85b835SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1858c85b835SJames Wright } 1868c85b835SJames Wright 1878c85b835SJames Wright /** 1888c85b835SJames Wright @brief Setup indirect projection of divergence of diffusive flux 1898c85b835SJames Wright 1900c373b74SJames Wright @param[in] honee `Honee` context 1918561fee2SJames Wright @param[in,out] diff_flux_proj Flux projection object to setup 1928c85b835SJames Wright **/ 193e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionSetup_Indirect(Honee honee, DivDiffFluxProjectionData diff_flux_proj) { 1940c373b74SJames Wright Ceed ceed = honee->ceed; 1958c85b835SJames Wright NodalProjectionData projection = diff_flux_proj->projection; 1968c85b835SJames Wright CeedBasis basis_diff_flux; 1978c85b835SJames Wright CeedElemRestriction elem_restr_diff_flux, elem_restr_qd; 1988c85b835SJames Wright CeedVector q_data; 1990880fbb6SJames Wright CeedInt q_data_size; 20036038bbcSJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)projection->dm); 2018c85b835SJames Wright 2028c85b835SJames Wright PetscFunctionBeginUser; 2030880fbb6SJames Wright { 2040880fbb6SJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 2050880fbb6SJames Wright DMLabel domain_label = NULL; 2068c85b835SJames Wright 2078c85b835SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, domain_label, label_value, height, dm_field, &elem_restr_diff_flux)); 2088c85b835SJames Wright PetscCall(CreateBasisFromPlex(ceed, projection->dm, domain_label, label_value, height, dm_field, &basis_diff_flux)); 209e3663b90SJames 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, 210e3663b90SJames Wright &q_data_size)); 2110880fbb6SJames Wright } 2128c85b835SJames Wright 21336038bbcSJames Wright { 2148c85b835SJames Wright CeedOperator op_rhs; 2158c85b835SJames Wright 21636038bbcSJames Wright PetscCheck(diff_flux_proj->CreateRHSOperator_Indirect, comm, PETSC_ERR_ARG_WRONGSTATE, 21736038bbcSJames Wright "Must define CreateRHSOperator_Indirect to use indirect div_diff_flux projection"); 218e3663b90SJames Wright PetscCall(diff_flux_proj->CreateRHSOperator_Indirect(honee, diff_flux_proj, &op_rhs)); 2190c373b74SJames Wright PetscCall(OperatorApplyContextCreate(honee->dm, projection->dm, ceed, op_rhs, NULL, NULL, NULL, NULL, &projection->l2_rhs_ctx)); 2208c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs)); 2218c85b835SJames Wright } 2228c85b835SJames Wright 2238c85b835SJames Wright { // -- Build Mass operator 2248c85b835SJames Wright CeedQFunction qf_mass; 2258c85b835SJames Wright CeedOperator op_mass; 2268c85b835SJames Wright 22764dd23feSJames Wright PetscCall(HoneeMassQFunctionCreate(ceed, projection->num_comp, q_data_size, &qf_mass)); 2288c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass)); 2298c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE)); 2308c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 2318c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE)); 2328c85b835SJames Wright 2338c85b835SJames Wright { // -- Setup KSP for L^2 projection 2348c85b835SJames Wright Mat mat_mass; 2358c85b835SJames Wright 2368c85b835SJames Wright PetscCall(MatCreateCeed(projection->dm, projection->dm, op_mass, NULL, &mat_mass)); 2378c85b835SJames Wright 2388c85b835SJames Wright PetscCall(KSPCreate(comm, &projection->ksp)); 2398c85b835SJames Wright PetscCall(KSPSetOptionsPrefix(projection->ksp, "div_diff_flux_projection_")); 2408c85b835SJames Wright { // lumped by default 2418c85b835SJames Wright PC pc; 2428c85b835SJames Wright PetscCall(KSPGetPC(projection->ksp, &pc)); 2438c85b835SJames Wright PetscCall(PCSetType(pc, PCJACOBI)); 2448c85b835SJames Wright PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 2458c85b835SJames Wright PetscCall(KSPSetType(projection->ksp, KSPPREONLY)); 2468c85b835SJames Wright } 2478c85b835SJames Wright PetscCall(KSPSetFromOptions_WithMatCeed(projection->ksp, mat_mass)); 2488c85b835SJames Wright PetscCall(MatDestroy(&mat_mass)); 2498c85b835SJames Wright } 2508c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 2518c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 2528c85b835SJames Wright } 2538c85b835SJames Wright 2548c85b835SJames Wright { // Create OperatorApplyContext to calculate divergence at quadrature points 25540b78511SJames Wright CeedQFunction qf_calc_divergence = NULL; 2568c85b835SJames Wright CeedOperator op_calc_divergence; 2570880fbb6SJames Wright CeedElemRestriction elem_restr_div_diff_flux = NULL; 2580880fbb6SJames Wright PetscInt dim; 2598c85b835SJames Wright 2600880fbb6SJames Wright PetscCall(DMGetDimension(projection->dm, &dim)); 2610880fbb6SJames Wright PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_div_diff_flux, NULL, NULL, NULL)); 2628c85b835SJames Wright 26340b78511SJames Wright switch (dim) { 26440b78511SJames Wright case 2: 26540b78511SJames Wright switch (diff_flux_proj->num_diff_flux_comps) { 26640b78511SJames Wright case 1: 26740b78511SJames Wright PetscCallCeed(ceed, 26840b78511SJames Wright CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux2D_1, ComputeDivDiffusiveFlux2D_1_loc, &qf_calc_divergence)); 26940b78511SJames Wright break; 27040b78511SJames Wright } 27140b78511SJames Wright break; 27240b78511SJames Wright case 3: 27340b78511SJames Wright switch (diff_flux_proj->num_diff_flux_comps) { 27440b78511SJames Wright case 1: 27540b78511SJames Wright PetscCallCeed(ceed, 27640b78511SJames Wright CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux3D_1, ComputeDivDiffusiveFlux3D_1_loc, &qf_calc_divergence)); 27740b78511SJames Wright break; 27840b78511SJames Wright case 4: 27940b78511SJames Wright PetscCallCeed(ceed, 28040b78511SJames Wright CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux3D_4, ComputeDivDiffusiveFlux3D_4_loc, &qf_calc_divergence)); 28140b78511SJames Wright break; 28240b78511SJames Wright } 28340b78511SJames Wright break; 28440b78511SJames Wright } 28540b78511SJames Wright PetscCheck(qf_calc_divergence, comm, PETSC_ERR_SUP, 28640b78511SJames Wright "QFunction for calculating divergence of diffusive flux does not exist for" 28740b78511SJames Wright " %" PetscInt_FMT " dimensional grid and %" PetscInt_FMT 28840b78511SJames Wright " number of components.\nA new qfunction can be easily added; see source code for pattern.", 28940b78511SJames Wright dim, diff_flux_proj->num_diff_flux_comps); 2908c85b835SJames Wright 2918c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_calc_divergence, "Grad F_diff", projection->num_comp * dim, CEED_EVAL_GRAD)); 2928c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_calc_divergence, "qdata", q_data_size, CEED_EVAL_NONE)); 29340b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_calc_divergence, "Div F_diff", diff_flux_proj->num_diff_flux_comps, CEED_EVAL_NONE)); 2948c85b835SJames Wright 2958c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_calc_divergence, NULL, NULL, &op_calc_divergence)); 2968c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "Grad F_diff", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE)); 2978c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 2988c85b835SJames Wright PetscCallCeed( 2998c85b835SJames Wright ceed, CeedOperatorSetField(op_calc_divergence, "Div F_diff", elem_restr_div_diff_flux, CEED_BASIS_NONE, diff_flux_proj->div_diff_flux_ceed)); 3008c85b835SJames Wright 30136038bbcSJames Wright PetscCall( 30236038bbcSJames Wright OperatorApplyContextCreate(projection->dm, NULL, ceed, op_calc_divergence, NULL, NULL, NULL, NULL, &diff_flux_proj->calc_div_diff_flux)); 3038c85b835SJames Wright 3040880fbb6SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_div_diff_flux)); 3058c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_calc_divergence)); 3068c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_calc_divergence)); 3078c85b835SJames Wright } 3088c85b835SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 3098c85b835SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 3108c85b835SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 3118c85b835SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux)); 3128c85b835SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3138c85b835SJames Wright } 3148c85b835SJames Wright 3158c85b835SJames Wright /** 3168c85b835SJames Wright @brief Setup projection of divergence of diffusive flux 3178c85b835SJames Wright 3180c373b74SJames Wright @param[in] honee `Honee` context 3198561fee2SJames Wright @param[in,out] diff_flux_proj Flux projection object to setup 3208c85b835SJames Wright **/ 321e3663b90SJames Wright PetscErrorCode DivDiffFluxProjectionSetup(Honee honee, DivDiffFluxProjectionData diff_flux_proj) { 3228c85b835SJames Wright PetscFunctionBeginUser; 3230c373b74SJames Wright switch (honee->app_ctx->divFdiffproj_method) { 3248c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_DIRECT: 325e3663b90SJames Wright PetscCall(DivDiffFluxProjectionSetup_Direct(honee, diff_flux_proj)); 3268c85b835SJames Wright break; 3278c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_INDIRECT: 328e3663b90SJames Wright PetscCall(DivDiffFluxProjectionSetup_Indirect(honee, diff_flux_proj)); 3298c85b835SJames Wright break; 3308c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_NONE: 3310c373b74SJames Wright SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s", 3320c373b74SJames Wright DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]); 3338c85b835SJames Wright break; 3348c85b835SJames Wright } 3358c85b835SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3368c85b835SJames Wright } 3378c85b835SJames Wright 3388c85b835SJames Wright /** 3398c85b835SJames Wright @brief Project the divergence of diffusive flux 3408c85b835SJames Wright 3418c85b835SJames Wright This implicitly sets the `CeedVector` input (`div_diff_flux_ceed`) to the divergence of diffusive flux. 3428c85b835SJames Wright 3438c85b835SJames Wright @param[in] diff_flux_proj `NodalProjectionData` for the projection 3448c85b835SJames Wright @param[in] Q_loc Localized solution vector 3458c85b835SJames Wright **/ 34636038bbcSJames Wright PetscErrorCode DivDiffFluxProjectionApply(DivDiffFluxProjectionData diff_flux_proj, Vec Q_loc) { 3478c85b835SJames Wright NodalProjectionData projection = diff_flux_proj->projection; 3488c85b835SJames Wright 3498c85b835SJames Wright PetscFunctionBeginUser; 350ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_DivDiffFluxProjection, Q_loc, 0, 0, 0)); 3518c85b835SJames Wright switch (diff_flux_proj->method) { 3528c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_DIRECT: { 353*c8cb0c9bSJames Wright Vec DivDiffFlux, RHS; 3548c85b835SJames Wright 3558c85b835SJames Wright PetscCall(DMGetGlobalVector(projection->dm, &DivDiffFlux)); 356*c8cb0c9bSJames Wright PetscCall(DMGetGlobalVector(projection->dm, &RHS)); 3578c85b835SJames Wright if (diff_flux_proj->ceed_vec_has_array) { 3588c85b835SJames Wright PetscCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->DivDiffFlux_loc)); 3598c85b835SJames Wright diff_flux_proj->ceed_vec_has_array = PETSC_FALSE; 3608c85b835SJames Wright } 361*c8cb0c9bSJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, RHS, projection->l2_rhs_ctx)); 3628c85b835SJames Wright PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_rhs_view")); 3638c85b835SJames Wright 364*c8cb0c9bSJames Wright { 365*c8cb0c9bSJames Wright // Run PCApply manually if using ksp_type preonly -pc_type jacobi 366*c8cb0c9bSJames Wright // This is to avoid an AllReduce call in KSPSolve_Preonly, which causes significant slowdowns for lumped mass matrix solves. 367*c8cb0c9bSJames Wright // See https://gitlab.com/petsc/petsc/-/merge_requests/8048 for more details and a possible fix 368*c8cb0c9bSJames Wright PC pc; 369*c8cb0c9bSJames Wright PetscBool ispreonly, isjacobi; 370*c8cb0c9bSJames Wright PetscCall(KSPGetPC(projection->ksp, &pc)); 371*c8cb0c9bSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)projection->ksp, KSPPREONLY, &ispreonly)); 372*c8cb0c9bSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &isjacobi)); 373*c8cb0c9bSJames Wright if (ispreonly && isjacobi) PetscCall(PCApply(pc, RHS, DivDiffFlux)); 374*c8cb0c9bSJames Wright else PetscCall(KSPSolve(projection->ksp, RHS, DivDiffFlux)); 375*c8cb0c9bSJames Wright } 3768c85b835SJames Wright PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_view")); 3778c85b835SJames Wright 3788c85b835SJames Wright PetscCall(DMGlobalToLocal(projection->dm, DivDiffFlux, INSERT_VALUES, diff_flux_proj->DivDiffFlux_loc)); 3798c85b835SJames Wright PetscCall(VecReadPetscToCeed(diff_flux_proj->DivDiffFlux_loc, &diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->div_diff_flux_ceed)); 3808c85b835SJames Wright diff_flux_proj->ceed_vec_has_array = PETSC_TRUE; 3818c85b835SJames Wright 382*c8cb0c9bSJames Wright PetscCall(DMRestoreGlobalVector(projection->dm, &RHS)); 3838c85b835SJames Wright PetscCall(DMRestoreGlobalVector(projection->dm, &DivDiffFlux)); 3848c85b835SJames Wright break; 3858c85b835SJames Wright } 3868c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_INDIRECT: { 387*c8cb0c9bSJames Wright Vec DiffFlux, RHS; 3888c85b835SJames Wright 3898c85b835SJames Wright PetscCall(DMGetGlobalVector(projection->dm, &DiffFlux)); 390*c8cb0c9bSJames Wright PetscCall(DMGetGlobalVector(projection->dm, &RHS)); 391*c8cb0c9bSJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, RHS, projection->l2_rhs_ctx)); 3928c85b835SJames Wright PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_rhs_view")); 3938c85b835SJames Wright 394*c8cb0c9bSJames Wright { 395*c8cb0c9bSJames Wright // Run PCApply manually if using -ksp_type preonly -pc_type jacobi 396*c8cb0c9bSJames Wright // This is to avoid an AllReduce call in KSPSolve_Preonly, which causes significant slowdowns for lumped mass matrix solves. 397*c8cb0c9bSJames Wright // See https://gitlab.com/petsc/petsc/-/merge_requests/8048 for more details and a possible fix 398*c8cb0c9bSJames Wright PC pc; 399*c8cb0c9bSJames Wright PetscBool ispreonly, isjacobi; 400*c8cb0c9bSJames Wright PetscCall(KSPGetPC(projection->ksp, &pc)); 401*c8cb0c9bSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)projection->ksp, KSPPREONLY, &ispreonly)); 402*c8cb0c9bSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &isjacobi)); 403*c8cb0c9bSJames Wright if (ispreonly && isjacobi) PetscCall(PCApply(pc, RHS, DiffFlux)); 404*c8cb0c9bSJames Wright else PetscCall(KSPSolve(projection->ksp, RHS, DiffFlux)); 405*c8cb0c9bSJames Wright } 4068c85b835SJames Wright PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_view")); 4078c85b835SJames Wright 4088c85b835SJames Wright PetscCall(ApplyCeedOperatorGlobalToLocal(DiffFlux, NULL, diff_flux_proj->calc_div_diff_flux)); 409*c8cb0c9bSJames Wright PetscCall(DMRestoreGlobalVector(projection->dm, &RHS)); 4108c85b835SJames Wright PetscCall(DMRestoreGlobalVector(projection->dm, &DiffFlux)); 4118c85b835SJames Wright } break; 4128c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_NONE: 4138c85b835SJames Wright SETERRQ(PetscObjectComm((PetscObject)projection->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s", 4148c85b835SJames Wright DivDiffFluxProjectionMethods[diff_flux_proj->method]); 4158c85b835SJames Wright break; 4168c85b835SJames Wright } 417ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_DivDiffFluxProjection, Q_loc, 0, 0, 0)); 4188c85b835SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 4198c85b835SJames Wright } 4208c85b835SJames Wright 4218c85b835SJames Wright /** 4228c85b835SJames Wright @brief Destroy `DivDiffFluxProjectionData` object 4238c85b835SJames Wright 4248c85b835SJames Wright @param[in,out] diff_flux_proj Object to destroy 4258c85b835SJames Wright **/ 4268c85b835SJames Wright PetscErrorCode DivDiffFluxProjectionDataDestroy(DivDiffFluxProjectionData diff_flux_proj) { 4278c85b835SJames Wright PetscFunctionBeginUser; 4288c85b835SJames Wright if (diff_flux_proj == NULL) PetscFunctionReturn(PETSC_SUCCESS); 4290880fbb6SJames Wright Ceed ceed = CeedVectorReturnCeed(diff_flux_proj->div_diff_flux_ceed); 4300880fbb6SJames Wright 4314ea616f4SJames Wright PetscCall(NodalProjectionDataDestroy(&diff_flux_proj->projection)); 4328c85b835SJames Wright PetscCall(OperatorApplyContextDestroy(diff_flux_proj->calc_div_diff_flux)); 4338c85b835SJames Wright if (diff_flux_proj->ceed_vec_has_array) { 4348c85b835SJames Wright PetscCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->DivDiffFlux_loc)); 4358c85b835SJames Wright diff_flux_proj->ceed_vec_has_array = PETSC_FALSE; 4368c85b835SJames Wright } 4370880fbb6SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&diff_flux_proj->div_diff_flux_ceed)); 4380880fbb6SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&diff_flux_proj->elem_restr_div_diff_flux)); 4390880fbb6SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&diff_flux_proj->basis_div_diff_flux)); 4408c85b835SJames Wright PetscCall(VecDestroy(&diff_flux_proj->DivDiffFlux_loc)); 4418c85b835SJames Wright PetscCall(PetscFree(diff_flux_proj)); 4428c85b835SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 4438c85b835SJames Wright } 444