xref: /honee/src/diff_flux_projection.c (revision 1af555e8a723c759a2e546a290a4d15c89159397)
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 /**
138c85b835SJames Wright   @brief Initialize projection of divergence of diffusive flux
148c85b835SJames Wright 
158c85b835SJames Wright   Creates underlying `DM` for the projection operation and creates the restriction and basis to use with the CeedOperator
168c85b835SJames Wright 
178c85b835SJames Wright   @param[in]  user                     `User` context
188c85b835SJames Wright   @param[out] elem_restr_div_diff_flux `CeedElemRestriction` of the divergence of diffusive flux vector
198c85b835SJames Wright   @param[out] basis_div_diff_flux      `CeedBasis` of the divergence of diffusive flux vector
208c85b835SJames Wright   @param[out] eval_mode_diff_flux      `CeedEvalMode` for the divergence of the diffusive flux
218c85b835SJames Wright **/
228c85b835SJames Wright PetscErrorCode DiffFluxProjectionInitialize(User user, CeedElemRestriction *elem_restr_div_diff_flux, CeedBasis *basis_div_diff_flux,
238c85b835SJames Wright                                             CeedEvalMode *eval_mode_diff_flux) {
248c85b835SJames Wright   PetscSection              section;
258c85b835SJames Wright   PetscInt                  label_value = 0, height = 0, dm_field = 0, dim;
268c85b835SJames Wright   DMLabel                   domain_label = NULL;
278c85b835SJames Wright   DivDiffFluxProjectionData diff_flux_proj;
288c85b835SJames Wright   NodalProjectionData       projection;
298c85b835SJames Wright 
308c85b835SJames Wright   PetscFunctionBeginUser;
318c85b835SJames Wright   PetscCall(PetscNew(&user->diff_flux_proj));
328c85b835SJames Wright   diff_flux_proj = user->diff_flux_proj;
338c85b835SJames Wright   PetscCall(PetscNew(&user->diff_flux_proj->projection));
348c85b835SJames Wright   projection                          = user->diff_flux_proj->projection;
358c85b835SJames Wright   diff_flux_proj->method              = user->app_ctx->divFdiffproj_method;
368c85b835SJames Wright   diff_flux_proj->num_diff_flux_comps = 4;
378c85b835SJames Wright 
388c85b835SJames Wright   PetscCall(DMClone(user->dm, &projection->dm));
398c85b835SJames Wright   PetscCall(DMSetMatrixPreallocateSkip(projection->dm, PETSC_TRUE));
408c85b835SJames Wright   PetscCall(DMGetDimension(projection->dm, &dim));
418c85b835SJames Wright   switch (diff_flux_proj->method) {
428c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_DIRECT: {
438c85b835SJames Wright       projection->num_comp = diff_flux_proj->num_diff_flux_comps;
448c85b835SJames Wright       PetscCall(PetscObjectSetName((PetscObject)projection->dm, "DivDiffFluxProj"));
458c85b835SJames Wright       PetscCall(
468c85b835SJames Wright           DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, user->app_ctx->degree, 1, user->app_ctx->q_extra, 1, &projection->num_comp, projection->dm));
478c85b835SJames Wright 
488c85b835SJames Wright       PetscCall(DMGetLocalSection(projection->dm, &section));
498c85b835SJames Wright       PetscCall(PetscSectionSetFieldName(section, 0, ""));
508c85b835SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 0, "DivDiffusiveFlux_MomentumX"));
518c85b835SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 1, "DivDiffusiveFlux_MomentumY"));
528c85b835SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 2, "DivDiffusiveFlux_MomentumZ"));
538c85b835SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 3, "DivDiffusiveFlux_Energy"));
548c85b835SJames Wright 
558c85b835SJames Wright       PetscCall(DMPlexCeedElemRestrictionCreate(user->ceed, projection->dm, domain_label, label_value, height, dm_field, elem_restr_div_diff_flux));
568c85b835SJames Wright       PetscCallCeed(user->ceed, CeedElemRestrictionCreateVector(*elem_restr_div_diff_flux, &diff_flux_proj->div_diff_flux_ceed, NULL));
578c85b835SJames Wright       PetscCall(CreateBasisFromPlex(user->ceed, projection->dm, domain_label, label_value, height, dm_field, basis_div_diff_flux));
588c85b835SJames Wright       *eval_mode_diff_flux = CEED_EVAL_INTERP;
59*1af555e8SJames Wright 
60*1af555e8SJames Wright       {  // Create face labels on projection->dm for boundary integrals
61*1af555e8SJames Wright         DMLabel  face_sets_label;
62*1af555e8SJames Wright         PetscInt num_face_set_values, *face_set_values;
63*1af555e8SJames Wright 
64*1af555e8SJames Wright         PetscCall(DMGetLabel(user->dm, "Face Sets", &face_sets_label));
65*1af555e8SJames Wright         PetscCall(DMLabelCreateGlobalValueArray(user->dm, face_sets_label, &num_face_set_values, &face_set_values));
66*1af555e8SJames Wright         for (PetscInt f = 0; f < num_face_set_values; f++) {
67*1af555e8SJames Wright           DMLabel face_orientation_label;
68*1af555e8SJames Wright           char   *face_orientation_label_name;
69*1af555e8SJames Wright 
70*1af555e8SJames Wright           PetscCall(DMPlexCreateFaceLabel(user->dm, face_set_values[f], &face_orientation_label_name));
71*1af555e8SJames Wright           PetscCall(DMGetLabel(user->dm, face_orientation_label_name, &face_orientation_label));
72*1af555e8SJames Wright           PetscCall(DMAddLabel(projection->dm, face_orientation_label));
73*1af555e8SJames Wright           PetscCall(PetscFree(face_orientation_label_name));
74*1af555e8SJames Wright         }
75*1af555e8SJames Wright       }
768c85b835SJames Wright     } break;
778c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_INDIRECT: {
788c85b835SJames Wright       projection->num_comp = diff_flux_proj->num_diff_flux_comps * dim;
798c85b835SJames Wright       PetscCall(PetscObjectSetName((PetscObject)projection->dm, "DiffFluxProj"));
808c85b835SJames Wright       PetscCall(
818c85b835SJames Wright           DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, user->app_ctx->degree, 1, user->app_ctx->q_extra, 1, &projection->num_comp, projection->dm));
828c85b835SJames Wright 
838c85b835SJames Wright       PetscCall(DMGetLocalSection(projection->dm, &section));
848c85b835SJames Wright       PetscCall(PetscSectionSetFieldName(section, 0, ""));
858c85b835SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 0, "DiffusiveFlux_MomentumXX"));
868c85b835SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 1, "DiffusiveFlux_MomentumXY"));
878c85b835SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 2, "DiffusiveFlux_MomentumXZ"));
888c85b835SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 3, "DiffusiveFlux_MomentumYX"));
898c85b835SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 4, "DiffusiveFlux_MomentumYY"));
908c85b835SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 5, "DiffusiveFlux_MomentumYZ"));
918c85b835SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 6, "DiffusiveFlux_MomentumZX"));
928c85b835SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 7, "DiffusiveFlux_MomentumZY"));
938c85b835SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 8, "DiffusiveFlux_MomentumZZ"));
948c85b835SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 9, "DiffusiveFlux_EnergyX"));
958c85b835SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 10, "DiffusiveFlux_EnergyY"));
968c85b835SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 11, "DiffusiveFlux_EnergyZ"));
978c85b835SJames Wright 
988c85b835SJames Wright       PetscCall(DMPlexCeedElemRestrictionQDataCreate(user->ceed, projection->dm, domain_label, label_value, height,
998c85b835SJames Wright                                                      diff_flux_proj->num_diff_flux_comps, elem_restr_div_diff_flux));
1008c85b835SJames Wright       PetscCallCeed(user->ceed, CeedElemRestrictionCreateVector(*elem_restr_div_diff_flux, &diff_flux_proj->div_diff_flux_ceed, NULL));
1018c85b835SJames Wright       *basis_div_diff_flux = CEED_BASIS_NONE;
1028c85b835SJames Wright       *eval_mode_diff_flux = CEED_EVAL_NONE;
1038c85b835SJames Wright     } break;
1048c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_NONE:
1058c85b835SJames Wright       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
1068c85b835SJames Wright               DivDiffFluxProjectionMethods[user->app_ctx->divFdiffproj_method]);
1078c85b835SJames Wright       break;
1088c85b835SJames Wright   }
1098c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1108c85b835SJames Wright };
1118c85b835SJames Wright 
1128c85b835SJames Wright /**
1138c85b835SJames Wright   @brief Setup direct projection of divergence of diffusive flux
1148c85b835SJames Wright 
1158c85b835SJames Wright   @param[in] ceed      `Ceed` context
1168c85b835SJames Wright   @param[in] user      `User` context
1178c85b835SJames Wright   @param[in] ceed_data `CeedData` context
1188c85b835SJames Wright   @param[in] problem   `ProblemData` context
1198c85b835SJames Wright **/
1208c85b835SJames Wright static PetscErrorCode DivDiffFluxProjectionSetup_Direct(Ceed ceed, User user, CeedData ceed_data, ProblemData problem) {
1218c85b835SJames Wright   DivDiffFluxProjectionData diff_flux_proj = user->diff_flux_proj;
1228c85b835SJames Wright   NodalProjectionData       projection     = diff_flux_proj->projection;
1238c85b835SJames Wright   CeedOperator              op_rhs;
1248c85b835SJames Wright   CeedBasis                 basis_diff_flux;
1258c85b835SJames Wright   CeedElemRestriction       elem_restr_diff_flux_volume, elem_restr_qd;
1268c85b835SJames Wright   CeedVector                q_data;
1278c85b835SJames Wright   CeedInt                   num_comp_q, q_data_size;
1288c85b835SJames Wright   PetscInt                  dim, label_value = 0;
1298c85b835SJames Wright   DMLabel                   domain_label = NULL;
1308c85b835SJames Wright 
1318c85b835SJames Wright   PetscFunctionBeginUser;
1328c85b835SJames Wright   // -- Get Pre-requisite things
1338c85b835SJames Wright   PetscCall(DMGetDimension(projection->dm, &dim));
1348c85b835SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q));
1358c85b835SJames Wright 
1368c85b835SJames Wright   {  // Get elem_restr_diff_flux and basis_diff_flux
1378c85b835SJames Wright     CeedOperator     *sub_ops;
1388c85b835SJames Wright     CeedOperatorField op_field;
1398c85b835SJames Wright     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
1408c85b835SJames Wright 
1418c85b835SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_ifunction, &sub_ops));
1428c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "div F_diff", &op_field));
1438c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_diff_flux_volume));
1448c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetBasis(op_field, &basis_diff_flux));
1458c85b835SJames Wright   }
1468c85b835SJames Wright   PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, &elem_restr_qd,
1478c85b835SJames Wright                      &q_data, &q_data_size));
1488c85b835SJames Wright 
1498c85b835SJames Wright   PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_rhs));
1508c85b835SJames Wright   {  // Add the volume integral CeedOperator
1518c85b835SJames Wright     CeedQFunction qf_rhs_volume;
1528c85b835SJames Wright     CeedOperator  op_rhs_volume;
1538c85b835SJames Wright 
1548c85b835SJames Wright     switch (user->phys->state_var) {
1558c85b835SJames Wright       case STATEVAR_PRIMITIVE:
1568c85b835SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_Prim, DivDiffusiveFluxVolumeRHS_Prim_loc, &qf_rhs_volume));
1578c85b835SJames Wright         break;
1588c85b835SJames Wright       case STATEVAR_CONSERVATIVE:
1598c85b835SJames Wright         PetscCallCeed(ceed,
1608c85b835SJames Wright                       CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_Conserv, DivDiffusiveFluxVolumeRHS_Conserv_loc, &qf_rhs_volume));
1618c85b835SJames Wright         break;
1628c85b835SJames Wright       case STATEVAR_ENTROPY:
1638c85b835SJames Wright         PetscCallCeed(ceed,
1648c85b835SJames Wright                       CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_Entropy, DivDiffusiveFluxVolumeRHS_Entropy_loc, &qf_rhs_volume));
1658c85b835SJames Wright         break;
1668c85b835SJames Wright     }
1678c85b835SJames Wright 
1688c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_volume, problem->apply_vol_ifunction.qfctx));
1698c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "q", num_comp_q, CEED_EVAL_INTERP));
1708c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
1718c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "qdata", q_data_size, CEED_EVAL_NONE));
1728c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_volume, "diffusive flux RHS", projection->num_comp * dim, CEED_EVAL_GRAD));
1738c85b835SJames Wright 
1748c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_volume, NULL, NULL, &op_rhs_volume));
1758c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
1768c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
1778c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
1788c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "diffusive flux RHS", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE));
1798c85b835SJames Wright 
1808c85b835SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_rhs, op_rhs_volume));
1818c85b835SJames Wright 
1828c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_volume));
1838c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_volume));
1848c85b835SJames Wright   }
1858c85b835SJames Wright 
1868c85b835SJames Wright   {  // Add the boundary integral CeedOperator
1878c85b835SJames Wright     CeedQFunction qf_rhs_boundary;
1888c85b835SJames Wright     DMLabel       face_sets_label;
1898c85b835SJames Wright     PetscInt      num_face_set_values, *face_set_values;
1908c85b835SJames Wright     CeedInt       q_data_size;
1918c85b835SJames Wright 
1928c85b835SJames Wright     // -- Build RHS operator
1938c85b835SJames Wright     switch (user->phys->state_var) {
1948c85b835SJames Wright       case STATEVAR_PRIMITIVE:
1958c85b835SJames Wright         PetscCallCeed(ceed,
1968c85b835SJames Wright                       CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_Prim, DivDiffusiveFluxBoundaryRHS_Prim_loc, &qf_rhs_boundary));
1978c85b835SJames Wright         break;
1988c85b835SJames Wright       case STATEVAR_CONSERVATIVE:
1998c85b835SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_Conserv, DivDiffusiveFluxBoundaryRHS_Conserv_loc,
2008c85b835SJames Wright                                                         &qf_rhs_boundary));
2018c85b835SJames Wright         break;
2028c85b835SJames Wright       case STATEVAR_ENTROPY:
2038c85b835SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_Entropy, DivDiffusiveFluxBoundaryRHS_Entropy_loc,
2048c85b835SJames Wright                                                         &qf_rhs_boundary));
2058c85b835SJames Wright         break;
2068c85b835SJames Wright     }
2078c85b835SJames Wright 
2088c85b835SJames Wright     PetscCall(QDataBoundaryGradientGetNumComponents(user->dm, &q_data_size));
2098c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_boundary, problem->apply_vol_ifunction.qfctx));
2108c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "q", num_comp_q, CEED_EVAL_INTERP));
2118c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
2128c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "qdata", q_data_size, CEED_EVAL_NONE));
2138c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_boundary, "diffusive flux RHS", projection->num_comp, CEED_EVAL_INTERP));
2148c85b835SJames Wright 
215*1af555e8SJames Wright     PetscCall(DMGetLabel(projection->dm, "Face Sets", &face_sets_label));
216*1af555e8SJames Wright     PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_sets_label, &num_face_set_values, &face_set_values));
2178c85b835SJames Wright     for (PetscInt f = 0; f < num_face_set_values; f++) {
2188c85b835SJames Wright       DMLabel  face_orientation_label;
2198c85b835SJames Wright       PetscInt num_orientations_values, *orientation_values;
2208c85b835SJames Wright 
2218c85b835SJames Wright       {
2228c85b835SJames Wright         char *face_orientation_label_name;
2238c85b835SJames Wright 
224*1af555e8SJames Wright         PetscCall(DMPlexCreateFaceLabel(projection->dm, face_set_values[f], &face_orientation_label_name));
225*1af555e8SJames Wright         PetscCall(DMGetLabel(projection->dm, face_orientation_label_name, &face_orientation_label));
2268c85b835SJames Wright         PetscCall(PetscFree(face_orientation_label_name));
2278c85b835SJames Wright       }
228*1af555e8SJames Wright       PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_orientation_label, &num_orientations_values, &orientation_values));
2298c85b835SJames Wright       for (PetscInt o = 0; o < num_orientations_values; o++) {
2308c85b835SJames Wright         CeedOperator        op_rhs_boundary;
2318c85b835SJames Wright         CeedBasis           basis_q, basis_diff_flux_boundary;
2328c85b835SJames Wright         CeedElemRestriction elem_restr_qdata, elem_restr_q, elem_restr_diff_flux_boundary;
2338c85b835SJames Wright         CeedVector          q_data;
2348c85b835SJames Wright         CeedInt             q_data_size;
2358c85b835SJames Wright         PetscInt            orientation = orientation_values[o], dm_field_q = 0, height_cell = 0, height_face = 1;
2368c85b835SJames Wright 
2378c85b835SJames Wright         PetscCall(DMPlexCeedElemRestrictionCreate(ceed, user->dm, face_orientation_label, orientation, height_cell, dm_field_q, &elem_restr_q));
2388c85b835SJames Wright         PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0,
2398c85b835SJames Wright                                                   &elem_restr_diff_flux_boundary));
2408c85b835SJames Wright         PetscCall(QDataBoundaryGradientGet(ceed, user->dm, face_orientation_label, orientation, ceed_data->x_coord, &elem_restr_qdata, &q_data,
2418c85b835SJames Wright                                            &q_data_size));
2428c85b835SJames Wright         PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, user->dm, face_orientation_label, orientation, orientation, dm_field_q, &basis_q));
2438c85b835SJames Wright         PetscCall(CreateBasisFromPlex(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, &basis_diff_flux_boundary));
2448c85b835SJames Wright 
2458c85b835SJames Wright         PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_boundary, NULL, NULL, &op_rhs_boundary));
2468c85b835SJames Wright         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
2478c85b835SJames Wright         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
2488c85b835SJames Wright         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "qdata", elem_restr_qdata, CEED_BASIS_NONE, q_data));
2498c85b835SJames Wright         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "diffusive flux RHS", elem_restr_diff_flux_boundary, basis_diff_flux_boundary,
2508c85b835SJames Wright                                                  CEED_VECTOR_ACTIVE));
2518c85b835SJames Wright 
2528c85b835SJames Wright         PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_rhs, op_rhs_boundary));
2538c85b835SJames Wright 
2548c85b835SJames Wright         PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_boundary));
2558c85b835SJames Wright         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qdata));
2568c85b835SJames Wright         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
2578c85b835SJames Wright         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_boundary));
2588c85b835SJames Wright         PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
2598c85b835SJames Wright         PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux_boundary));
2608c85b835SJames Wright         PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
2618c85b835SJames Wright       }
2628c85b835SJames Wright       PetscCall(PetscFree(orientation_values));
2638c85b835SJames Wright     }
2648c85b835SJames Wright     PetscCall(PetscFree(face_set_values));
2658c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_boundary));
2668c85b835SJames Wright   }
2678c85b835SJames Wright 
2688c85b835SJames Wright   PetscCall(DMCreateLocalVector(projection->dm, &diff_flux_proj->DivDiffFlux_loc));
2698c85b835SJames Wright   diff_flux_proj->ceed_vec_has_array = PETSC_FALSE;
2708c85b835SJames Wright   PetscCall(
2718c85b835SJames Wright       OperatorApplyContextCreate(user->dm, projection->dm, ceed, op_rhs, NULL, NULL, NULL, diff_flux_proj->DivDiffFlux_loc, &projection->l2_rhs_ctx));
2728c85b835SJames Wright 
2738c85b835SJames Wright   {  // -- Build Mass operator
2748c85b835SJames Wright     CeedQFunction qf_mass;
2758c85b835SJames Wright     CeedOperator  op_mass;
2768c85b835SJames Wright 
2778c85b835SJames Wright     PetscCall(CreateMassQFunction(ceed, projection->num_comp, q_data_size, &qf_mass));
2788c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
2798c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE));
2808c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
2818c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE));
2828c85b835SJames Wright 
2838c85b835SJames Wright     {  // -- Setup KSP for L^2 projection
2848c85b835SJames Wright       Mat      mat_mass;
2858c85b835SJames Wright       MPI_Comm comm = PetscObjectComm((PetscObject)projection->dm);
2868c85b835SJames Wright 
2878c85b835SJames Wright       PetscCall(MatCreateCeed(projection->dm, projection->dm, op_mass, NULL, &mat_mass));
2888c85b835SJames Wright 
2898c85b835SJames Wright       PetscCall(KSPCreate(comm, &projection->ksp));
2908c85b835SJames Wright       PetscCall(KSPSetOptionsPrefix(projection->ksp, "div_diff_flux_projection_"));
2918c85b835SJames Wright       {  // lumped by default
2928c85b835SJames Wright         PC pc;
2938c85b835SJames Wright         PetscCall(KSPGetPC(projection->ksp, &pc));
2948c85b835SJames Wright         PetscCall(PCSetType(pc, PCJACOBI));
2958c85b835SJames Wright         PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
2968c85b835SJames Wright         PetscCall(KSPSetType(projection->ksp, KSPPREONLY));
2978c85b835SJames Wright       }
2988c85b835SJames Wright       PetscCall(KSPSetFromOptions_WithMatCeed(projection->ksp, mat_mass));
2998c85b835SJames Wright       PetscCall(MatDestroy(&mat_mass));
3008c85b835SJames Wright     }
3018c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
3028c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
3038c85b835SJames Wright   }
3048c85b835SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
3058c85b835SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
3068c85b835SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
3078c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3088c85b835SJames Wright }
3098c85b835SJames Wright 
3108c85b835SJames Wright /**
3118c85b835SJames Wright   @brief Setup indirect projection of divergence of diffusive flux
3128c85b835SJames Wright 
3138c85b835SJames Wright   @param[in]     ceed      `Ceed` context
3148c85b835SJames Wright   @param[in,out] user      `User` context
3158c85b835SJames Wright   @param[in]     ceed_data `CeedData` context
3168c85b835SJames Wright   @param[in]     problem   `ProblemData` context
3178c85b835SJames Wright **/
3188c85b835SJames Wright static PetscErrorCode DivDiffFluxProjectionSetup_Indirect(Ceed ceed, User user, CeedData ceed_data, ProblemData problem) {
3198c85b835SJames Wright   DivDiffFluxProjectionData diff_flux_proj = user->diff_flux_proj;
3208c85b835SJames Wright   NodalProjectionData       projection     = diff_flux_proj->projection;
3218c85b835SJames Wright   CeedBasis                 basis_diff_flux;
3228c85b835SJames Wright   CeedElemRestriction       elem_restr_diff_flux, elem_restr_qd;
3238c85b835SJames Wright   CeedVector                q_data;
3248c85b835SJames Wright   CeedInt                   num_comp_q, q_data_size;
3258c85b835SJames Wright   PetscInt                  dim;
3268c85b835SJames Wright   PetscInt                  label_value = 0, height = 0, dm_field = 0;
3278c85b835SJames Wright   DMLabel                   domain_label = NULL;
3288c85b835SJames Wright 
3298c85b835SJames Wright   PetscFunctionBeginUser;
3308c85b835SJames Wright   PetscCall(DMGetDimension(projection->dm, &dim));
3318c85b835SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q));
3328c85b835SJames Wright 
3338c85b835SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, domain_label, label_value, height, dm_field, &elem_restr_diff_flux));
3348c85b835SJames Wright   PetscCall(CreateBasisFromPlex(ceed, projection->dm, domain_label, label_value, height, dm_field, &basis_diff_flux));
3358c85b835SJames Wright   PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, &elem_restr_qd,
3368c85b835SJames Wright                      &q_data, &q_data_size));
3378c85b835SJames Wright 
3388c85b835SJames Wright   {  // Create RHS CeedOperator for L^2 projection
3398c85b835SJames Wright     CeedQFunction qf_rhs;
3408c85b835SJames Wright     CeedOperator  op_rhs;
3418c85b835SJames Wright 
3428c85b835SJames Wright     switch (user->phys->state_var) {
3438c85b835SJames Wright       case STATEVAR_PRIMITIVE:
3448c85b835SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_Prim, DiffusiveFluxRHS_Prim_loc, &qf_rhs));
3458c85b835SJames Wright         break;
3468c85b835SJames Wright       case STATEVAR_CONSERVATIVE:
3478c85b835SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_Conserv, DiffusiveFluxRHS_Conserv_loc, &qf_rhs));
3488c85b835SJames Wright         break;
3498c85b835SJames Wright       case STATEVAR_ENTROPY:
3508c85b835SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_Entropy, DiffusiveFluxRHS_Entropy_loc, &qf_rhs));
3518c85b835SJames Wright         break;
3528c85b835SJames Wright     }
3538c85b835SJames Wright 
3548c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs, problem->apply_vol_ifunction.qfctx));
3558c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "q", num_comp_q, CEED_EVAL_INTERP));
3568c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
3578c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "qdata", q_data_size, CEED_EVAL_NONE));
3588c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs, "F_diff RHS", projection->num_comp, CEED_EVAL_INTERP));
3598c85b835SJames Wright 
3608c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, &op_rhs));
3618c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
3628c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
3638c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
3648c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs, "F_diff RHS", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
3658c85b835SJames Wright 
3668c85b835SJames Wright     PetscCall(OperatorApplyContextCreate(user->dm, projection->dm, ceed, op_rhs, NULL, NULL, NULL, NULL, &projection->l2_rhs_ctx));
3678c85b835SJames Wright 
3688c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs));
3698c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
3708c85b835SJames Wright   }
3718c85b835SJames Wright 
3728c85b835SJames Wright   {  // -- Build Mass operator
3738c85b835SJames Wright     CeedQFunction qf_mass;
3748c85b835SJames Wright     CeedOperator  op_mass;
3758c85b835SJames Wright 
3768c85b835SJames Wright     PetscCall(CreateMassQFunction(ceed, projection->num_comp, q_data_size, &qf_mass));
3778c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
3788c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
3798c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
3808c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
3818c85b835SJames Wright 
3828c85b835SJames Wright     {  // -- Setup KSP for L^2 projection
3838c85b835SJames Wright       Mat      mat_mass;
3848c85b835SJames Wright       MPI_Comm comm = PetscObjectComm((PetscObject)projection->dm);
3858c85b835SJames Wright 
3868c85b835SJames Wright       PetscCall(MatCreateCeed(projection->dm, projection->dm, op_mass, NULL, &mat_mass));
3878c85b835SJames Wright 
3888c85b835SJames Wright       PetscCall(KSPCreate(comm, &projection->ksp));
3898c85b835SJames Wright       PetscCall(KSPSetOptionsPrefix(projection->ksp, "div_diff_flux_projection_"));
3908c85b835SJames Wright       {  // lumped by default
3918c85b835SJames Wright         PC pc;
3928c85b835SJames Wright         PetscCall(KSPGetPC(projection->ksp, &pc));
3938c85b835SJames Wright         PetscCall(PCSetType(pc, PCJACOBI));
3948c85b835SJames Wright         PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
3958c85b835SJames Wright         PetscCall(KSPSetType(projection->ksp, KSPPREONLY));
3968c85b835SJames Wright       }
3978c85b835SJames Wright       PetscCall(KSPSetFromOptions_WithMatCeed(projection->ksp, mat_mass));
3988c85b835SJames Wright       PetscCall(MatDestroy(&mat_mass));
3998c85b835SJames Wright     }
4008c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
4018c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
4028c85b835SJames Wright   }
4038c85b835SJames Wright 
4048c85b835SJames Wright   {  // Create OperatorApplyContext to calculate divergence at quadrature points
4058c85b835SJames Wright     CeedQFunction       qf_calc_divergence;
4068c85b835SJames Wright     CeedOperator        op_calc_divergence;
4078c85b835SJames Wright     CeedElemRestriction elem_restr_div_diff_flux;
4088c85b835SJames Wright 
4098c85b835SJames Wright     {  // Get elem_restr_div_diff_flux
4108c85b835SJames Wright       CeedOperator     *sub_ops;
4118c85b835SJames Wright       CeedOperatorField op_field;
4128c85b835SJames Wright       PetscInt          sub_op_index = 0;  // will be 0 for the volume op
4138c85b835SJames Wright 
4148c85b835SJames Wright       PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_ifunction, &sub_ops));
4158c85b835SJames Wright       PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "div F_diff", &op_field));
4168c85b835SJames Wright       PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_div_diff_flux));
4178c85b835SJames Wright     }
4188c85b835SJames Wright 
4198c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux3D_4, ComputeDivDiffusiveFlux3D_4_loc, &qf_calc_divergence));
4208c85b835SJames Wright 
4218c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_calc_divergence, "Grad F_diff", projection->num_comp * dim, CEED_EVAL_GRAD));
4228c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_calc_divergence, "qdata", q_data_size, CEED_EVAL_NONE));
4238c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_calc_divergence, "Div F_diff", 4, CEED_EVAL_NONE));
4248c85b835SJames Wright 
4258c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_calc_divergence, NULL, NULL, &op_calc_divergence));
4268c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "Grad F_diff", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
4278c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
4288c85b835SJames Wright     PetscCallCeed(
4298c85b835SJames Wright         ceed, CeedOperatorSetField(op_calc_divergence, "Div F_diff", elem_restr_div_diff_flux, CEED_BASIS_NONE, diff_flux_proj->div_diff_flux_ceed));
4308c85b835SJames Wright 
4318c85b835SJames Wright     PetscCall(OperatorApplyContextCreate(projection->dm, NULL, ceed, op_calc_divergence, NULL, NULL, NULL, NULL,
4328c85b835SJames Wright                                          &user->diff_flux_proj->calc_div_diff_flux));
4338c85b835SJames Wright 
4348c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_calc_divergence));
4358c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_calc_divergence));
4368c85b835SJames Wright   }
4378c85b835SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
4388c85b835SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
4398c85b835SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
4408c85b835SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux));
4418c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4428c85b835SJames Wright }
4438c85b835SJames Wright 
4448c85b835SJames Wright /**
4458c85b835SJames Wright   @brief Setup projection of divergence of diffusive flux
4468c85b835SJames Wright 
4478c85b835SJames Wright   @param[in]     ceed      `Ceed` context
4488c85b835SJames Wright   @param[in,out] user      `User` context
4498c85b835SJames Wright   @param[in]     ceed_data `CeedData` context
4508c85b835SJames Wright   @param[in]     problem   `ProblemData` context
4518c85b835SJames Wright **/
4528c85b835SJames Wright PetscErrorCode DivDiffFluxProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem) {
4538c85b835SJames Wright   PetscFunctionBeginUser;
4548c85b835SJames Wright   switch (user->app_ctx->divFdiffproj_method) {
4558c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_DIRECT:
4568c85b835SJames Wright       PetscCall(DivDiffFluxProjectionSetup_Direct(ceed, user, ceed_data, problem));
4578c85b835SJames Wright       break;
4588c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_INDIRECT:
4598c85b835SJames Wright       PetscCall(DivDiffFluxProjectionSetup_Indirect(ceed, user, ceed_data, problem));
4608c85b835SJames Wright       break;
4618c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_NONE:
4628c85b835SJames Wright       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
4638c85b835SJames Wright               DivDiffFluxProjectionMethods[user->app_ctx->divFdiffproj_method]);
4648c85b835SJames Wright       break;
4658c85b835SJames Wright   }
4668c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4678c85b835SJames Wright }
4688c85b835SJames Wright 
4698c85b835SJames Wright /**
4708c85b835SJames Wright   @brief Project the divergence of diffusive flux
4718c85b835SJames Wright 
4728c85b835SJames Wright   This implicitly sets the `CeedVector` input (`div_diff_flux_ceed`) to the divergence of diffusive flux.
4738c85b835SJames Wright 
4748c85b835SJames Wright   @param[in]  diff_flux_proj `NodalProjectionData` for the projection
4758c85b835SJames Wright   @param[in]  Q_loc          Localized solution vector
4768c85b835SJames Wright **/
4778c85b835SJames Wright PetscErrorCode DiffFluxProjectionApply(DivDiffFluxProjectionData diff_flux_proj, Vec Q_loc) {
4788c85b835SJames Wright   NodalProjectionData projection = diff_flux_proj->projection;
4798c85b835SJames Wright 
4808c85b835SJames Wright   PetscFunctionBeginUser;
4818c85b835SJames Wright   PetscCall(PetscLogEventBegin(FLUIDS_DivDiffFluxProjection, Q_loc, 0, 0, 0));
4828c85b835SJames Wright   switch (diff_flux_proj->method) {
4838c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_DIRECT: {
4848c85b835SJames Wright       Vec DivDiffFlux;
4858c85b835SJames Wright 
4868c85b835SJames Wright       PetscCall(DMGetGlobalVector(projection->dm, &DivDiffFlux));
4878c85b835SJames Wright       if (diff_flux_proj->ceed_vec_has_array) {
4888c85b835SJames Wright         PetscCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->DivDiffFlux_loc));
4898c85b835SJames Wright         diff_flux_proj->ceed_vec_has_array = PETSC_FALSE;
4908c85b835SJames Wright       }
4918c85b835SJames Wright       PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, DivDiffFlux, projection->l2_rhs_ctx));
4928c85b835SJames Wright       PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_rhs_view"));
4938c85b835SJames Wright 
4948c85b835SJames Wright       PetscCall(KSPSolve(projection->ksp, DivDiffFlux, DivDiffFlux));
4958c85b835SJames Wright       PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_view"));
4968c85b835SJames Wright 
4978c85b835SJames Wright       PetscCall(DMGlobalToLocal(projection->dm, DivDiffFlux, INSERT_VALUES, diff_flux_proj->DivDiffFlux_loc));
4988c85b835SJames Wright       PetscCall(VecReadPetscToCeed(diff_flux_proj->DivDiffFlux_loc, &diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->div_diff_flux_ceed));
4998c85b835SJames Wright       diff_flux_proj->ceed_vec_has_array = PETSC_TRUE;
5008c85b835SJames Wright 
5018c85b835SJames Wright       PetscCall(DMRestoreGlobalVector(projection->dm, &DivDiffFlux));
5028c85b835SJames Wright       break;
5038c85b835SJames Wright     }
5048c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_INDIRECT: {
5058c85b835SJames Wright       Vec DiffFlux;
5068c85b835SJames Wright 
5078c85b835SJames Wright       PetscCall(DMGetGlobalVector(projection->dm, &DiffFlux));
5088c85b835SJames Wright       PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, DiffFlux, projection->l2_rhs_ctx));
5098c85b835SJames Wright       PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_rhs_view"));
5108c85b835SJames Wright 
5118c85b835SJames Wright       PetscCall(KSPSolve(projection->ksp, DiffFlux, DiffFlux));
5128c85b835SJames Wright       PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_view"));
5138c85b835SJames Wright 
5148c85b835SJames Wright       PetscCall(ApplyCeedOperatorGlobalToLocal(DiffFlux, NULL, diff_flux_proj->calc_div_diff_flux));
5158c85b835SJames Wright       PetscCall(DMRestoreGlobalVector(projection->dm, &DiffFlux));
5168c85b835SJames Wright     } break;
5178c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_NONE:
5188c85b835SJames Wright       SETERRQ(PetscObjectComm((PetscObject)projection->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
5198c85b835SJames Wright               DivDiffFluxProjectionMethods[diff_flux_proj->method]);
5208c85b835SJames Wright       break;
5218c85b835SJames Wright   }
5228c85b835SJames Wright   PetscCall(PetscLogEventEnd(FLUIDS_DivDiffFluxProjection, Q_loc, 0, 0, 0));
5238c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
5248c85b835SJames Wright }
5258c85b835SJames Wright 
5268c85b835SJames Wright /**
5278c85b835SJames Wright   @brief Destroy `DivDiffFluxProjectionData` object
5288c85b835SJames Wright 
5298c85b835SJames Wright   @param[in,out] diff_flux_proj Object to destroy
5308c85b835SJames Wright **/
5318c85b835SJames Wright PetscErrorCode DivDiffFluxProjectionDataDestroy(DivDiffFluxProjectionData diff_flux_proj) {
5328c85b835SJames Wright   PetscFunctionBeginUser;
5338c85b835SJames Wright   if (diff_flux_proj == NULL) PetscFunctionReturn(PETSC_SUCCESS);
5348c85b835SJames Wright   PetscCall(NodalProjectionDataDestroy(diff_flux_proj->projection));
5358c85b835SJames Wright   PetscCall(OperatorApplyContextDestroy(diff_flux_proj->calc_div_diff_flux));
5368c85b835SJames Wright   if (diff_flux_proj->ceed_vec_has_array) {
5378c85b835SJames Wright     PetscCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->DivDiffFlux_loc));
5388c85b835SJames Wright     diff_flux_proj->ceed_vec_has_array = PETSC_FALSE;
5398c85b835SJames Wright   }
5408c85b835SJames Wright   PetscCall(CeedVectorDestroy(&diff_flux_proj->div_diff_flux_ceed));
5418c85b835SJames Wright   PetscCall(VecDestroy(&diff_flux_proj->DivDiffFlux_loc));
5428c85b835SJames Wright   PetscCall(PetscFree(diff_flux_proj));
5438c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
5448c85b835SJames Wright }
545