xref: /honee/src/diff_flux_projection.c (revision d9e6962164416f0bf8fc41cae34da27fa35469f8)
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
18*d9e69621SJames 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)
20*d9e69621SJames Wright   @param[out] diff_flux_proj      The `DivDiffFluxProjectionData` object created, or set to `NULL` if `divFdiffproj_method = DIV_DIFF_FLUX_PROJ_NONE`
218c85b835SJames Wright **/
22*d9e69621SJames Wright PetscErrorCode DivDiffFluxProjectionCreate(Honee honee, DivDiffFluxProjectionMethod divFdiffproj_method, PetscInt num_diff_flux_comps,
23*d9e69621SJames 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;
26*d9e69621SJames Wright   DivDiffFluxProjectionData diff_flux_proj_;
278c85b835SJames Wright   NodalProjectionData       projection;
288c85b835SJames Wright 
298c85b835SJames Wright   PetscFunctionBeginUser;
30*d9e69621SJames Wright   if (divFdiffproj_method == DIV_DIFF_FLUX_PROJ_NONE) {
31*d9e69621SJames Wright     *diff_flux_proj = NULL;
3236038bbcSJames Wright     PetscFunctionReturn(PETSC_SUCCESS);
3336038bbcSJames Wright   }
34*d9e69621SJames Wright   PetscCall(PetscNew(&diff_flux_proj_));
35*d9e69621SJames Wright   PetscCall(PetscNew(&diff_flux_proj_->projection));
36*d9e69621SJames Wright   projection                           = diff_flux_proj_->projection;
37*d9e69621SJames Wright   diff_flux_proj_->method              = divFdiffproj_method;
38*d9e69621SJames Wright   diff_flux_proj_->num_diff_flux_comps = num_diff_flux_comps;
398c85b835SJames Wright 
400c373b74SJames Wright   PetscCall(DMClone(honee->dm, &projection->dm));
418c85b835SJames Wright   PetscCall(DMSetMatrixPreallocateSkip(projection->dm, PETSC_TRUE));
428c85b835SJames Wright   PetscCall(DMGetDimension(projection->dm, &dim));
43*d9e69621SJames Wright   switch (diff_flux_proj_->method) {
448c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_DIRECT: {
45*d9e69621SJames Wright       projection->num_comp = diff_flux_proj_->num_diff_flux_comps;
468c85b835SJames Wright       PetscCall(PetscObjectSetName((PetscObject)projection->dm, "DivDiffFluxProj"));
4736038bbcSJames Wright       PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &projection->num_comp, projection->dm));
488c85b835SJames Wright 
490c373b74SJames Wright       PetscCall(DMPlexCeedElemRestrictionCreate(honee->ceed, projection->dm, domain_label, label_value, height, dm_field,
50*d9e69621SJames Wright                                                 &diff_flux_proj_->elem_restr_div_diff_flux));
510c373b74SJames Wright       PetscCallCeed(honee->ceed,
52*d9e69621SJames Wright                     CeedElemRestrictionCreateVector(diff_flux_proj_->elem_restr_div_diff_flux, &diff_flux_proj_->div_diff_flux_ceed, NULL));
53*d9e69621SJames Wright       PetscCall(CreateBasisFromPlex(honee->ceed, projection->dm, domain_label, label_value, height, dm_field, &diff_flux_proj_->basis_div_diff_flux));
54*d9e69621SJames Wright       diff_flux_proj_->eval_mode_div_diff_flux = CEED_EVAL_INTERP;
551af555e8SJames Wright 
561af555e8SJames Wright       {  // Create face labels on projection->dm for boundary integrals
571af555e8SJames Wright         DMLabel  face_sets_label;
581af555e8SJames Wright         PetscInt num_face_set_values, *face_set_values;
591af555e8SJames Wright 
600c373b74SJames Wright         PetscCall(DMGetLabel(honee->dm, "Face Sets", &face_sets_label));
610c373b74SJames Wright         PetscCall(DMLabelCreateGlobalValueArray(honee->dm, face_sets_label, &num_face_set_values, &face_set_values));
621af555e8SJames Wright         for (PetscInt f = 0; f < num_face_set_values; f++) {
631af555e8SJames Wright           DMLabel face_orientation_label;
641af555e8SJames Wright           char   *face_orientation_label_name;
651af555e8SJames Wright 
660c373b74SJames Wright           PetscCall(DMPlexCreateFaceLabel(honee->dm, face_set_values[f], &face_orientation_label_name));
670c373b74SJames Wright           PetscCall(DMGetLabel(honee->dm, face_orientation_label_name, &face_orientation_label));
681af555e8SJames Wright           PetscCall(DMAddLabel(projection->dm, face_orientation_label));
691af555e8SJames Wright           PetscCall(PetscFree(face_orientation_label_name));
701af555e8SJames Wright         }
7136038bbcSJames Wright         PetscCall(PetscFree(face_set_values));
721af555e8SJames Wright       }
738c85b835SJames Wright     } break;
748c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_INDIRECT: {
75*d9e69621SJames Wright       projection->num_comp = diff_flux_proj_->num_diff_flux_comps * dim;
768c85b835SJames Wright       PetscCall(PetscObjectSetName((PetscObject)projection->dm, "DiffFluxProj"));
7736038bbcSJames Wright       PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &projection->num_comp, projection->dm));
788c85b835SJames Wright 
790c373b74SJames Wright       PetscCall(DMPlexCeedElemRestrictionQDataCreate(honee->ceed, projection->dm, domain_label, label_value, height,
80*d9e69621SJames Wright                                                      diff_flux_proj_->num_diff_flux_comps, &diff_flux_proj_->elem_restr_div_diff_flux));
810c373b74SJames Wright       PetscCallCeed(honee->ceed,
82*d9e69621SJames Wright                     CeedElemRestrictionCreateVector(diff_flux_proj_->elem_restr_div_diff_flux, &diff_flux_proj_->div_diff_flux_ceed, NULL));
83*d9e69621SJames Wright       diff_flux_proj_->basis_div_diff_flux     = CEED_BASIS_NONE;
84*d9e69621SJames Wright       diff_flux_proj_->eval_mode_div_diff_flux = CEED_EVAL_NONE;
858c85b835SJames Wright     } break;
868c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_NONE:
87*d9e69621SJames Wright       SETERRQ(honee->comm, PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
88*d9e69621SJames Wright               DivDiffFluxProjectionMethods[divFdiffproj_method]);
898c85b835SJames Wright       break;
908c85b835SJames Wright   }
91*d9e69621SJames Wright   *diff_flux_proj = diff_flux_proj_;
928c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
938c85b835SJames Wright };
948c85b835SJames Wright 
958c85b835SJames Wright /**
960880fbb6SJames Wright   @brief Return the objects required for the Divergence of Diffusive flux to be read by a `CeedOperator`
970880fbb6SJames Wright 
980880fbb6SJames Wright   @param[in]  diff_flux_proj Projection object
990880fbb6SJames Wright   @param[out] elem_restr     Element restriction for the divergence of diffusive flux, or `NULL`
1000880fbb6SJames Wright   @param[out] basis          Basis for the divergence of diffusive flux, or `NULL`
1010880fbb6SJames Wright   @param[out] vector         Vector for the divergence of diffusive flux, or `NULL`
1020880fbb6SJames Wright   @param[out] eval_mode      Eval mode for the divergence of diffusive flux, or `NULL`
1030880fbb6SJames Wright **/
1040880fbb6SJames Wright PetscErrorCode DivDiffFluxProjectionGetOperatorFieldData(DivDiffFluxProjectionData diff_flux_proj, CeedElemRestriction *elem_restr, CeedBasis *basis,
1050880fbb6SJames Wright                                                          CeedVector *vector, CeedEvalMode *eval_mode) {
1060880fbb6SJames Wright   Ceed ceed = CeedVectorReturnCeed(diff_flux_proj->div_diff_flux_ceed);
1070880fbb6SJames Wright 
1080880fbb6SJames Wright   PetscFunctionBeginUser;
1090880fbb6SJames Wright   if (elem_restr) PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(diff_flux_proj->elem_restr_div_diff_flux, elem_restr));
1100880fbb6SJames Wright   if (basis) PetscCallCeed(ceed, CeedBasisReferenceCopy(diff_flux_proj->basis_div_diff_flux, basis));
1110880fbb6SJames Wright   if (vector) PetscCallCeed(ceed, CeedVectorReferenceCopy(diff_flux_proj->div_diff_flux_ceed, vector));
1120880fbb6SJames Wright   if (eval_mode) *eval_mode = diff_flux_proj->eval_mode_div_diff_flux;
1130880fbb6SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1140880fbb6SJames Wright }
1150880fbb6SJames Wright 
1160880fbb6SJames Wright /**
1178c85b835SJames Wright   @brief Setup direct projection of divergence of diffusive flux
1188c85b835SJames Wright 
1190c373b74SJames Wright   @param[in]     honee          `Honee` context
1208561fee2SJames Wright   @param[in,out] diff_flux_proj Flux projection object to setup
1218c85b835SJames Wright **/
122e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionSetup_Direct(Honee honee, DivDiffFluxProjectionData diff_flux_proj) {
1230c373b74SJames Wright   Ceed                ceed       = honee->ceed;
1248c85b835SJames Wright   NodalProjectionData projection = diff_flux_proj->projection;
12536038bbcSJames Wright   MPI_Comm            comm       = PetscObjectComm((PetscObject)projection->dm);
1268c85b835SJames Wright 
1278c85b835SJames Wright   PetscFunctionBeginUser;
12836038bbcSJames Wright   {  // Create Projection RHS OperatorApplyContext
12936038bbcSJames Wright     CeedOperator op_rhs;
1308c85b835SJames Wright 
13136038bbcSJames Wright     PetscCheck(diff_flux_proj->CreateRHSOperator_Direct, comm, PETSC_ERR_ARG_WRONGSTATE,
13236038bbcSJames Wright                "Must define CreateRHSOperator_Direct to use indirect div_diff_flux projection");
133e3663b90SJames Wright     PetscCall(diff_flux_proj->CreateRHSOperator_Direct(honee, diff_flux_proj, &op_rhs));
1348c85b835SJames Wright     PetscCall(DMCreateLocalVector(projection->dm, &diff_flux_proj->DivDiffFlux_loc));
1358c85b835SJames Wright     diff_flux_proj->ceed_vec_has_array = PETSC_FALSE;
1360c373b74SJames Wright     PetscCall(OperatorApplyContextCreate(honee->dm, projection->dm, ceed, op_rhs, NULL, NULL, NULL, diff_flux_proj->DivDiffFlux_loc,
13736038bbcSJames Wright                                          &projection->l2_rhs_ctx));
13836038bbcSJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
13936038bbcSJames Wright   }
1408c85b835SJames Wright 
1418c85b835SJames Wright   {  // -- Build Mass operator
1428c85b835SJames Wright     CeedQFunction       qf_mass;
1438c85b835SJames Wright     CeedOperator        op_mass;
144590d9cddSJames Wright     CeedBasis           basis_div_diff_flux             = NULL;
145590d9cddSJames Wright     CeedElemRestriction elem_restr_div_diff_flux_volume = NULL, elem_restr_qd;
1460880fbb6SJames Wright     CeedVector          q_data;
1470880fbb6SJames Wright     CeedInt             q_data_size;
1480880fbb6SJames Wright     PetscInt            label_value  = 0;
1490880fbb6SJames Wright     DMLabel             domain_label = NULL;
1500880fbb6SJames Wright 
151590d9cddSJames Wright     PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_div_diff_flux_volume, &basis_div_diff_flux, NULL, NULL));
152e3663b90SJames 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,
153e3663b90SJames Wright                        &q_data_size));
1548c85b835SJames Wright 
15564dd23feSJames Wright     PetscCall(HoneeMassQFunctionCreate(ceed, projection->num_comp, q_data_size, &qf_mass));
1568c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
157590d9cddSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_div_diff_flux_volume, basis_div_diff_flux, CEED_VECTOR_ACTIVE));
1588c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
159590d9cddSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_div_diff_flux_volume, basis_div_diff_flux, CEED_VECTOR_ACTIVE));
1608c85b835SJames Wright 
1618c85b835SJames Wright     {  // -- Setup KSP for L^2 projection
1628c85b835SJames Wright       Mat mat_mass;
1638c85b835SJames Wright 
1648c85b835SJames Wright       PetscCall(MatCreateCeed(projection->dm, projection->dm, op_mass, NULL, &mat_mass));
1658c85b835SJames Wright 
1668c85b835SJames Wright       PetscCall(KSPCreate(comm, &projection->ksp));
1678c85b835SJames Wright       PetscCall(KSPSetOptionsPrefix(projection->ksp, "div_diff_flux_projection_"));
1688c85b835SJames Wright       {  // lumped by default
1698c85b835SJames Wright         PC pc;
1708c85b835SJames Wright         PetscCall(KSPGetPC(projection->ksp, &pc));
1718c85b835SJames Wright         PetscCall(PCSetType(pc, PCJACOBI));
1728c85b835SJames Wright         PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
1738c85b835SJames Wright         PetscCall(KSPSetType(projection->ksp, KSPPREONLY));
1748c85b835SJames Wright       }
1758c85b835SJames Wright       PetscCall(KSPSetFromOptions_WithMatCeed(projection->ksp, mat_mass));
1768c85b835SJames Wright       PetscCall(MatDestroy(&mat_mass));
1778c85b835SJames Wright     }
1780880fbb6SJames Wright 
179590d9cddSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_div_diff_flux_volume));
180590d9cddSJames Wright     PetscCallCeed(ceed, CeedBasisDestroy(&basis_div_diff_flux));
1810880fbb6SJames Wright     PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
1820880fbb6SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
1838c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
1848c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
1858c85b835SJames Wright   }
1868c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1878c85b835SJames Wright }
1888c85b835SJames Wright 
1898c85b835SJames Wright /**
1908c85b835SJames Wright   @brief Setup indirect projection of divergence of diffusive flux
1918c85b835SJames Wright 
1920c373b74SJames Wright   @param[in]     honee          `Honee` context
1938561fee2SJames Wright   @param[in,out] diff_flux_proj Flux projection object to setup
1948c85b835SJames Wright **/
195e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionSetup_Indirect(Honee honee, DivDiffFluxProjectionData diff_flux_proj) {
1960c373b74SJames Wright   Ceed                ceed       = honee->ceed;
1978c85b835SJames Wright   NodalProjectionData projection = diff_flux_proj->projection;
1988c85b835SJames Wright   CeedBasis           basis_diff_flux;
1998c85b835SJames Wright   CeedElemRestriction elem_restr_diff_flux, elem_restr_qd;
2008c85b835SJames Wright   CeedVector          q_data;
2010880fbb6SJames Wright   CeedInt             q_data_size;
20236038bbcSJames Wright   MPI_Comm            comm = PetscObjectComm((PetscObject)projection->dm);
2038c85b835SJames Wright 
2048c85b835SJames Wright   PetscFunctionBeginUser;
2050880fbb6SJames Wright   {
2060880fbb6SJames Wright     PetscInt label_value = 0, height = 0, dm_field = 0;
2070880fbb6SJames Wright     DMLabel  domain_label = NULL;
2088c85b835SJames Wright 
2098c85b835SJames Wright     PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, domain_label, label_value, height, dm_field, &elem_restr_diff_flux));
2108c85b835SJames Wright     PetscCall(CreateBasisFromPlex(ceed, projection->dm, domain_label, label_value, height, dm_field, &basis_diff_flux));
211e3663b90SJames 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,
212e3663b90SJames Wright                        &q_data_size));
2130880fbb6SJames Wright   }
2148c85b835SJames Wright 
21536038bbcSJames Wright   {
2168c85b835SJames Wright     CeedOperator op_rhs;
2178c85b835SJames Wright 
21836038bbcSJames Wright     PetscCheck(diff_flux_proj->CreateRHSOperator_Indirect, comm, PETSC_ERR_ARG_WRONGSTATE,
21936038bbcSJames Wright                "Must define CreateRHSOperator_Indirect to use indirect div_diff_flux projection");
220e3663b90SJames Wright     PetscCall(diff_flux_proj->CreateRHSOperator_Indirect(honee, diff_flux_proj, &op_rhs));
2210c373b74SJames Wright     PetscCall(OperatorApplyContextCreate(honee->dm, projection->dm, ceed, op_rhs, NULL, NULL, NULL, NULL, &projection->l2_rhs_ctx));
2228c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
2238c85b835SJames Wright   }
2248c85b835SJames Wright 
2258c85b835SJames Wright   {  // -- Build Mass operator
2268c85b835SJames Wright     CeedQFunction qf_mass;
2278c85b835SJames Wright     CeedOperator  op_mass;
2288c85b835SJames Wright 
22964dd23feSJames Wright     PetscCall(HoneeMassQFunctionCreate(ceed, projection->num_comp, q_data_size, &qf_mass));
2308c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
2318c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
2328c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
2338c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
2348c85b835SJames Wright 
2358c85b835SJames Wright     {  // -- Setup KSP for L^2 projection
2368c85b835SJames Wright       Mat mat_mass;
2378c85b835SJames Wright 
2388c85b835SJames Wright       PetscCall(MatCreateCeed(projection->dm, projection->dm, op_mass, NULL, &mat_mass));
2398c85b835SJames Wright 
2408c85b835SJames Wright       PetscCall(KSPCreate(comm, &projection->ksp));
2418c85b835SJames Wright       PetscCall(KSPSetOptionsPrefix(projection->ksp, "div_diff_flux_projection_"));
2428c85b835SJames Wright       {  // lumped by default
2438c85b835SJames Wright         PC pc;
2448c85b835SJames Wright         PetscCall(KSPGetPC(projection->ksp, &pc));
2458c85b835SJames Wright         PetscCall(PCSetType(pc, PCJACOBI));
2468c85b835SJames Wright         PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
2478c85b835SJames Wright         PetscCall(KSPSetType(projection->ksp, KSPPREONLY));
2488c85b835SJames Wright       }
2498c85b835SJames Wright       PetscCall(KSPSetFromOptions_WithMatCeed(projection->ksp, mat_mass));
2508c85b835SJames Wright       PetscCall(MatDestroy(&mat_mass));
2518c85b835SJames Wright     }
2528c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
2538c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
2548c85b835SJames Wright   }
2558c85b835SJames Wright 
2568c85b835SJames Wright   {  // Create OperatorApplyContext to calculate divergence at quadrature points
25740b78511SJames Wright     CeedQFunction       qf_calc_divergence = NULL;
2588c85b835SJames Wright     CeedOperator        op_calc_divergence;
2590880fbb6SJames Wright     CeedElemRestriction elem_restr_div_diff_flux = NULL;
2600880fbb6SJames Wright     PetscInt            dim;
2618c85b835SJames Wright 
2620880fbb6SJames Wright     PetscCall(DMGetDimension(projection->dm, &dim));
2630880fbb6SJames Wright     PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_div_diff_flux, NULL, NULL, NULL));
2648c85b835SJames Wright 
26540b78511SJames Wright     switch (dim) {
26640b78511SJames Wright       case 2:
26740b78511SJames Wright         switch (diff_flux_proj->num_diff_flux_comps) {
26840b78511SJames Wright           case 1:
26940b78511SJames Wright             PetscCallCeed(ceed,
27040b78511SJames Wright                           CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux2D_1, ComputeDivDiffusiveFlux2D_1_loc, &qf_calc_divergence));
27140b78511SJames Wright             break;
27240b78511SJames Wright         }
27340b78511SJames Wright         break;
27440b78511SJames Wright       case 3:
27540b78511SJames Wright         switch (diff_flux_proj->num_diff_flux_comps) {
27640b78511SJames Wright           case 1:
27740b78511SJames Wright             PetscCallCeed(ceed,
27840b78511SJames Wright                           CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux3D_1, ComputeDivDiffusiveFlux3D_1_loc, &qf_calc_divergence));
27940b78511SJames Wright             break;
28040b78511SJames Wright           case 4:
28140b78511SJames Wright             PetscCallCeed(ceed,
28240b78511SJames Wright                           CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux3D_4, ComputeDivDiffusiveFlux3D_4_loc, &qf_calc_divergence));
28340b78511SJames Wright             break;
28440b78511SJames Wright         }
28540b78511SJames Wright         break;
28640b78511SJames Wright     }
28740b78511SJames Wright     PetscCheck(qf_calc_divergence, comm, PETSC_ERR_SUP,
28840b78511SJames Wright                "QFunction for calculating divergence of diffusive flux does not exist for"
28940b78511SJames Wright                " %" PetscInt_FMT " dimensional grid and %" PetscInt_FMT
29040b78511SJames Wright                " number of components.\nA new qfunction can be easily added; see source code for pattern.",
29140b78511SJames Wright                dim, diff_flux_proj->num_diff_flux_comps);
2928c85b835SJames Wright 
2938c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_calc_divergence, "Grad F_diff", projection->num_comp * dim, CEED_EVAL_GRAD));
2948c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_calc_divergence, "qdata", q_data_size, CEED_EVAL_NONE));
29540b78511SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_calc_divergence, "Div F_diff", diff_flux_proj->num_diff_flux_comps, CEED_EVAL_NONE));
2968c85b835SJames Wright 
2978c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_calc_divergence, NULL, NULL, &op_calc_divergence));
2988c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "Grad F_diff", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
2998c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
3008c85b835SJames Wright     PetscCallCeed(
3018c85b835SJames Wright         ceed, CeedOperatorSetField(op_calc_divergence, "Div F_diff", elem_restr_div_diff_flux, CEED_BASIS_NONE, diff_flux_proj->div_diff_flux_ceed));
3028c85b835SJames Wright 
30336038bbcSJames Wright     PetscCall(
30436038bbcSJames Wright         OperatorApplyContextCreate(projection->dm, NULL, ceed, op_calc_divergence, NULL, NULL, NULL, NULL, &diff_flux_proj->calc_div_diff_flux));
3058c85b835SJames Wright 
3060880fbb6SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_div_diff_flux));
3078c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_calc_divergence));
3088c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_calc_divergence));
3098c85b835SJames Wright   }
3108c85b835SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
3118c85b835SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
3128c85b835SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
3138c85b835SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux));
3148c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3158c85b835SJames Wright }
3168c85b835SJames Wright 
3178c85b835SJames Wright /**
3188c85b835SJames Wright   @brief Setup projection of divergence of diffusive flux
3198c85b835SJames Wright 
3200c373b74SJames Wright   @param[in]     honee          `Honee` context
3218561fee2SJames Wright   @param[in,out] diff_flux_proj Flux projection object to setup
3228c85b835SJames Wright **/
323e3663b90SJames Wright PetscErrorCode DivDiffFluxProjectionSetup(Honee honee, DivDiffFluxProjectionData diff_flux_proj) {
3248c85b835SJames Wright   PetscFunctionBeginUser;
3250c373b74SJames Wright   switch (honee->app_ctx->divFdiffproj_method) {
3268c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_DIRECT:
327e3663b90SJames Wright       PetscCall(DivDiffFluxProjectionSetup_Direct(honee, diff_flux_proj));
3288c85b835SJames Wright       break;
3298c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_INDIRECT:
330e3663b90SJames Wright       PetscCall(DivDiffFluxProjectionSetup_Indirect(honee, diff_flux_proj));
3318c85b835SJames Wright       break;
3328c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_NONE:
3330c373b74SJames Wright       SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
3340c373b74SJames Wright               DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]);
3358c85b835SJames Wright       break;
3368c85b835SJames Wright   }
3378c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3388c85b835SJames Wright }
3398c85b835SJames Wright 
3408c85b835SJames Wright /**
3418c85b835SJames Wright   @brief Project the divergence of diffusive flux
3428c85b835SJames Wright 
3438c85b835SJames Wright   This implicitly sets the `CeedVector` input (`div_diff_flux_ceed`) to the divergence of diffusive flux.
3448c85b835SJames Wright 
3458c85b835SJames Wright   @param[in]  diff_flux_proj `NodalProjectionData` for the projection
3468c85b835SJames Wright   @param[in]  Q_loc          Localized solution vector
3478c85b835SJames Wright **/
34836038bbcSJames Wright PetscErrorCode DivDiffFluxProjectionApply(DivDiffFluxProjectionData diff_flux_proj, Vec Q_loc) {
3498c85b835SJames Wright   NodalProjectionData projection = diff_flux_proj->projection;
3508c85b835SJames Wright 
3518c85b835SJames Wright   PetscFunctionBeginUser;
352ea615d4cSJames Wright   PetscCall(PetscLogEventBegin(HONEE_DivDiffFluxProjection, Q_loc, 0, 0, 0));
3538c85b835SJames Wright   switch (diff_flux_proj->method) {
3548c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_DIRECT: {
355c8cb0c9bSJames Wright       Vec DivDiffFlux, RHS;
3568c85b835SJames Wright 
3578c85b835SJames Wright       PetscCall(DMGetGlobalVector(projection->dm, &DivDiffFlux));
358c8cb0c9bSJames Wright       PetscCall(DMGetGlobalVector(projection->dm, &RHS));
3598c85b835SJames Wright       if (diff_flux_proj->ceed_vec_has_array) {
3608c85b835SJames Wright         PetscCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->DivDiffFlux_loc));
3618c85b835SJames Wright         diff_flux_proj->ceed_vec_has_array = PETSC_FALSE;
3628c85b835SJames Wright       }
363c8cb0c9bSJames Wright       PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, RHS, projection->l2_rhs_ctx));
3648c85b835SJames Wright       PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_rhs_view"));
3658c85b835SJames Wright 
366c8cb0c9bSJames Wright       {
367c8cb0c9bSJames Wright         // Run PCApply manually if using ksp_type preonly -pc_type jacobi
368c8cb0c9bSJames Wright         // This is to avoid an AllReduce call in KSPSolve_Preonly, which causes significant slowdowns for lumped mass matrix solves.
369c8cb0c9bSJames Wright         // See https://gitlab.com/petsc/petsc/-/merge_requests/8048 for more details and a possible fix
370c8cb0c9bSJames Wright         PC        pc;
371c8cb0c9bSJames Wright         PetscBool ispreonly, isjacobi;
372c8cb0c9bSJames Wright         PetscCall(KSPGetPC(projection->ksp, &pc));
373c8cb0c9bSJames Wright         PetscCall(PetscObjectTypeCompare((PetscObject)projection->ksp, KSPPREONLY, &ispreonly));
374c8cb0c9bSJames Wright         PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &isjacobi));
375c8cb0c9bSJames Wright         if (ispreonly && isjacobi) PetscCall(PCApply(pc, RHS, DivDiffFlux));
376c8cb0c9bSJames Wright         else PetscCall(KSPSolve(projection->ksp, RHS, DivDiffFlux));
377c8cb0c9bSJames Wright       }
3788c85b835SJames Wright       PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_view"));
3798c85b835SJames Wright 
3808c85b835SJames Wright       PetscCall(DMGlobalToLocal(projection->dm, DivDiffFlux, INSERT_VALUES, diff_flux_proj->DivDiffFlux_loc));
3818c85b835SJames Wright       PetscCall(VecReadPetscToCeed(diff_flux_proj->DivDiffFlux_loc, &diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->div_diff_flux_ceed));
3828c85b835SJames Wright       diff_flux_proj->ceed_vec_has_array = PETSC_TRUE;
3838c85b835SJames Wright 
384c8cb0c9bSJames Wright       PetscCall(DMRestoreGlobalVector(projection->dm, &RHS));
3858c85b835SJames Wright       PetscCall(DMRestoreGlobalVector(projection->dm, &DivDiffFlux));
3868c85b835SJames Wright       break;
3878c85b835SJames Wright     }
3888c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_INDIRECT: {
389c8cb0c9bSJames Wright       Vec DiffFlux, RHS;
3908c85b835SJames Wright 
3918c85b835SJames Wright       PetscCall(DMGetGlobalVector(projection->dm, &DiffFlux));
392c8cb0c9bSJames Wright       PetscCall(DMGetGlobalVector(projection->dm, &RHS));
393c8cb0c9bSJames Wright       PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, RHS, projection->l2_rhs_ctx));
3948c85b835SJames Wright       PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_rhs_view"));
3958c85b835SJames Wright 
396c8cb0c9bSJames Wright       {
397c8cb0c9bSJames Wright         // Run PCApply manually if using -ksp_type preonly -pc_type jacobi
398c8cb0c9bSJames Wright         // This is to avoid an AllReduce call in KSPSolve_Preonly, which causes significant slowdowns for lumped mass matrix solves.
399c8cb0c9bSJames Wright         // See https://gitlab.com/petsc/petsc/-/merge_requests/8048 for more details and a possible fix
400c8cb0c9bSJames Wright         PC        pc;
401c8cb0c9bSJames Wright         PetscBool ispreonly, isjacobi;
402c8cb0c9bSJames Wright         PetscCall(KSPGetPC(projection->ksp, &pc));
403c8cb0c9bSJames Wright         PetscCall(PetscObjectTypeCompare((PetscObject)projection->ksp, KSPPREONLY, &ispreonly));
404c8cb0c9bSJames Wright         PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &isjacobi));
405c8cb0c9bSJames Wright         if (ispreonly && isjacobi) PetscCall(PCApply(pc, RHS, DiffFlux));
406c8cb0c9bSJames Wright         else PetscCall(KSPSolve(projection->ksp, RHS, DiffFlux));
407c8cb0c9bSJames Wright       }
4088c85b835SJames Wright       PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_view"));
4098c85b835SJames Wright 
4108c85b835SJames Wright       PetscCall(ApplyCeedOperatorGlobalToLocal(DiffFlux, NULL, diff_flux_proj->calc_div_diff_flux));
411c8cb0c9bSJames Wright       PetscCall(DMRestoreGlobalVector(projection->dm, &RHS));
4128c85b835SJames Wright       PetscCall(DMRestoreGlobalVector(projection->dm, &DiffFlux));
4138c85b835SJames Wright     } break;
4148c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_NONE:
4158c85b835SJames Wright       SETERRQ(PetscObjectComm((PetscObject)projection->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
4168c85b835SJames Wright               DivDiffFluxProjectionMethods[diff_flux_proj->method]);
4178c85b835SJames Wright       break;
4188c85b835SJames Wright   }
419ea615d4cSJames Wright   PetscCall(PetscLogEventEnd(HONEE_DivDiffFluxProjection, Q_loc, 0, 0, 0));
4208c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4218c85b835SJames Wright }
4228c85b835SJames Wright 
4238c85b835SJames Wright /**
4248c85b835SJames Wright   @brief Destroy `DivDiffFluxProjectionData` object
4258c85b835SJames Wright 
4268c85b835SJames Wright   @param[in,out] diff_flux_proj Object to destroy
4278c85b835SJames Wright **/
4288c85b835SJames Wright PetscErrorCode DivDiffFluxProjectionDataDestroy(DivDiffFluxProjectionData diff_flux_proj) {
4298c85b835SJames Wright   PetscFunctionBeginUser;
4308c85b835SJames Wright   if (diff_flux_proj == NULL) PetscFunctionReturn(PETSC_SUCCESS);
4310880fbb6SJames Wright   Ceed ceed = CeedVectorReturnCeed(diff_flux_proj->div_diff_flux_ceed);
4320880fbb6SJames Wright 
4334ea616f4SJames Wright   PetscCall(NodalProjectionDataDestroy(&diff_flux_proj->projection));
4348c85b835SJames Wright   PetscCall(OperatorApplyContextDestroy(diff_flux_proj->calc_div_diff_flux));
4358c85b835SJames Wright   if (diff_flux_proj->ceed_vec_has_array) {
4368c85b835SJames Wright     PetscCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->DivDiffFlux_loc));
4378c85b835SJames Wright     diff_flux_proj->ceed_vec_has_array = PETSC_FALSE;
4388c85b835SJames Wright   }
4390880fbb6SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&diff_flux_proj->div_diff_flux_ceed));
4400880fbb6SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&diff_flux_proj->elem_restr_div_diff_flux));
4410880fbb6SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&diff_flux_proj->basis_div_diff_flux));
4428c85b835SJames Wright   PetscCall(VecDestroy(&diff_flux_proj->DivDiffFlux_loc));
4438c85b835SJames Wright   PetscCall(PetscFree(diff_flux_proj));
4448c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4458c85b835SJames Wright }
446