xref: /honee/src/diff_flux_projection.c (revision 64dd23fed3bf7ffd926f9c8ee5c8cd1891ca06e3)
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 
128c85b835SJames Wright /**
130c373b74SJames Wright   @brief Create `DivDiffFluxProjectionData` for solution DM in `honee`
148c85b835SJames Wright 
150c373b74SJames Wright   @param[in]  honee               `Honee` context
1636038bbcSJames Wright   @param[in]  num_diff_flux_comps Number of components that makes up the diffusive flux (e.g. 1 for scalar advection-diffusion)
1736038bbcSJames Wright   @param[out] pdiff_flux_proj     The `DivDiffFluxProjectionData` object created, or `NULL` if not created
188c85b835SJames Wright **/
190c373b74SJames Wright PetscErrorCode DivDiffFluxProjectionCreate(Honee honee, PetscInt num_diff_flux_comps, DivDiffFluxProjectionData *pdiff_flux_proj) {
200c373b74SJames Wright   PetscInt                  label_value = 0, height = 0, dm_field = 0, dim, degree = honee->app_ctx->degree, q_extra = honee->app_ctx->q_extra;
218c85b835SJames Wright   DMLabel                   domain_label = NULL;
228c85b835SJames Wright   DivDiffFluxProjectionData diff_flux_proj;
238c85b835SJames Wright   NodalProjectionData       projection;
248c85b835SJames Wright 
258c85b835SJames Wright   PetscFunctionBeginUser;
260c373b74SJames Wright   if (honee->app_ctx->divFdiffproj_method == DIV_DIFF_FLUX_PROJ_NONE) {
2736038bbcSJames Wright     *pdiff_flux_proj = NULL;
2836038bbcSJames Wright     PetscFunctionReturn(PETSC_SUCCESS);
2936038bbcSJames Wright   }
3036038bbcSJames Wright   PetscCall(PetscNew(pdiff_flux_proj));
3136038bbcSJames Wright   diff_flux_proj = *pdiff_flux_proj;
320c373b74SJames Wright   PetscCall(PetscNew(&honee->diff_flux_proj->projection));
330c373b74SJames Wright   projection                          = honee->diff_flux_proj->projection;
340c373b74SJames Wright   diff_flux_proj->method              = honee->app_ctx->divFdiffproj_method;
3536038bbcSJames Wright   diff_flux_proj->num_diff_flux_comps = num_diff_flux_comps;
368c85b835SJames Wright 
370c373b74SJames Wright   PetscCall(DMClone(honee->dm, &projection->dm));
388c85b835SJames Wright   PetscCall(DMSetMatrixPreallocateSkip(projection->dm, PETSC_TRUE));
398c85b835SJames Wright   PetscCall(DMGetDimension(projection->dm, &dim));
408c85b835SJames Wright   switch (diff_flux_proj->method) {
418c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_DIRECT: {
428c85b835SJames Wright       projection->num_comp = diff_flux_proj->num_diff_flux_comps;
438c85b835SJames Wright       PetscCall(PetscObjectSetName((PetscObject)projection->dm, "DivDiffFluxProj"));
4436038bbcSJames Wright       PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &projection->num_comp, projection->dm));
458c85b835SJames Wright 
460c373b74SJames Wright       PetscCall(DMPlexCeedElemRestrictionCreate(honee->ceed, projection->dm, domain_label, label_value, height, dm_field,
4736038bbcSJames Wright                                                 &diff_flux_proj->elem_restr_div_diff_flux));
480c373b74SJames Wright       PetscCallCeed(honee->ceed,
490c373b74SJames Wright                     CeedElemRestrictionCreateVector(diff_flux_proj->elem_restr_div_diff_flux, &diff_flux_proj->div_diff_flux_ceed, NULL));
500c373b74SJames Wright       PetscCall(CreateBasisFromPlex(honee->ceed, projection->dm, domain_label, label_value, height, dm_field, &diff_flux_proj->basis_div_diff_flux));
5136038bbcSJames Wright       diff_flux_proj->eval_mode_div_diff_flux = CEED_EVAL_INTERP;
521af555e8SJames Wright 
531af555e8SJames Wright       {  // Create face labels on projection->dm for boundary integrals
541af555e8SJames Wright         DMLabel  face_sets_label;
551af555e8SJames Wright         PetscInt num_face_set_values, *face_set_values;
561af555e8SJames Wright 
570c373b74SJames Wright         PetscCall(DMGetLabel(honee->dm, "Face Sets", &face_sets_label));
580c373b74SJames Wright         PetscCall(DMLabelCreateGlobalValueArray(honee->dm, face_sets_label, &num_face_set_values, &face_set_values));
591af555e8SJames Wright         for (PetscInt f = 0; f < num_face_set_values; f++) {
601af555e8SJames Wright           DMLabel face_orientation_label;
611af555e8SJames Wright           char   *face_orientation_label_name;
621af555e8SJames Wright 
630c373b74SJames Wright           PetscCall(DMPlexCreateFaceLabel(honee->dm, face_set_values[f], &face_orientation_label_name));
640c373b74SJames Wright           PetscCall(DMGetLabel(honee->dm, face_orientation_label_name, &face_orientation_label));
651af555e8SJames Wright           PetscCall(DMAddLabel(projection->dm, face_orientation_label));
661af555e8SJames Wright           PetscCall(PetscFree(face_orientation_label_name));
671af555e8SJames Wright         }
6836038bbcSJames Wright         PetscCall(PetscFree(face_set_values));
691af555e8SJames Wright       }
708c85b835SJames Wright     } break;
718c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_INDIRECT: {
728c85b835SJames Wright       projection->num_comp = diff_flux_proj->num_diff_flux_comps * dim;
738c85b835SJames Wright       PetscCall(PetscObjectSetName((PetscObject)projection->dm, "DiffFluxProj"));
7436038bbcSJames Wright       PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &projection->num_comp, projection->dm));
758c85b835SJames Wright 
760c373b74SJames Wright       PetscCall(DMPlexCeedElemRestrictionQDataCreate(honee->ceed, projection->dm, domain_label, label_value, height,
7736038bbcSJames Wright                                                      diff_flux_proj->num_diff_flux_comps, &diff_flux_proj->elem_restr_div_diff_flux));
780c373b74SJames Wright       PetscCallCeed(honee->ceed,
790c373b74SJames Wright                     CeedElemRestrictionCreateVector(diff_flux_proj->elem_restr_div_diff_flux, &diff_flux_proj->div_diff_flux_ceed, NULL));
8036038bbcSJames Wright       diff_flux_proj->basis_div_diff_flux     = CEED_BASIS_NONE;
8136038bbcSJames Wright       diff_flux_proj->eval_mode_div_diff_flux = CEED_EVAL_NONE;
828c85b835SJames Wright     } break;
838c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_NONE:
840c373b74SJames Wright       SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
850c373b74SJames Wright               DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]);
868c85b835SJames Wright       break;
878c85b835SJames Wright   }
888c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
898c85b835SJames Wright };
908c85b835SJames Wright 
918c85b835SJames Wright /**
920880fbb6SJames Wright   @brief Return the objects required for the Divergence of Diffusive flux to be read by a `CeedOperator`
930880fbb6SJames Wright 
940880fbb6SJames Wright   @param[in]  diff_flux_proj Projection object
950880fbb6SJames Wright   @param[out] elem_restr     Element restriction for the divergence of diffusive flux, or `NULL`
960880fbb6SJames Wright   @param[out] basis          Basis for the divergence of diffusive flux, or `NULL`
970880fbb6SJames Wright   @param[out] vector         Vector for the divergence of diffusive flux, or `NULL`
980880fbb6SJames Wright   @param[out] eval_mode      Eval mode for the divergence of diffusive flux, or `NULL`
990880fbb6SJames Wright **/
1000880fbb6SJames Wright PetscErrorCode DivDiffFluxProjectionGetOperatorFieldData(DivDiffFluxProjectionData diff_flux_proj, CeedElemRestriction *elem_restr, CeedBasis *basis,
1010880fbb6SJames Wright                                                          CeedVector *vector, CeedEvalMode *eval_mode) {
1020880fbb6SJames Wright   Ceed ceed = CeedVectorReturnCeed(diff_flux_proj->div_diff_flux_ceed);
1030880fbb6SJames Wright 
1040880fbb6SJames Wright   PetscFunctionBeginUser;
1050880fbb6SJames Wright   if (elem_restr) PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(diff_flux_proj->elem_restr_div_diff_flux, elem_restr));
1060880fbb6SJames Wright   if (basis) PetscCallCeed(ceed, CeedBasisReferenceCopy(diff_flux_proj->basis_div_diff_flux, basis));
1070880fbb6SJames Wright   if (vector) PetscCallCeed(ceed, CeedVectorReferenceCopy(diff_flux_proj->div_diff_flux_ceed, vector));
1080880fbb6SJames Wright   if (eval_mode) *eval_mode = diff_flux_proj->eval_mode_div_diff_flux;
1090880fbb6SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1100880fbb6SJames Wright }
1110880fbb6SJames Wright 
1120880fbb6SJames Wright /**
1138c85b835SJames Wright   @brief Setup direct projection of divergence of diffusive flux
1148c85b835SJames Wright 
1150c373b74SJames Wright   @param[in]     honee          `Honee` context
1168561fee2SJames Wright   @param[in,out] diff_flux_proj Flux projection object to setup
1178c85b835SJames Wright **/
118e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionSetup_Direct(Honee honee, DivDiffFluxProjectionData diff_flux_proj) {
1190c373b74SJames Wright   Ceed                ceed       = honee->ceed;
1208c85b835SJames Wright   NodalProjectionData projection = diff_flux_proj->projection;
12136038bbcSJames Wright   MPI_Comm            comm       = PetscObjectComm((PetscObject)projection->dm);
1228c85b835SJames Wright 
1238c85b835SJames Wright   PetscFunctionBeginUser;
12436038bbcSJames Wright   {  // Create Projection RHS OperatorApplyContext
12536038bbcSJames Wright     CeedOperator op_rhs;
1268c85b835SJames Wright 
12736038bbcSJames Wright     PetscCheck(diff_flux_proj->CreateRHSOperator_Direct, comm, PETSC_ERR_ARG_WRONGSTATE,
12836038bbcSJames Wright                "Must define CreateRHSOperator_Direct to use indirect div_diff_flux projection");
129e3663b90SJames Wright     PetscCall(diff_flux_proj->CreateRHSOperator_Direct(honee, diff_flux_proj, &op_rhs));
1308c85b835SJames Wright     PetscCall(DMCreateLocalVector(projection->dm, &diff_flux_proj->DivDiffFlux_loc));
1318c85b835SJames Wright     diff_flux_proj->ceed_vec_has_array = PETSC_FALSE;
1320c373b74SJames Wright     PetscCall(OperatorApplyContextCreate(honee->dm, projection->dm, ceed, op_rhs, NULL, NULL, NULL, diff_flux_proj->DivDiffFlux_loc,
13336038bbcSJames Wright                                          &projection->l2_rhs_ctx));
13436038bbcSJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
13536038bbcSJames Wright   }
1368c85b835SJames Wright 
1378c85b835SJames Wright   {  // -- Build Mass operator
1388c85b835SJames Wright     CeedQFunction       qf_mass;
1398c85b835SJames Wright     CeedOperator        op_mass;
140590d9cddSJames Wright     CeedBasis           basis_div_diff_flux             = NULL;
141590d9cddSJames Wright     CeedElemRestriction elem_restr_div_diff_flux_volume = NULL, elem_restr_qd;
1420880fbb6SJames Wright     CeedVector          q_data;
1430880fbb6SJames Wright     CeedInt             q_data_size;
1440880fbb6SJames Wright     PetscInt            label_value  = 0;
1450880fbb6SJames Wright     DMLabel             domain_label = NULL;
1460880fbb6SJames Wright 
147590d9cddSJames Wright     PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_div_diff_flux_volume, &basis_div_diff_flux, NULL, NULL));
148e3663b90SJames 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,
149e3663b90SJames Wright                        &q_data_size));
1508c85b835SJames Wright 
151*64dd23feSJames Wright     PetscCall(HoneeMassQFunctionCreate(ceed, projection->num_comp, q_data_size, &qf_mass));
1528c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
153590d9cddSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_div_diff_flux_volume, basis_div_diff_flux, CEED_VECTOR_ACTIVE));
1548c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
155590d9cddSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_div_diff_flux_volume, basis_div_diff_flux, CEED_VECTOR_ACTIVE));
1568c85b835SJames Wright 
1578c85b835SJames Wright     {  // -- Setup KSP for L^2 projection
1588c85b835SJames Wright       Mat mat_mass;
1598c85b835SJames Wright 
1608c85b835SJames Wright       PetscCall(MatCreateCeed(projection->dm, projection->dm, op_mass, NULL, &mat_mass));
1618c85b835SJames Wright 
1628c85b835SJames Wright       PetscCall(KSPCreate(comm, &projection->ksp));
1638c85b835SJames Wright       PetscCall(KSPSetOptionsPrefix(projection->ksp, "div_diff_flux_projection_"));
1648c85b835SJames Wright       {  // lumped by default
1658c85b835SJames Wright         PC pc;
1668c85b835SJames Wright         PetscCall(KSPGetPC(projection->ksp, &pc));
1678c85b835SJames Wright         PetscCall(PCSetType(pc, PCJACOBI));
1688c85b835SJames Wright         PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
1698c85b835SJames Wright         PetscCall(KSPSetType(projection->ksp, KSPPREONLY));
1708c85b835SJames Wright       }
1718c85b835SJames Wright       PetscCall(KSPSetFromOptions_WithMatCeed(projection->ksp, mat_mass));
1728c85b835SJames Wright       PetscCall(MatDestroy(&mat_mass));
1738c85b835SJames Wright     }
1740880fbb6SJames Wright 
175590d9cddSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_div_diff_flux_volume));
176590d9cddSJames Wright     PetscCallCeed(ceed, CeedBasisDestroy(&basis_div_diff_flux));
1770880fbb6SJames Wright     PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
1780880fbb6SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
1798c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
1808c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
1818c85b835SJames Wright   }
1828c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1838c85b835SJames Wright }
1848c85b835SJames Wright 
1858c85b835SJames Wright /**
1868c85b835SJames Wright   @brief Setup indirect projection of divergence of diffusive flux
1878c85b835SJames Wright 
1880c373b74SJames Wright   @param[in]     honee          `Honee` context
1898561fee2SJames Wright   @param[in,out] diff_flux_proj Flux projection object to setup
1908c85b835SJames Wright **/
191e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionSetup_Indirect(Honee honee, DivDiffFluxProjectionData diff_flux_proj) {
1920c373b74SJames Wright   Ceed                ceed       = honee->ceed;
1938c85b835SJames Wright   NodalProjectionData projection = diff_flux_proj->projection;
1948c85b835SJames Wright   CeedBasis           basis_diff_flux;
1958c85b835SJames Wright   CeedElemRestriction elem_restr_diff_flux, elem_restr_qd;
1968c85b835SJames Wright   CeedVector          q_data;
1970880fbb6SJames Wright   CeedInt             q_data_size;
19836038bbcSJames Wright   MPI_Comm            comm = PetscObjectComm((PetscObject)projection->dm);
1998c85b835SJames Wright 
2008c85b835SJames Wright   PetscFunctionBeginUser;
2010880fbb6SJames Wright   {
2020880fbb6SJames Wright     PetscInt label_value = 0, height = 0, dm_field = 0;
2030880fbb6SJames Wright     DMLabel  domain_label = NULL;
2048c85b835SJames Wright 
2058c85b835SJames Wright     PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, domain_label, label_value, height, dm_field, &elem_restr_diff_flux));
2068c85b835SJames Wright     PetscCall(CreateBasisFromPlex(ceed, projection->dm, domain_label, label_value, height, dm_field, &basis_diff_flux));
207e3663b90SJames 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,
208e3663b90SJames Wright                        &q_data_size));
2090880fbb6SJames Wright   }
2108c85b835SJames Wright 
21136038bbcSJames Wright   {
2128c85b835SJames Wright     CeedOperator op_rhs;
2138c85b835SJames Wright 
21436038bbcSJames Wright     PetscCheck(diff_flux_proj->CreateRHSOperator_Indirect, comm, PETSC_ERR_ARG_WRONGSTATE,
21536038bbcSJames Wright                "Must define CreateRHSOperator_Indirect to use indirect div_diff_flux projection");
216e3663b90SJames Wright     PetscCall(diff_flux_proj->CreateRHSOperator_Indirect(honee, diff_flux_proj, &op_rhs));
2170c373b74SJames Wright     PetscCall(OperatorApplyContextCreate(honee->dm, projection->dm, ceed, op_rhs, NULL, NULL, NULL, NULL, &projection->l2_rhs_ctx));
2188c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
2198c85b835SJames Wright   }
2208c85b835SJames Wright 
2218c85b835SJames Wright   {  // -- Build Mass operator
2228c85b835SJames Wright     CeedQFunction qf_mass;
2238c85b835SJames Wright     CeedOperator  op_mass;
2248c85b835SJames Wright 
225*64dd23feSJames Wright     PetscCall(HoneeMassQFunctionCreate(ceed, projection->num_comp, q_data_size, &qf_mass));
2268c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
2278c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
2288c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
2298c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
2308c85b835SJames Wright 
2318c85b835SJames Wright     {  // -- Setup KSP for L^2 projection
2328c85b835SJames Wright       Mat mat_mass;
2338c85b835SJames Wright 
2348c85b835SJames Wright       PetscCall(MatCreateCeed(projection->dm, projection->dm, op_mass, NULL, &mat_mass));
2358c85b835SJames Wright 
2368c85b835SJames Wright       PetscCall(KSPCreate(comm, &projection->ksp));
2378c85b835SJames Wright       PetscCall(KSPSetOptionsPrefix(projection->ksp, "div_diff_flux_projection_"));
2388c85b835SJames Wright       {  // lumped by default
2398c85b835SJames Wright         PC pc;
2408c85b835SJames Wright         PetscCall(KSPGetPC(projection->ksp, &pc));
2418c85b835SJames Wright         PetscCall(PCSetType(pc, PCJACOBI));
2428c85b835SJames Wright         PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
2438c85b835SJames Wright         PetscCall(KSPSetType(projection->ksp, KSPPREONLY));
2448c85b835SJames Wright       }
2458c85b835SJames Wright       PetscCall(KSPSetFromOptions_WithMatCeed(projection->ksp, mat_mass));
2468c85b835SJames Wright       PetscCall(MatDestroy(&mat_mass));
2478c85b835SJames Wright     }
2488c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
2498c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
2508c85b835SJames Wright   }
2518c85b835SJames Wright 
2528c85b835SJames Wright   {  // Create OperatorApplyContext to calculate divergence at quadrature points
25340b78511SJames Wright     CeedQFunction       qf_calc_divergence = NULL;
2548c85b835SJames Wright     CeedOperator        op_calc_divergence;
2550880fbb6SJames Wright     CeedElemRestriction elem_restr_div_diff_flux = NULL;
2560880fbb6SJames Wright     PetscInt            dim;
2578c85b835SJames Wright 
2580880fbb6SJames Wright     PetscCall(DMGetDimension(projection->dm, &dim));
2590880fbb6SJames Wright     PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_div_diff_flux, NULL, NULL, NULL));
2608c85b835SJames Wright 
26140b78511SJames Wright     switch (dim) {
26240b78511SJames Wright       case 2:
26340b78511SJames Wright         switch (diff_flux_proj->num_diff_flux_comps) {
26440b78511SJames Wright           case 1:
26540b78511SJames Wright             PetscCallCeed(ceed,
26640b78511SJames Wright                           CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux2D_1, ComputeDivDiffusiveFlux2D_1_loc, &qf_calc_divergence));
26740b78511SJames Wright             break;
26840b78511SJames Wright         }
26940b78511SJames Wright         break;
27040b78511SJames Wright       case 3:
27140b78511SJames Wright         switch (diff_flux_proj->num_diff_flux_comps) {
27240b78511SJames Wright           case 1:
27340b78511SJames Wright             PetscCallCeed(ceed,
27440b78511SJames Wright                           CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux3D_1, ComputeDivDiffusiveFlux3D_1_loc, &qf_calc_divergence));
27540b78511SJames Wright             break;
27640b78511SJames Wright           case 4:
27740b78511SJames Wright             PetscCallCeed(ceed,
27840b78511SJames Wright                           CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux3D_4, ComputeDivDiffusiveFlux3D_4_loc, &qf_calc_divergence));
27940b78511SJames Wright             break;
28040b78511SJames Wright         }
28140b78511SJames Wright         break;
28240b78511SJames Wright     }
28340b78511SJames Wright     PetscCheck(qf_calc_divergence, comm, PETSC_ERR_SUP,
28440b78511SJames Wright                "QFunction for calculating divergence of diffusive flux does not exist for"
28540b78511SJames Wright                " %" PetscInt_FMT " dimensional grid and %" PetscInt_FMT
28640b78511SJames Wright                " number of components.\nA new qfunction can be easily added; see source code for pattern.",
28740b78511SJames Wright                dim, diff_flux_proj->num_diff_flux_comps);
2888c85b835SJames Wright 
2898c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_calc_divergence, "Grad F_diff", projection->num_comp * dim, CEED_EVAL_GRAD));
2908c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_calc_divergence, "qdata", q_data_size, CEED_EVAL_NONE));
29140b78511SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_calc_divergence, "Div F_diff", diff_flux_proj->num_diff_flux_comps, CEED_EVAL_NONE));
2928c85b835SJames Wright 
2938c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_calc_divergence, NULL, NULL, &op_calc_divergence));
2948c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "Grad F_diff", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
2958c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
2968c85b835SJames Wright     PetscCallCeed(
2978c85b835SJames Wright         ceed, CeedOperatorSetField(op_calc_divergence, "Div F_diff", elem_restr_div_diff_flux, CEED_BASIS_NONE, diff_flux_proj->div_diff_flux_ceed));
2988c85b835SJames Wright 
29936038bbcSJames Wright     PetscCall(
30036038bbcSJames Wright         OperatorApplyContextCreate(projection->dm, NULL, ceed, op_calc_divergence, NULL, NULL, NULL, NULL, &diff_flux_proj->calc_div_diff_flux));
3018c85b835SJames Wright 
3020880fbb6SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_div_diff_flux));
3038c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_calc_divergence));
3048c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_calc_divergence));
3058c85b835SJames Wright   }
3068c85b835SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
3078c85b835SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
3088c85b835SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
3098c85b835SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux));
3108c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3118c85b835SJames Wright }
3128c85b835SJames Wright 
3138c85b835SJames Wright /**
3148c85b835SJames Wright   @brief Setup projection of divergence of diffusive flux
3158c85b835SJames Wright 
3160c373b74SJames Wright   @param[in]     honee          `Honee` context
3178561fee2SJames Wright   @param[in,out] diff_flux_proj Flux projection object to setup
3188c85b835SJames Wright **/
319e3663b90SJames Wright PetscErrorCode DivDiffFluxProjectionSetup(Honee honee, DivDiffFluxProjectionData diff_flux_proj) {
3208c85b835SJames Wright   PetscFunctionBeginUser;
3210c373b74SJames Wright   switch (honee->app_ctx->divFdiffproj_method) {
3228c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_DIRECT:
323e3663b90SJames Wright       PetscCall(DivDiffFluxProjectionSetup_Direct(honee, diff_flux_proj));
3248c85b835SJames Wright       break;
3258c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_INDIRECT:
326e3663b90SJames Wright       PetscCall(DivDiffFluxProjectionSetup_Indirect(honee, diff_flux_proj));
3278c85b835SJames Wright       break;
3288c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_NONE:
3290c373b74SJames Wright       SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
3300c373b74SJames Wright               DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]);
3318c85b835SJames Wright       break;
3328c85b835SJames Wright   }
3338c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3348c85b835SJames Wright }
3358c85b835SJames Wright 
3368c85b835SJames Wright /**
3378c85b835SJames Wright   @brief Project the divergence of diffusive flux
3388c85b835SJames Wright 
3398c85b835SJames Wright   This implicitly sets the `CeedVector` input (`div_diff_flux_ceed`) to the divergence of diffusive flux.
3408c85b835SJames Wright 
3418c85b835SJames Wright   @param[in]  diff_flux_proj `NodalProjectionData` for the projection
3428c85b835SJames Wright   @param[in]  Q_loc          Localized solution vector
3438c85b835SJames Wright **/
34436038bbcSJames Wright PetscErrorCode DivDiffFluxProjectionApply(DivDiffFluxProjectionData diff_flux_proj, Vec Q_loc) {
3458c85b835SJames Wright   NodalProjectionData projection = diff_flux_proj->projection;
3468c85b835SJames Wright 
3478c85b835SJames Wright   PetscFunctionBeginUser;
3488c85b835SJames Wright   PetscCall(PetscLogEventBegin(FLUIDS_DivDiffFluxProjection, Q_loc, 0, 0, 0));
3498c85b835SJames Wright   switch (diff_flux_proj->method) {
3508c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_DIRECT: {
3518c85b835SJames Wright       Vec DivDiffFlux;
3528c85b835SJames Wright 
3538c85b835SJames Wright       PetscCall(DMGetGlobalVector(projection->dm, &DivDiffFlux));
3548c85b835SJames Wright       if (diff_flux_proj->ceed_vec_has_array) {
3558c85b835SJames Wright         PetscCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->DivDiffFlux_loc));
3568c85b835SJames Wright         diff_flux_proj->ceed_vec_has_array = PETSC_FALSE;
3578c85b835SJames Wright       }
3588c85b835SJames Wright       PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, DivDiffFlux, projection->l2_rhs_ctx));
3598c85b835SJames Wright       PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_rhs_view"));
3608c85b835SJames Wright 
3618c85b835SJames Wright       PetscCall(KSPSolve(projection->ksp, DivDiffFlux, DivDiffFlux));
3628c85b835SJames Wright       PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_view"));
3638c85b835SJames Wright 
3648c85b835SJames Wright       PetscCall(DMGlobalToLocal(projection->dm, DivDiffFlux, INSERT_VALUES, diff_flux_proj->DivDiffFlux_loc));
3658c85b835SJames Wright       PetscCall(VecReadPetscToCeed(diff_flux_proj->DivDiffFlux_loc, &diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->div_diff_flux_ceed));
3668c85b835SJames Wright       diff_flux_proj->ceed_vec_has_array = PETSC_TRUE;
3678c85b835SJames Wright 
3688c85b835SJames Wright       PetscCall(DMRestoreGlobalVector(projection->dm, &DivDiffFlux));
3698c85b835SJames Wright       break;
3708c85b835SJames Wright     }
3718c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_INDIRECT: {
3728c85b835SJames Wright       Vec DiffFlux;
3738c85b835SJames Wright 
3748c85b835SJames Wright       PetscCall(DMGetGlobalVector(projection->dm, &DiffFlux));
3758c85b835SJames Wright       PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, DiffFlux, projection->l2_rhs_ctx));
3768c85b835SJames Wright       PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_rhs_view"));
3778c85b835SJames Wright 
3788c85b835SJames Wright       PetscCall(KSPSolve(projection->ksp, DiffFlux, DiffFlux));
3798c85b835SJames Wright       PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_view"));
3808c85b835SJames Wright 
3818c85b835SJames Wright       PetscCall(ApplyCeedOperatorGlobalToLocal(DiffFlux, NULL, diff_flux_proj->calc_div_diff_flux));
3828c85b835SJames Wright       PetscCall(DMRestoreGlobalVector(projection->dm, &DiffFlux));
3838c85b835SJames Wright     } break;
3848c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_NONE:
3858c85b835SJames Wright       SETERRQ(PetscObjectComm((PetscObject)projection->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
3868c85b835SJames Wright               DivDiffFluxProjectionMethods[diff_flux_proj->method]);
3878c85b835SJames Wright       break;
3888c85b835SJames Wright   }
3898c85b835SJames Wright   PetscCall(PetscLogEventEnd(FLUIDS_DivDiffFluxProjection, Q_loc, 0, 0, 0));
3908c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3918c85b835SJames Wright }
3928c85b835SJames Wright 
3938c85b835SJames Wright /**
3948c85b835SJames Wright   @brief Destroy `DivDiffFluxProjectionData` object
3958c85b835SJames Wright 
3968c85b835SJames Wright   @param[in,out] diff_flux_proj Object to destroy
3978c85b835SJames Wright **/
3988c85b835SJames Wright PetscErrorCode DivDiffFluxProjectionDataDestroy(DivDiffFluxProjectionData diff_flux_proj) {
3998c85b835SJames Wright   PetscFunctionBeginUser;
4008c85b835SJames Wright   if (diff_flux_proj == NULL) PetscFunctionReturn(PETSC_SUCCESS);
4010880fbb6SJames Wright   Ceed ceed = CeedVectorReturnCeed(diff_flux_proj->div_diff_flux_ceed);
4020880fbb6SJames Wright 
4038c85b835SJames Wright   PetscCall(NodalProjectionDataDestroy(diff_flux_proj->projection));
4048c85b835SJames Wright   PetscCall(OperatorApplyContextDestroy(diff_flux_proj->calc_div_diff_flux));
4058c85b835SJames Wright   if (diff_flux_proj->ceed_vec_has_array) {
4068c85b835SJames Wright     PetscCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->DivDiffFlux_loc));
4078c85b835SJames Wright     diff_flux_proj->ceed_vec_has_array = PETSC_FALSE;
4088c85b835SJames Wright   }
4090880fbb6SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&diff_flux_proj->div_diff_flux_ceed));
4100880fbb6SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&diff_flux_proj->elem_restr_div_diff_flux));
4110880fbb6SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&diff_flux_proj->basis_div_diff_flux));
4128c85b835SJames Wright   PetscCall(VecDestroy(&diff_flux_proj->DivDiffFlux_loc));
4138c85b835SJames Wright   PetscCall(PetscFree(diff_flux_proj));
4148c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4158c85b835SJames Wright }
416