xref: /honee/src/diff_flux_projection.c (revision 8561fee26e1c8104caa6d917d7ecea939d705f29)
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 /**
1336038bbcSJames Wright   @brief Create `DivDiffFluxProjectionData` for solution DM in `user`
148c85b835SJames Wright 
158c85b835SJames Wright   @param[in]  user                `User` 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 **/
1936038bbcSJames Wright PetscErrorCode DivDiffFluxProjectionCreate(User user, PetscInt num_diff_flux_comps, DivDiffFluxProjectionData *pdiff_flux_proj) {
2036038bbcSJames Wright   PetscInt                  label_value = 0, height = 0, dm_field = 0, dim, degree = user->app_ctx->degree, q_extra = user->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;
2636038bbcSJames Wright   if (user->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;
328c85b835SJames Wright   PetscCall(PetscNew(&user->diff_flux_proj->projection));
338c85b835SJames Wright   projection                          = user->diff_flux_proj->projection;
348c85b835SJames Wright   diff_flux_proj->method              = user->app_ctx->divFdiffproj_method;
3536038bbcSJames Wright   diff_flux_proj->num_diff_flux_comps = num_diff_flux_comps;
368c85b835SJames Wright 
378c85b835SJames Wright   PetscCall(DMClone(user->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 
4636038bbcSJames Wright       PetscCall(DMPlexCeedElemRestrictionCreate(user->ceed, projection->dm, domain_label, label_value, height, dm_field,
4736038bbcSJames Wright                                                 &diff_flux_proj->elem_restr_div_diff_flux));
4836038bbcSJames Wright       PetscCallCeed(user->ceed, CeedElemRestrictionCreateVector(diff_flux_proj->elem_restr_div_diff_flux, &diff_flux_proj->div_diff_flux_ceed, NULL));
4936038bbcSJames Wright       PetscCall(CreateBasisFromPlex(user->ceed, projection->dm, domain_label, label_value, height, dm_field, &diff_flux_proj->basis_div_diff_flux));
5036038bbcSJames Wright       diff_flux_proj->eval_mode_div_diff_flux = CEED_EVAL_INTERP;
511af555e8SJames Wright 
521af555e8SJames Wright       {  // Create face labels on projection->dm for boundary integrals
531af555e8SJames Wright         DMLabel  face_sets_label;
541af555e8SJames Wright         PetscInt num_face_set_values, *face_set_values;
551af555e8SJames Wright 
561af555e8SJames Wright         PetscCall(DMGetLabel(user->dm, "Face Sets", &face_sets_label));
571af555e8SJames Wright         PetscCall(DMLabelCreateGlobalValueArray(user->dm, face_sets_label, &num_face_set_values, &face_set_values));
581af555e8SJames Wright         for (PetscInt f = 0; f < num_face_set_values; f++) {
591af555e8SJames Wright           DMLabel face_orientation_label;
601af555e8SJames Wright           char   *face_orientation_label_name;
611af555e8SJames Wright 
621af555e8SJames Wright           PetscCall(DMPlexCreateFaceLabel(user->dm, face_set_values[f], &face_orientation_label_name));
631af555e8SJames Wright           PetscCall(DMGetLabel(user->dm, face_orientation_label_name, &face_orientation_label));
641af555e8SJames Wright           PetscCall(DMAddLabel(projection->dm, face_orientation_label));
651af555e8SJames Wright           PetscCall(PetscFree(face_orientation_label_name));
661af555e8SJames Wright         }
6736038bbcSJames Wright         PetscCall(PetscFree(face_set_values));
681af555e8SJames Wright       }
698c85b835SJames Wright     } break;
708c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_INDIRECT: {
718c85b835SJames Wright       projection->num_comp = diff_flux_proj->num_diff_flux_comps * dim;
728c85b835SJames Wright       PetscCall(PetscObjectSetName((PetscObject)projection->dm, "DiffFluxProj"));
7336038bbcSJames Wright       PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &projection->num_comp, projection->dm));
748c85b835SJames Wright 
758c85b835SJames Wright       PetscCall(DMPlexCeedElemRestrictionQDataCreate(user->ceed, projection->dm, domain_label, label_value, height,
7636038bbcSJames Wright                                                      diff_flux_proj->num_diff_flux_comps, &diff_flux_proj->elem_restr_div_diff_flux));
7736038bbcSJames Wright       PetscCallCeed(user->ceed, CeedElemRestrictionCreateVector(diff_flux_proj->elem_restr_div_diff_flux, &diff_flux_proj->div_diff_flux_ceed, NULL));
7836038bbcSJames Wright       diff_flux_proj->basis_div_diff_flux     = CEED_BASIS_NONE;
7936038bbcSJames Wright       diff_flux_proj->eval_mode_div_diff_flux = CEED_EVAL_NONE;
808c85b835SJames Wright     } break;
818c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_NONE:
828c85b835SJames Wright       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
838c85b835SJames Wright               DivDiffFluxProjectionMethods[user->app_ctx->divFdiffproj_method]);
848c85b835SJames Wright       break;
858c85b835SJames Wright   }
868c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
878c85b835SJames Wright };
888c85b835SJames Wright 
898c85b835SJames Wright /**
900880fbb6SJames Wright   @brief Return the objects required for the Divergence of Diffusive flux to be read by a `CeedOperator`
910880fbb6SJames Wright 
920880fbb6SJames Wright   @param[in]  diff_flux_proj Projection object
930880fbb6SJames Wright   @param[out] elem_restr     Element restriction for the divergence of diffusive flux, or `NULL`
940880fbb6SJames Wright   @param[out] basis          Basis for the divergence of diffusive flux, or `NULL`
950880fbb6SJames Wright   @param[out] vector         Vector for the divergence of diffusive flux, or `NULL`
960880fbb6SJames Wright   @param[out] eval_mode      Eval mode for the divergence of diffusive flux, or `NULL`
970880fbb6SJames Wright **/
980880fbb6SJames Wright PetscErrorCode DivDiffFluxProjectionGetOperatorFieldData(DivDiffFluxProjectionData diff_flux_proj, CeedElemRestriction *elem_restr, CeedBasis *basis,
990880fbb6SJames Wright                                                          CeedVector *vector, CeedEvalMode *eval_mode) {
1000880fbb6SJames Wright   Ceed ceed = CeedVectorReturnCeed(diff_flux_proj->div_diff_flux_ceed);
1010880fbb6SJames Wright 
1020880fbb6SJames Wright   PetscFunctionBeginUser;
1030880fbb6SJames Wright   if (elem_restr) PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(diff_flux_proj->elem_restr_div_diff_flux, elem_restr));
1040880fbb6SJames Wright   if (basis) PetscCallCeed(ceed, CeedBasisReferenceCopy(diff_flux_proj->basis_div_diff_flux, basis));
1050880fbb6SJames Wright   if (vector) PetscCallCeed(ceed, CeedVectorReferenceCopy(diff_flux_proj->div_diff_flux_ceed, vector));
1060880fbb6SJames Wright   if (eval_mode) *eval_mode = diff_flux_proj->eval_mode_div_diff_flux;
1070880fbb6SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1080880fbb6SJames Wright }
1090880fbb6SJames Wright 
1100880fbb6SJames Wright /**
1118c85b835SJames Wright   @brief Setup direct projection of divergence of diffusive flux
1128c85b835SJames Wright 
1138c85b835SJames Wright   @param[in]     user           `User` context
1148c85b835SJames Wright   @param[in]     ceed_data      `CeedData` context
115*8561fee2SJames Wright   @param[in,out] diff_flux_proj Flux projection object to setup
1168c85b835SJames Wright **/
117*8561fee2SJames Wright static PetscErrorCode DivDiffFluxProjectionSetup_Direct(User user, CeedData ceed_data, DivDiffFluxProjectionData diff_flux_proj) {
118*8561fee2SJames Wright   Ceed                ceed       = user->ceed;
1198c85b835SJames Wright   NodalProjectionData projection = diff_flux_proj->projection;
12036038bbcSJames Wright   MPI_Comm            comm       = PetscObjectComm((PetscObject)projection->dm);
1218c85b835SJames Wright 
1228c85b835SJames Wright   PetscFunctionBeginUser;
12336038bbcSJames Wright   {  // Create Projection RHS OperatorApplyContext
12436038bbcSJames Wright     CeedOperator op_rhs;
1258c85b835SJames Wright 
12636038bbcSJames Wright     PetscCheck(diff_flux_proj->CreateRHSOperator_Direct, comm, PETSC_ERR_ARG_WRONGSTATE,
12736038bbcSJames Wright                "Must define CreateRHSOperator_Direct to use indirect div_diff_flux projection");
12836038bbcSJames Wright     PetscCall(diff_flux_proj->CreateRHSOperator_Direct(user, ceed_data, diff_flux_proj, &op_rhs));
1298c85b835SJames Wright     PetscCall(DMCreateLocalVector(projection->dm, &diff_flux_proj->DivDiffFlux_loc));
1308c85b835SJames Wright     diff_flux_proj->ceed_vec_has_array = PETSC_FALSE;
13136038bbcSJames Wright     PetscCall(OperatorApplyContextCreate(user->dm, projection->dm, ceed, op_rhs, NULL, NULL, NULL, diff_flux_proj->DivDiffFlux_loc,
13236038bbcSJames Wright                                          &projection->l2_rhs_ctx));
13336038bbcSJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
13436038bbcSJames Wright   }
1358c85b835SJames Wright 
1368c85b835SJames Wright   {  // -- Build Mass operator
1378c85b835SJames Wright     CeedQFunction       qf_mass;
1388c85b835SJames Wright     CeedOperator        op_mass;
139590d9cddSJames Wright     CeedBasis           basis_div_diff_flux             = NULL;
140590d9cddSJames Wright     CeedElemRestriction elem_restr_div_diff_flux_volume = NULL, elem_restr_qd;
1410880fbb6SJames Wright     CeedVector          q_data;
1420880fbb6SJames Wright     CeedInt             q_data_size;
1430880fbb6SJames Wright     PetscInt            label_value  = 0;
1440880fbb6SJames Wright     DMLabel             domain_label = NULL;
1450880fbb6SJames Wright 
146590d9cddSJames Wright     PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_div_diff_flux_volume, &basis_div_diff_flux, NULL, NULL));
1470880fbb6SJames Wright     PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord,
1480880fbb6SJames Wright                        &elem_restr_qd, &q_data, &q_data_size));
1498c85b835SJames Wright 
1508c85b835SJames Wright     PetscCall(CreateMassQFunction(ceed, projection->num_comp, q_data_size, &qf_mass));
1518c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
152590d9cddSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_div_diff_flux_volume, basis_div_diff_flux, CEED_VECTOR_ACTIVE));
1538c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
154590d9cddSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_div_diff_flux_volume, basis_div_diff_flux, CEED_VECTOR_ACTIVE));
1558c85b835SJames Wright 
1568c85b835SJames Wright     {  // -- Setup KSP for L^2 projection
1578c85b835SJames Wright       Mat mat_mass;
1588c85b835SJames Wright 
1598c85b835SJames Wright       PetscCall(MatCreateCeed(projection->dm, projection->dm, op_mass, NULL, &mat_mass));
1608c85b835SJames Wright 
1618c85b835SJames Wright       PetscCall(KSPCreate(comm, &projection->ksp));
1628c85b835SJames Wright       PetscCall(KSPSetOptionsPrefix(projection->ksp, "div_diff_flux_projection_"));
1638c85b835SJames Wright       {  // lumped by default
1648c85b835SJames Wright         PC pc;
1658c85b835SJames Wright         PetscCall(KSPGetPC(projection->ksp, &pc));
1668c85b835SJames Wright         PetscCall(PCSetType(pc, PCJACOBI));
1678c85b835SJames Wright         PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
1688c85b835SJames Wright         PetscCall(KSPSetType(projection->ksp, KSPPREONLY));
1698c85b835SJames Wright       }
1708c85b835SJames Wright       PetscCall(KSPSetFromOptions_WithMatCeed(projection->ksp, mat_mass));
1718c85b835SJames Wright       PetscCall(MatDestroy(&mat_mass));
1728c85b835SJames Wright     }
1730880fbb6SJames Wright 
174590d9cddSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_div_diff_flux_volume));
175590d9cddSJames Wright     PetscCallCeed(ceed, CeedBasisDestroy(&basis_div_diff_flux));
1760880fbb6SJames Wright     PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
1770880fbb6SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
1788c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
1798c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
1808c85b835SJames Wright   }
1818c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1828c85b835SJames Wright }
1838c85b835SJames Wright 
1848c85b835SJames Wright /**
1858c85b835SJames Wright   @brief Setup indirect projection of divergence of diffusive flux
1868c85b835SJames Wright 
187*8561fee2SJames Wright   @param[in]     user           `User` context
1888c85b835SJames Wright   @param[in]     ceed_data      `CeedData` context
189*8561fee2SJames Wright   @param[in,out] diff_flux_proj Flux projection object to setup
1908c85b835SJames Wright **/
191*8561fee2SJames Wright static PetscErrorCode DivDiffFluxProjectionSetup_Indirect(User user, CeedData ceed_data, DivDiffFluxProjectionData diff_flux_proj) {
192*8561fee2SJames Wright   Ceed                ceed       = user->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));
2070880fbb6SJames Wright     PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord,
2080880fbb6SJames Wright                        &elem_restr_qd, &q_data, &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");
21636038bbcSJames Wright     PetscCall(diff_flux_proj->CreateRHSOperator_Indirect(user, ceed_data, diff_flux_proj, &op_rhs));
2178c85b835SJames Wright     PetscCall(OperatorApplyContextCreate(user->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 
2258c85b835SJames Wright     PetscCall(CreateMassQFunction(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
2538c85b835SJames Wright     CeedQFunction       qf_calc_divergence;
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 
2618c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux3D_4, ComputeDivDiffusiveFlux3D_4_loc, &qf_calc_divergence));
2628c85b835SJames Wright 
2638c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_calc_divergence, "Grad F_diff", projection->num_comp * dim, CEED_EVAL_GRAD));
2648c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_calc_divergence, "qdata", q_data_size, CEED_EVAL_NONE));
2658c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_calc_divergence, "Div F_diff", 4, CEED_EVAL_NONE));
2668c85b835SJames Wright 
2678c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_calc_divergence, NULL, NULL, &op_calc_divergence));
2688c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "Grad F_diff", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
2698c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
2708c85b835SJames Wright     PetscCallCeed(
2718c85b835SJames Wright         ceed, CeedOperatorSetField(op_calc_divergence, "Div F_diff", elem_restr_div_diff_flux, CEED_BASIS_NONE, diff_flux_proj->div_diff_flux_ceed));
2728c85b835SJames Wright 
27336038bbcSJames Wright     PetscCall(
27436038bbcSJames Wright         OperatorApplyContextCreate(projection->dm, NULL, ceed, op_calc_divergence, NULL, NULL, NULL, NULL, &diff_flux_proj->calc_div_diff_flux));
2758c85b835SJames Wright 
2760880fbb6SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_div_diff_flux));
2778c85b835SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_calc_divergence));
2788c85b835SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_calc_divergence));
2798c85b835SJames Wright   }
2808c85b835SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
2818c85b835SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
2828c85b835SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
2838c85b835SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux));
2848c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2858c85b835SJames Wright }
2868c85b835SJames Wright 
2878c85b835SJames Wright /**
2888c85b835SJames Wright   @brief Setup projection of divergence of diffusive flux
2898c85b835SJames Wright 
290*8561fee2SJames Wright   @param[in]     user           `User` context
2918c85b835SJames Wright   @param[in]     ceed_data      `CeedData` context
292*8561fee2SJames Wright   @param[in,out] diff_flux_proj Flux projection object to setup
2938c85b835SJames Wright **/
294*8561fee2SJames Wright PetscErrorCode DivDiffFluxProjectionSetup(User user, CeedData ceed_data, DivDiffFluxProjectionData diff_flux_proj) {
2958c85b835SJames Wright   PetscFunctionBeginUser;
2968c85b835SJames Wright   switch (user->app_ctx->divFdiffproj_method) {
2978c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_DIRECT:
298*8561fee2SJames Wright       PetscCall(DivDiffFluxProjectionSetup_Direct(user, ceed_data, diff_flux_proj));
2998c85b835SJames Wright       break;
3008c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_INDIRECT:
301*8561fee2SJames Wright       PetscCall(DivDiffFluxProjectionSetup_Indirect(user, ceed_data, diff_flux_proj));
3028c85b835SJames Wright       break;
3038c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_NONE:
3048c85b835SJames Wright       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
3058c85b835SJames Wright               DivDiffFluxProjectionMethods[user->app_ctx->divFdiffproj_method]);
3068c85b835SJames Wright       break;
3078c85b835SJames Wright   }
3088c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3098c85b835SJames Wright }
3108c85b835SJames Wright 
3118c85b835SJames Wright /**
3128c85b835SJames Wright   @brief Project the divergence of diffusive flux
3138c85b835SJames Wright 
3148c85b835SJames Wright   This implicitly sets the `CeedVector` input (`div_diff_flux_ceed`) to the divergence of diffusive flux.
3158c85b835SJames Wright 
3168c85b835SJames Wright   @param[in]  diff_flux_proj `NodalProjectionData` for the projection
3178c85b835SJames Wright   @param[in]  Q_loc          Localized solution vector
3188c85b835SJames Wright **/
31936038bbcSJames Wright PetscErrorCode DivDiffFluxProjectionApply(DivDiffFluxProjectionData diff_flux_proj, Vec Q_loc) {
3208c85b835SJames Wright   NodalProjectionData projection = diff_flux_proj->projection;
3218c85b835SJames Wright 
3228c85b835SJames Wright   PetscFunctionBeginUser;
3238c85b835SJames Wright   PetscCall(PetscLogEventBegin(FLUIDS_DivDiffFluxProjection, Q_loc, 0, 0, 0));
3248c85b835SJames Wright   switch (diff_flux_proj->method) {
3258c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_DIRECT: {
3268c85b835SJames Wright       Vec DivDiffFlux;
3278c85b835SJames Wright 
3288c85b835SJames Wright       PetscCall(DMGetGlobalVector(projection->dm, &DivDiffFlux));
3298c85b835SJames Wright       if (diff_flux_proj->ceed_vec_has_array) {
3308c85b835SJames Wright         PetscCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->DivDiffFlux_loc));
3318c85b835SJames Wright         diff_flux_proj->ceed_vec_has_array = PETSC_FALSE;
3328c85b835SJames Wright       }
3338c85b835SJames Wright       PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, DivDiffFlux, projection->l2_rhs_ctx));
3348c85b835SJames Wright       PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_rhs_view"));
3358c85b835SJames Wright 
3368c85b835SJames Wright       PetscCall(KSPSolve(projection->ksp, DivDiffFlux, DivDiffFlux));
3378c85b835SJames Wright       PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_view"));
3388c85b835SJames Wright 
3398c85b835SJames Wright       PetscCall(DMGlobalToLocal(projection->dm, DivDiffFlux, INSERT_VALUES, diff_flux_proj->DivDiffFlux_loc));
3408c85b835SJames Wright       PetscCall(VecReadPetscToCeed(diff_flux_proj->DivDiffFlux_loc, &diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->div_diff_flux_ceed));
3418c85b835SJames Wright       diff_flux_proj->ceed_vec_has_array = PETSC_TRUE;
3428c85b835SJames Wright 
3438c85b835SJames Wright       PetscCall(DMRestoreGlobalVector(projection->dm, &DivDiffFlux));
3448c85b835SJames Wright       break;
3458c85b835SJames Wright     }
3468c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_INDIRECT: {
3478c85b835SJames Wright       Vec DiffFlux;
3488c85b835SJames Wright 
3498c85b835SJames Wright       PetscCall(DMGetGlobalVector(projection->dm, &DiffFlux));
3508c85b835SJames Wright       PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, DiffFlux, projection->l2_rhs_ctx));
3518c85b835SJames Wright       PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_rhs_view"));
3528c85b835SJames Wright 
3538c85b835SJames Wright       PetscCall(KSPSolve(projection->ksp, DiffFlux, DiffFlux));
3548c85b835SJames Wright       PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_view"));
3558c85b835SJames Wright 
3568c85b835SJames Wright       PetscCall(ApplyCeedOperatorGlobalToLocal(DiffFlux, NULL, diff_flux_proj->calc_div_diff_flux));
3578c85b835SJames Wright       PetscCall(DMRestoreGlobalVector(projection->dm, &DiffFlux));
3588c85b835SJames Wright     } break;
3598c85b835SJames Wright     case DIV_DIFF_FLUX_PROJ_NONE:
3608c85b835SJames Wright       SETERRQ(PetscObjectComm((PetscObject)projection->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
3618c85b835SJames Wright               DivDiffFluxProjectionMethods[diff_flux_proj->method]);
3628c85b835SJames Wright       break;
3638c85b835SJames Wright   }
3648c85b835SJames Wright   PetscCall(PetscLogEventEnd(FLUIDS_DivDiffFluxProjection, Q_loc, 0, 0, 0));
3658c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3668c85b835SJames Wright }
3678c85b835SJames Wright 
3688c85b835SJames Wright /**
3698c85b835SJames Wright   @brief Destroy `DivDiffFluxProjectionData` object
3708c85b835SJames Wright 
3718c85b835SJames Wright   @param[in,out] diff_flux_proj Object to destroy
3728c85b835SJames Wright **/
3738c85b835SJames Wright PetscErrorCode DivDiffFluxProjectionDataDestroy(DivDiffFluxProjectionData diff_flux_proj) {
3748c85b835SJames Wright   PetscFunctionBeginUser;
3758c85b835SJames Wright   if (diff_flux_proj == NULL) PetscFunctionReturn(PETSC_SUCCESS);
3760880fbb6SJames Wright   Ceed ceed = CeedVectorReturnCeed(diff_flux_proj->div_diff_flux_ceed);
3770880fbb6SJames Wright 
3788c85b835SJames Wright   PetscCall(NodalProjectionDataDestroy(diff_flux_proj->projection));
3798c85b835SJames Wright   PetscCall(OperatorApplyContextDestroy(diff_flux_proj->calc_div_diff_flux));
3808c85b835SJames Wright   if (diff_flux_proj->ceed_vec_has_array) {
3818c85b835SJames Wright     PetscCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->DivDiffFlux_loc));
3828c85b835SJames Wright     diff_flux_proj->ceed_vec_has_array = PETSC_FALSE;
3838c85b835SJames Wright   }
3840880fbb6SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&diff_flux_proj->div_diff_flux_ceed));
3850880fbb6SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&diff_flux_proj->elem_restr_div_diff_flux));
3860880fbb6SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&diff_flux_proj->basis_div_diff_flux));
3878c85b835SJames Wright   PetscCall(VecDestroy(&diff_flux_proj->DivDiffFlux_loc));
3888c85b835SJames Wright   PetscCall(PetscFree(diff_flux_proj));
3898c85b835SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3908c85b835SJames Wright }
391