xref: /honee/src/diff_flux_projection.c (revision 6b9fd99336306b4ee2da369fb0ba126183859527)
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));
55*6b9fd993SJames Wright       PetscCall(DMPlexCeedBasisCreate(honee->ceed, projection->dm, domain_label, label_value, height, dm_field,
56*6b9fd993SJames Wright                                       &diff_flux_proj_->basis_div_diff_flux));
57d9e69621SJames Wright       diff_flux_proj_->eval_mode_div_diff_flux = CEED_EVAL_INTERP;
581af555e8SJames Wright 
591af555e8SJames Wright       {  // Create face labels on projection->dm for boundary integrals
601af555e8SJames Wright         DMLabel  face_sets_label;
611af555e8SJames Wright         PetscInt num_face_set_values, *face_set_values;
621af555e8SJames Wright 
630c373b74SJames Wright         PetscCall(DMGetLabel(honee->dm, "Face Sets", &face_sets_label));
640c373b74SJames Wright         PetscCall(DMLabelCreateGlobalValueArray(honee->dm, face_sets_label, &num_face_set_values, &face_set_values));
651af555e8SJames Wright         for (PetscInt f = 0; f < num_face_set_values; f++) {
661af555e8SJames Wright           DMLabel face_orientation_label;
671af555e8SJames Wright           char   *face_orientation_label_name;
681af555e8SJames Wright 
690c373b74SJames Wright           PetscCall(DMPlexCreateFaceLabel(honee->dm, face_set_values[f], &face_orientation_label_name));
700c373b74SJames Wright           PetscCall(DMGetLabel(honee->dm, face_orientation_label_name, &face_orientation_label));
711af555e8SJames Wright           PetscCall(DMAddLabel(projection->dm, face_orientation_label));
721af555e8SJames Wright           PetscCall(PetscFree(face_orientation_label_name));
731af555e8SJames Wright         }
7436038bbcSJames Wright         PetscCall(PetscFree(face_set_values));
751af555e8SJames Wright       }
768c85b835SJames Wright     } break;
778c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_INDIRECT: {
78d9e69621SJames Wright       projection->num_comp = diff_flux_proj_->num_diff_flux_comps * dim;
798c85b835SJames Wright       PetscCall(PetscObjectSetName((PetscObject)projection->dm, "DiffFluxProj"));
8036038bbcSJames Wright       PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &projection->num_comp, projection->dm));
818c85b835SJames Wright 
820c373b74SJames Wright       PetscCall(DMPlexCeedElemRestrictionQDataCreate(honee->ceed, projection->dm, domain_label, label_value, height,
83d9e69621SJames Wright                                                      diff_flux_proj_->num_diff_flux_comps, &diff_flux_proj_->elem_restr_div_diff_flux));
840c373b74SJames Wright       PetscCallCeed(honee->ceed,
85d9e69621SJames Wright                     CeedElemRestrictionCreateVector(diff_flux_proj_->elem_restr_div_diff_flux, &diff_flux_proj_->div_diff_flux_ceed, NULL));
86d9e69621SJames Wright       diff_flux_proj_->basis_div_diff_flux     = CEED_BASIS_NONE;
87d9e69621SJames Wright       diff_flux_proj_->eval_mode_div_diff_flux = CEED_EVAL_NONE;
888c85b835SJames Wright     } break;
898c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_NONE:
90d9e69621SJames Wright       SETERRQ(honee->comm, PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
91d9e69621SJames Wright               DivDiffFluxProjectionMethods[divFdiffproj_method]);
928c85b835SJames Wright       break;
938c85b835SJames Wright   }
94d9e69621SJames Wright   *diff_flux_proj = diff_flux_proj_;
958c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
968c85b835SJames Wright };
978c85b835SJames Wright 
988c85b835SJames Wright /**
990880fbb6SJames Wright   @brief Return the objects required for the Divergence of Diffusive flux to be read by a `CeedOperator`
1000880fbb6SJames Wright 
1010880fbb6SJames Wright   @param[in]  diff_flux_proj Projection object
1020880fbb6SJames Wright   @param[out] elem_restr     Element restriction for the divergence of diffusive flux, or `NULL`
1030880fbb6SJames Wright   @param[out] basis          Basis for the divergence of diffusive flux, or `NULL`
1040880fbb6SJames Wright   @param[out] vector         Vector for the divergence of diffusive flux, or `NULL`
1050880fbb6SJames Wright   @param[out] eval_mode      Eval mode for the divergence of diffusive flux, or `NULL`
1060880fbb6SJames Wright **/
1070880fbb6SJames Wright PetscErrorCode DivDiffFluxProjectionGetOperatorFieldData(DivDiffFluxProjectionData diff_flux_proj, CeedElemRestriction *elem_restr, CeedBasis *basis,
1080880fbb6SJames Wright                                                          CeedVector *vector, CeedEvalMode *eval_mode) {
1090880fbb6SJames Wright   Ceed ceed = CeedVectorReturnCeed(diff_flux_proj->div_diff_flux_ceed);
1100880fbb6SJames Wright 
1110880fbb6SJames Wright   PetscFunctionBeginUser;
1120880fbb6SJames Wright   if (elem_restr) PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(diff_flux_proj->elem_restr_div_diff_flux, elem_restr));
1130880fbb6SJames Wright   if (basis) PetscCallCeed(ceed, CeedBasisReferenceCopy(diff_flux_proj->basis_div_diff_flux, basis));
1140880fbb6SJames Wright   if (vector) PetscCallCeed(ceed, CeedVectorReferenceCopy(diff_flux_proj->div_diff_flux_ceed, vector));
1150880fbb6SJames Wright   if (eval_mode) *eval_mode = diff_flux_proj->eval_mode_div_diff_flux;
1160880fbb6SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1170880fbb6SJames Wright }
1180880fbb6SJames Wright 
1190880fbb6SJames Wright /**
1208c85b835SJames Wright   @brief Setup direct projection of divergence of diffusive flux
1218c85b835SJames Wright 
1220c373b74SJames Wright   @param[in]     honee          `Honee` context
1238561fee2SJames Wright   @param[in,out] diff_flux_proj Flux projection object to setup
1248c85b835SJames Wright **/
125e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionSetup_Direct(Honee honee, DivDiffFluxProjectionData diff_flux_proj) {
1260c373b74SJames Wright   Ceed                ceed       = honee->ceed;
1278c85b835SJames Wright   NodalProjectionData projection = diff_flux_proj->projection;
12836038bbcSJames Wright   MPI_Comm            comm       = PetscObjectComm((PetscObject)projection->dm);
1298c85b835SJames Wright 
1308c85b835SJames Wright   PetscFunctionBeginUser;
13136038bbcSJames Wright   {  // Create Projection RHS OperatorApplyContext
13236038bbcSJames Wright     CeedOperator op_rhs;
1338c85b835SJames Wright 
13436038bbcSJames Wright     PetscCheck(diff_flux_proj->CreateRHSOperator_Direct, comm, PETSC_ERR_ARG_WRONGSTATE,
13536038bbcSJames Wright                "Must define CreateRHSOperator_Direct to use indirect div_diff_flux projection");
136e3663b90SJames Wright     PetscCall(diff_flux_proj->CreateRHSOperator_Direct(honee, diff_flux_proj, &op_rhs));
1378c85b835SJames Wright     PetscCall(DMCreateLocalVector(projection->dm, &diff_flux_proj->DivDiffFlux_loc));
1388c85b835SJames Wright     diff_flux_proj->ceed_vec_has_array = PETSC_FALSE;
1390c373b74SJames Wright     PetscCall(OperatorApplyContextCreate(honee->dm, projection->dm, ceed, op_rhs, NULL, NULL, NULL, diff_flux_proj->DivDiffFlux_loc,
14036038bbcSJames Wright                                          &projection->l2_rhs_ctx));
14136038bbcSJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
14236038bbcSJames Wright   }
1438c85b835SJames Wright 
1448c85b835SJames Wright   {  // -- Build Mass operator
1458c85b835SJames Wright     CeedQFunction       qf_mass;
1468c85b835SJames Wright     CeedOperator        op_mass;
147590d9cddSJames Wright     CeedBasis           basis_div_diff_flux             = NULL;
148590d9cddSJames Wright     CeedElemRestriction elem_restr_div_diff_flux_volume = NULL, elem_restr_qd;
1490880fbb6SJames Wright     CeedVector          q_data;
1500880fbb6SJames Wright     CeedInt             q_data_size;
1510880fbb6SJames Wright     PetscInt            label_value  = 0;
1520880fbb6SJames Wright     DMLabel             domain_label = NULL;
1530880fbb6SJames Wright 
154590d9cddSJames Wright     PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_div_diff_flux_volume, &basis_div_diff_flux, NULL, NULL));
155e3663b90SJames 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,
156e3663b90SJames Wright                        &q_data_size));
1578c85b835SJames Wright 
15864dd23feSJames Wright     PetscCall(HoneeMassQFunctionCreate(ceed, projection->num_comp, q_data_size, &qf_mass));
1598c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
160590d9cddSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_div_diff_flux_volume, basis_div_diff_flux, CEED_VECTOR_ACTIVE));
1618c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
162590d9cddSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_div_diff_flux_volume, basis_div_diff_flux, CEED_VECTOR_ACTIVE));
1638c85b835SJames Wright 
1648c85b835SJames Wright     {  // -- Setup KSP for L^2 projection
1658c85b835SJames Wright       Mat mat_mass;
1668c85b835SJames Wright 
1678c85b835SJames Wright       PetscCall(MatCreateCeed(projection->dm, projection->dm, op_mass, NULL, &mat_mass));
1688c85b835SJames Wright 
1698c85b835SJames Wright       PetscCall(KSPCreate(comm, &projection->ksp));
1708c85b835SJames Wright       PetscCall(KSPSetOptionsPrefix(projection->ksp, "div_diff_flux_projection_"));
1718c85b835SJames Wright       {  // lumped by default
1728c85b835SJames Wright         PC pc;
1738c85b835SJames Wright         PetscCall(KSPGetPC(projection->ksp, &pc));
1748c85b835SJames Wright         PetscCall(PCSetType(pc, PCJACOBI));
1758c85b835SJames Wright         PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
1768c85b835SJames Wright         PetscCall(KSPSetType(projection->ksp, KSPPREONLY));
1778c85b835SJames Wright       }
1788c85b835SJames Wright       PetscCall(KSPSetFromOptions_WithMatCeed(projection->ksp, mat_mass));
1798c85b835SJames Wright       PetscCall(MatDestroy(&mat_mass));
1808c85b835SJames Wright     }
1810880fbb6SJames Wright 
182590d9cddSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_div_diff_flux_volume));
183590d9cddSJames Wright     PetscCallCeed(ceed, CeedBasisDestroy(&basis_div_diff_flux));
1840880fbb6SJames Wright     PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
1850880fbb6SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
1868c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
1878c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
1888c85b835SJames Wright   }
1898c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1908c85b835SJames Wright }
1918c85b835SJames Wright 
1928c85b835SJames Wright /**
1938c85b835SJames Wright   @brief Setup indirect projection of divergence of diffusive flux
1948c85b835SJames Wright 
1950c373b74SJames Wright   @param[in]     honee          `Honee` context
1968561fee2SJames Wright   @param[in,out] diff_flux_proj Flux projection object to setup
1978c85b835SJames Wright **/
198e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionSetup_Indirect(Honee honee, DivDiffFluxProjectionData diff_flux_proj) {
1990c373b74SJames Wright   Ceed                ceed       = honee->ceed;
2008c85b835SJames Wright   NodalProjectionData projection = diff_flux_proj->projection;
2018c85b835SJames Wright   CeedBasis           basis_diff_flux;
2028c85b835SJames Wright   CeedElemRestriction elem_restr_diff_flux, elem_restr_qd;
2038c85b835SJames Wright   CeedVector          q_data;
2040880fbb6SJames Wright   CeedInt             q_data_size;
20536038bbcSJames Wright   MPI_Comm            comm = PetscObjectComm((PetscObject)projection->dm);
2068c85b835SJames Wright 
2078c85b835SJames Wright   PetscFunctionBeginUser;
2080880fbb6SJames Wright   {
2090880fbb6SJames Wright     PetscInt label_value = 0, height = 0, dm_field = 0;
2100880fbb6SJames Wright     DMLabel  domain_label = NULL;
2118c85b835SJames Wright 
2128c85b835SJames Wright     PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, domain_label, label_value, height, dm_field, &elem_restr_diff_flux));
213*6b9fd993SJames Wright     PetscCall(DMPlexCeedBasisCreate(ceed, projection->dm, domain_label, label_value, height, dm_field, &basis_diff_flux));
214e3663b90SJames 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,
215e3663b90SJames Wright                        &q_data_size));
2160880fbb6SJames Wright   }
2178c85b835SJames Wright 
21836038bbcSJames Wright   {
2198c85b835SJames Wright     CeedOperator op_rhs;
2208c85b835SJames Wright 
22136038bbcSJames Wright     PetscCheck(diff_flux_proj->CreateRHSOperator_Indirect, comm, PETSC_ERR_ARG_WRONGSTATE,
22236038bbcSJames Wright                "Must define CreateRHSOperator_Indirect to use indirect div_diff_flux projection");
223e3663b90SJames Wright     PetscCall(diff_flux_proj->CreateRHSOperator_Indirect(honee, diff_flux_proj, &op_rhs));
2240c373b74SJames Wright     PetscCall(OperatorApplyContextCreate(honee->dm, projection->dm, ceed, op_rhs, NULL, NULL, NULL, NULL, &projection->l2_rhs_ctx));
2258c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
2268c85b835SJames Wright   }
2278c85b835SJames Wright 
2288c85b835SJames Wright   {  // -- Build Mass operator
2298c85b835SJames Wright     CeedQFunction qf_mass;
2308c85b835SJames Wright     CeedOperator  op_mass;
2318c85b835SJames Wright 
23264dd23feSJames Wright     PetscCall(HoneeMassQFunctionCreate(ceed, projection->num_comp, q_data_size, &qf_mass));
2338c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
2348c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
2358c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
2368c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
2378c85b835SJames Wright 
2388c85b835SJames Wright     {  // -- Setup KSP for L^2 projection
2398c85b835SJames Wright       Mat mat_mass;
2408c85b835SJames Wright 
2418c85b835SJames Wright       PetscCall(MatCreateCeed(projection->dm, projection->dm, op_mass, NULL, &mat_mass));
2428c85b835SJames Wright 
2438c85b835SJames Wright       PetscCall(KSPCreate(comm, &projection->ksp));
2448c85b835SJames Wright       PetscCall(KSPSetOptionsPrefix(projection->ksp, "div_diff_flux_projection_"));
2458c85b835SJames Wright       {  // lumped by default
2468c85b835SJames Wright         PC pc;
2478c85b835SJames Wright         PetscCall(KSPGetPC(projection->ksp, &pc));
2488c85b835SJames Wright         PetscCall(PCSetType(pc, PCJACOBI));
2498c85b835SJames Wright         PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
2508c85b835SJames Wright         PetscCall(KSPSetType(projection->ksp, KSPPREONLY));
2518c85b835SJames Wright       }
2528c85b835SJames Wright       PetscCall(KSPSetFromOptions_WithMatCeed(projection->ksp, mat_mass));
2538c85b835SJames Wright       PetscCall(MatDestroy(&mat_mass));
2548c85b835SJames Wright     }
2558c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
2568c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
2578c85b835SJames Wright   }
2588c85b835SJames Wright 
2598c85b835SJames Wright   {  // Create OperatorApplyContext to calculate divergence at quadrature points
26040b78511SJames Wright     CeedQFunction       qf_calc_divergence = NULL;
2618c85b835SJames Wright     CeedOperator        op_calc_divergence;
2620880fbb6SJames Wright     CeedElemRestriction elem_restr_div_diff_flux = NULL;
2630880fbb6SJames Wright     PetscInt            dim;
2648c85b835SJames Wright 
2650880fbb6SJames Wright     PetscCall(DMGetDimension(projection->dm, &dim));
2660880fbb6SJames Wright     PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_div_diff_flux, NULL, NULL, NULL));
2678c85b835SJames Wright 
26840b78511SJames Wright     switch (dim) {
26940b78511SJames Wright       case 2:
27040b78511SJames Wright         switch (diff_flux_proj->num_diff_flux_comps) {
27140b78511SJames Wright           case 1:
27240b78511SJames Wright             PetscCallCeed(ceed,
27340b78511SJames Wright                           CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux2D_1, ComputeDivDiffusiveFlux2D_1_loc, &qf_calc_divergence));
27440b78511SJames Wright             break;
27540b78511SJames Wright         }
27640b78511SJames Wright         break;
27740b78511SJames Wright       case 3:
27840b78511SJames Wright         switch (diff_flux_proj->num_diff_flux_comps) {
27940b78511SJames Wright           case 1:
28040b78511SJames Wright             PetscCallCeed(ceed,
28140b78511SJames Wright                           CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux3D_1, ComputeDivDiffusiveFlux3D_1_loc, &qf_calc_divergence));
28240b78511SJames Wright             break;
28340b78511SJames Wright           case 4:
28440b78511SJames Wright             PetscCallCeed(ceed,
28540b78511SJames Wright                           CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux3D_4, ComputeDivDiffusiveFlux3D_4_loc, &qf_calc_divergence));
28640b78511SJames Wright             break;
28740b78511SJames Wright         }
28840b78511SJames Wright         break;
28940b78511SJames Wright     }
29040b78511SJames Wright     PetscCheck(qf_calc_divergence, comm, PETSC_ERR_SUP,
29140b78511SJames Wright                "QFunction for calculating divergence of diffusive flux does not exist for"
29240b78511SJames Wright                " %" PetscInt_FMT " dimensional grid and %" PetscInt_FMT
29340b78511SJames Wright                " number of components.\nA new qfunction can be easily added; see source code for pattern.",
29440b78511SJames Wright                dim, diff_flux_proj->num_diff_flux_comps);
2958c85b835SJames Wright 
2968c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_calc_divergence, "Grad F_diff", projection->num_comp * dim, CEED_EVAL_GRAD));
2978c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_calc_divergence, "qdata", q_data_size, CEED_EVAL_NONE));
29840b78511SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_calc_divergence, "Div F_diff", diff_flux_proj->num_diff_flux_comps, CEED_EVAL_NONE));
2998c85b835SJames Wright 
3008c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_calc_divergence, NULL, NULL, &op_calc_divergence));
3018c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "Grad F_diff", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
3028c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
3039eadbee4SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "Div F_diff", elem_restr_div_diff_flux, CEED_BASIS_NONE,
3049eadbee4SJames Wright                                              diff_flux_proj->div_diff_flux_ceed));
3058c85b835SJames Wright 
30665fa3c0aSJames Wright     PetscCall(OperatorApplyContextCreate(projection->dm, NULL, ceed, op_calc_divergence, NULL, CEED_VECTOR_NONE, NULL, NULL,
30765fa3c0aSJames Wright                                          &diff_flux_proj->calc_div_diff_flux));
3088c85b835SJames Wright 
3090880fbb6SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_div_diff_flux));
3108c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_calc_divergence));
3118c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_calc_divergence));
3128c85b835SJames Wright   }
3138c85b835SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
3148c85b835SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
3158c85b835SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
3168c85b835SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux));
3178c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3188c85b835SJames Wright }
3198c85b835SJames Wright 
3208c85b835SJames Wright /**
3218c85b835SJames Wright   @brief Setup projection of divergence of diffusive flux
3228c85b835SJames Wright 
3230c373b74SJames Wright   @param[in]     honee          `Honee` context
3248561fee2SJames Wright   @param[in,out] diff_flux_proj Flux projection object to setup
3258c85b835SJames Wright **/
326e3663b90SJames Wright PetscErrorCode DivDiffFluxProjectionSetup(Honee honee, DivDiffFluxProjectionData diff_flux_proj) {
3278c85b835SJames Wright   PetscFunctionBeginUser;
3280c373b74SJames Wright   switch (honee->app_ctx->divFdiffproj_method) {
3298c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_DIRECT:
330e3663b90SJames Wright       PetscCall(DivDiffFluxProjectionSetup_Direct(honee, diff_flux_proj));
3318c85b835SJames Wright       break;
3328c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_INDIRECT:
333e3663b90SJames Wright       PetscCall(DivDiffFluxProjectionSetup_Indirect(honee, diff_flux_proj));
3348c85b835SJames Wright       break;
3358c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_NONE:
3360c373b74SJames Wright       SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
3370c373b74SJames Wright               DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]);
3388c85b835SJames Wright       break;
3398c85b835SJames Wright   }
3408c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3418c85b835SJames Wright }
3428c85b835SJames Wright 
3438c85b835SJames Wright /**
3448c85b835SJames Wright   @brief Project the divergence of diffusive flux
3458c85b835SJames Wright 
3468c85b835SJames Wright   This implicitly sets the `CeedVector` input (`div_diff_flux_ceed`) to the divergence of diffusive flux.
3478c85b835SJames Wright 
3488c85b835SJames Wright   @param[in]  diff_flux_proj `NodalProjectionData` for the projection
3498c85b835SJames Wright   @param[in]  Q_loc          Localized solution vector
3508c85b835SJames Wright **/
35136038bbcSJames Wright PetscErrorCode DivDiffFluxProjectionApply(DivDiffFluxProjectionData diff_flux_proj, Vec Q_loc) {
3528c85b835SJames Wright   NodalProjectionData projection = diff_flux_proj->projection;
3538c85b835SJames Wright 
3548c85b835SJames Wright   PetscFunctionBeginUser;
355ea615d4cSJames Wright   PetscCall(PetscLogEventBegin(HONEE_DivDiffFluxProjection, Q_loc, 0, 0, 0));
3568c85b835SJames Wright   switch (diff_flux_proj->method) {
3578c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_DIRECT: {
358c8cb0c9bSJames Wright       Vec DivDiffFlux, RHS;
3598c85b835SJames Wright 
3608c85b835SJames Wright       PetscCall(DMGetGlobalVector(projection->dm, &DivDiffFlux));
361c8cb0c9bSJames Wright       PetscCall(DMGetGlobalVector(projection->dm, &RHS));
3628c85b835SJames Wright       if (diff_flux_proj->ceed_vec_has_array) {
3638c85b835SJames Wright         PetscCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->DivDiffFlux_loc));
3648c85b835SJames Wright         diff_flux_proj->ceed_vec_has_array = PETSC_FALSE;
3658c85b835SJames Wright       }
366c8cb0c9bSJames Wright       PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, RHS, projection->l2_rhs_ctx));
3678c85b835SJames Wright       PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_rhs_view"));
3688c85b835SJames Wright 
369c8cb0c9bSJames Wright       {
370c8cb0c9bSJames Wright         // Run PCApply manually if using ksp_type preonly -pc_type jacobi
371c8cb0c9bSJames Wright         // This is to avoid an AllReduce call in KSPSolve_Preonly, which causes significant slowdowns for lumped mass matrix solves.
372c8cb0c9bSJames Wright         // See https://gitlab.com/petsc/petsc/-/merge_requests/8048 for more details and a possible fix
373c8cb0c9bSJames Wright         PC        pc;
374c8cb0c9bSJames Wright         PetscBool ispreonly, isjacobi;
375c8cb0c9bSJames Wright         PetscCall(KSPGetPC(projection->ksp, &pc));
376c8cb0c9bSJames Wright         PetscCall(PetscObjectTypeCompare((PetscObject)projection->ksp, KSPPREONLY, &ispreonly));
377c8cb0c9bSJames Wright         PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &isjacobi));
378c8cb0c9bSJames Wright         if (ispreonly && isjacobi) PetscCall(PCApply(pc, RHS, DivDiffFlux));
379c8cb0c9bSJames Wright         else PetscCall(KSPSolve(projection->ksp, RHS, DivDiffFlux));
380c8cb0c9bSJames Wright       }
3818c85b835SJames Wright       PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_view"));
3828c85b835SJames Wright 
3838c85b835SJames Wright       PetscCall(DMGlobalToLocal(projection->dm, DivDiffFlux, INSERT_VALUES, diff_flux_proj->DivDiffFlux_loc));
3848c85b835SJames Wright       PetscCall(VecReadPetscToCeed(diff_flux_proj->DivDiffFlux_loc, &diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->div_diff_flux_ceed));
3858c85b835SJames Wright       diff_flux_proj->ceed_vec_has_array = PETSC_TRUE;
3868c85b835SJames Wright 
387c8cb0c9bSJames Wright       PetscCall(DMRestoreGlobalVector(projection->dm, &RHS));
3888c85b835SJames Wright       PetscCall(DMRestoreGlobalVector(projection->dm, &DivDiffFlux));
3898c85b835SJames Wright       break;
3908c85b835SJames Wright     }
3918c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_INDIRECT: {
392c8cb0c9bSJames Wright       Vec DiffFlux, RHS;
3938c85b835SJames Wright 
3948c85b835SJames Wright       PetscCall(DMGetGlobalVector(projection->dm, &DiffFlux));
395c8cb0c9bSJames Wright       PetscCall(DMGetGlobalVector(projection->dm, &RHS));
396c8cb0c9bSJames Wright       PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, RHS, projection->l2_rhs_ctx));
3978c85b835SJames Wright       PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_rhs_view"));
3988c85b835SJames Wright 
399c8cb0c9bSJames Wright       {
400c8cb0c9bSJames Wright         // Run PCApply manually if using -ksp_type preonly -pc_type jacobi
401c8cb0c9bSJames Wright         // This is to avoid an AllReduce call in KSPSolve_Preonly, which causes significant slowdowns for lumped mass matrix solves.
402c8cb0c9bSJames Wright         // See https://gitlab.com/petsc/petsc/-/merge_requests/8048 for more details and a possible fix
403c8cb0c9bSJames Wright         PC        pc;
404c8cb0c9bSJames Wright         PetscBool ispreonly, isjacobi;
405c8cb0c9bSJames Wright         PetscCall(KSPGetPC(projection->ksp, &pc));
406c8cb0c9bSJames Wright         PetscCall(PetscObjectTypeCompare((PetscObject)projection->ksp, KSPPREONLY, &ispreonly));
407c8cb0c9bSJames Wright         PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &isjacobi));
408c8cb0c9bSJames Wright         if (ispreonly && isjacobi) PetscCall(PCApply(pc, RHS, DiffFlux));
409c8cb0c9bSJames Wright         else PetscCall(KSPSolve(projection->ksp, RHS, DiffFlux));
410c8cb0c9bSJames Wright       }
4118c85b835SJames Wright       PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_view"));
4128c85b835SJames Wright 
4138c85b835SJames Wright       PetscCall(ApplyCeedOperatorGlobalToLocal(DiffFlux, NULL, diff_flux_proj->calc_div_diff_flux));
414c8cb0c9bSJames Wright       PetscCall(DMRestoreGlobalVector(projection->dm, &RHS));
4158c85b835SJames Wright       PetscCall(DMRestoreGlobalVector(projection->dm, &DiffFlux));
4168c85b835SJames Wright     } break;
4178c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_NONE:
4188c85b835SJames Wright       SETERRQ(PetscObjectComm((PetscObject)projection->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
4198c85b835SJames Wright               DivDiffFluxProjectionMethods[diff_flux_proj->method]);
4208c85b835SJames Wright       break;
4218c85b835SJames Wright   }
422ea615d4cSJames Wright   PetscCall(PetscLogEventEnd(HONEE_DivDiffFluxProjection, Q_loc, 0, 0, 0));
4238c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4248c85b835SJames Wright }
4258c85b835SJames Wright 
4268c85b835SJames Wright /**
4278c85b835SJames Wright   @brief Destroy `DivDiffFluxProjectionData` object
4288c85b835SJames Wright 
4298c85b835SJames Wright   @param[in,out] diff_flux_proj Object to destroy
4308c85b835SJames Wright **/
4318c85b835SJames Wright PetscErrorCode DivDiffFluxProjectionDataDestroy(DivDiffFluxProjectionData diff_flux_proj) {
4328c85b835SJames Wright   PetscFunctionBeginUser;
4338c85b835SJames Wright   if (diff_flux_proj == NULL) PetscFunctionReturn(PETSC_SUCCESS);
4340880fbb6SJames Wright   Ceed ceed = CeedVectorReturnCeed(diff_flux_proj->div_diff_flux_ceed);
4350880fbb6SJames Wright 
4364ea616f4SJames Wright   PetscCall(NodalProjectionDataDestroy(&diff_flux_proj->projection));
4378c85b835SJames Wright   PetscCall(OperatorApplyContextDestroy(diff_flux_proj->calc_div_diff_flux));
4388c85b835SJames Wright   if (diff_flux_proj->ceed_vec_has_array) {
4398c85b835SJames Wright     PetscCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->DivDiffFlux_loc));
4408c85b835SJames Wright     diff_flux_proj->ceed_vec_has_array = PETSC_FALSE;
4418c85b835SJames Wright   }
4420880fbb6SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&diff_flux_proj->div_diff_flux_ceed));
4430880fbb6SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&diff_flux_proj->elem_restr_div_diff_flux));
4440880fbb6SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&diff_flux_proj->basis_div_diff_flux));
4458c85b835SJames Wright   PetscCall(VecDestroy(&diff_flux_proj->DivDiffFlux_loc));
4468c85b835SJames Wright   PetscCall(PetscFree(diff_flux_proj));
4478c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4488c85b835SJames Wright }
449