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 /** 90*0880fbb6SJames Wright @brief Return the objects required for the Divergence of Diffusive flux to be read by a `CeedOperator` 91*0880fbb6SJames Wright 92*0880fbb6SJames Wright @param[in] diff_flux_proj Projection object 93*0880fbb6SJames Wright @param[out] elem_restr Element restriction for the divergence of diffusive flux, or `NULL` 94*0880fbb6SJames Wright @param[out] basis Basis for the divergence of diffusive flux, or `NULL` 95*0880fbb6SJames Wright @param[out] vector Vector for the divergence of diffusive flux, or `NULL` 96*0880fbb6SJames Wright @param[out] eval_mode Eval mode for the divergence of diffusive flux, or `NULL` 97*0880fbb6SJames Wright **/ 98*0880fbb6SJames Wright PetscErrorCode DivDiffFluxProjectionGetOperatorFieldData(DivDiffFluxProjectionData diff_flux_proj, CeedElemRestriction *elem_restr, CeedBasis *basis, 99*0880fbb6SJames Wright CeedVector *vector, CeedEvalMode *eval_mode) { 100*0880fbb6SJames Wright Ceed ceed = CeedVectorReturnCeed(diff_flux_proj->div_diff_flux_ceed); 101*0880fbb6SJames Wright 102*0880fbb6SJames Wright PetscFunctionBeginUser; 103*0880fbb6SJames Wright if (elem_restr) PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(diff_flux_proj->elem_restr_div_diff_flux, elem_restr)); 104*0880fbb6SJames Wright if (basis) PetscCallCeed(ceed, CeedBasisReferenceCopy(diff_flux_proj->basis_div_diff_flux, basis)); 105*0880fbb6SJames Wright if (vector) PetscCallCeed(ceed, CeedVectorReferenceCopy(diff_flux_proj->div_diff_flux_ceed, vector)); 106*0880fbb6SJames Wright if (eval_mode) *eval_mode = diff_flux_proj->eval_mode_div_diff_flux; 107*0880fbb6SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 108*0880fbb6SJames Wright } 109*0880fbb6SJames Wright 110*0880fbb6SJames Wright /** 1118c85b835SJames Wright @brief Setup direct projection of divergence of diffusive flux 1128c85b835SJames Wright 1138c85b835SJames Wright @param[in] ceed `Ceed` context 1148c85b835SJames Wright @param[in] user `User` context 1158c85b835SJames Wright @param[in] ceed_data `CeedData` context 1168c85b835SJames Wright @param[in] problem `ProblemData` context 1178c85b835SJames Wright **/ 1188c85b835SJames Wright static PetscErrorCode DivDiffFluxProjectionSetup_Direct(Ceed ceed, User user, CeedData ceed_data, ProblemData problem) { 1198c85b835SJames Wright DivDiffFluxProjectionData diff_flux_proj = user->diff_flux_proj; 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"); 12936038bbcSJames Wright PetscCall(diff_flux_proj->CreateRHSOperator_Direct(user, ceed_data, 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; 13236038bbcSJames Wright PetscCall(OperatorApplyContextCreate(user->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; 140*0880fbb6SJames Wright CeedBasis basis_diff_flux = NULL; 141*0880fbb6SJames Wright CeedElemRestriction elem_restr_diff_flux_volume = NULL, elem_restr_qd; 142*0880fbb6SJames Wright CeedVector q_data; 143*0880fbb6SJames Wright CeedInt q_data_size; 144*0880fbb6SJames Wright PetscInt label_value = 0; 145*0880fbb6SJames Wright DMLabel domain_label = NULL; 146*0880fbb6SJames Wright 147*0880fbb6SJames Wright PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_diff_flux_volume, &basis_diff_flux, NULL, NULL)); 148*0880fbb6SJames Wright PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, 149*0880fbb6SJames Wright &elem_restr_qd, &q_data, &q_data_size)); 1508c85b835SJames Wright 1518c85b835SJames Wright PetscCall(CreateMassQFunction(ceed, projection->num_comp, q_data_size, &qf_mass)); 1528c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass)); 1538c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE)); 1548c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 1558c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_diff_flux_volume, basis_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 } 174*0880fbb6SJames Wright 175*0880fbb6SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_volume)); 176*0880fbb6SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 177*0880fbb6SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 178*0880fbb6SJames 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 1888c85b835SJames Wright @param[in] ceed `Ceed` context 1898c85b835SJames Wright @param[in,out] user `User` context 1908c85b835SJames Wright @param[in] ceed_data `CeedData` context 1918c85b835SJames Wright @param[in] problem `ProblemData` context 1928c85b835SJames Wright **/ 1938c85b835SJames Wright static PetscErrorCode DivDiffFluxProjectionSetup_Indirect(Ceed ceed, User user, CeedData ceed_data, ProblemData problem) { 1948c85b835SJames Wright DivDiffFluxProjectionData diff_flux_proj = user->diff_flux_proj; 1958c85b835SJames Wright NodalProjectionData projection = diff_flux_proj->projection; 1968c85b835SJames Wright CeedBasis basis_diff_flux; 1978c85b835SJames Wright CeedElemRestriction elem_restr_diff_flux, elem_restr_qd; 1988c85b835SJames Wright CeedVector q_data; 199*0880fbb6SJames Wright CeedInt q_data_size; 20036038bbcSJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)projection->dm); 2018c85b835SJames Wright 2028c85b835SJames Wright PetscFunctionBeginUser; 203*0880fbb6SJames Wright { 204*0880fbb6SJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 205*0880fbb6SJames Wright DMLabel domain_label = NULL; 2068c85b835SJames Wright 2078c85b835SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, domain_label, label_value, height, dm_field, &elem_restr_diff_flux)); 2088c85b835SJames Wright PetscCall(CreateBasisFromPlex(ceed, projection->dm, domain_label, label_value, height, dm_field, &basis_diff_flux)); 209*0880fbb6SJames Wright PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, 210*0880fbb6SJames Wright &elem_restr_qd, &q_data, &q_data_size)); 211*0880fbb6SJames Wright } 2128c85b835SJames Wright 21336038bbcSJames Wright { 2148c85b835SJames Wright CeedOperator op_rhs; 2158c85b835SJames Wright 21636038bbcSJames Wright PetscCheck(diff_flux_proj->CreateRHSOperator_Indirect, comm, PETSC_ERR_ARG_WRONGSTATE, 21736038bbcSJames Wright "Must define CreateRHSOperator_Indirect to use indirect div_diff_flux projection"); 21836038bbcSJames Wright PetscCall(diff_flux_proj->CreateRHSOperator_Indirect(user, ceed_data, diff_flux_proj, &op_rhs)); 2198c85b835SJames Wright PetscCall(OperatorApplyContextCreate(user->dm, projection->dm, ceed, op_rhs, NULL, NULL, NULL, NULL, &projection->l2_rhs_ctx)); 2208c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs)); 2218c85b835SJames Wright } 2228c85b835SJames Wright 2238c85b835SJames Wright { // -- Build Mass operator 2248c85b835SJames Wright CeedQFunction qf_mass; 2258c85b835SJames Wright CeedOperator op_mass; 2268c85b835SJames Wright 2278c85b835SJames Wright PetscCall(CreateMassQFunction(ceed, projection->num_comp, q_data_size, &qf_mass)); 2288c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass)); 2298c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE)); 2308c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 2318c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE)); 2328c85b835SJames Wright 2338c85b835SJames Wright { // -- Setup KSP for L^2 projection 2348c85b835SJames Wright Mat mat_mass; 2358c85b835SJames Wright 2368c85b835SJames Wright PetscCall(MatCreateCeed(projection->dm, projection->dm, op_mass, NULL, &mat_mass)); 2378c85b835SJames Wright 2388c85b835SJames Wright PetscCall(KSPCreate(comm, &projection->ksp)); 2398c85b835SJames Wright PetscCall(KSPSetOptionsPrefix(projection->ksp, "div_diff_flux_projection_")); 2408c85b835SJames Wright { // lumped by default 2418c85b835SJames Wright PC pc; 2428c85b835SJames Wright PetscCall(KSPGetPC(projection->ksp, &pc)); 2438c85b835SJames Wright PetscCall(PCSetType(pc, PCJACOBI)); 2448c85b835SJames Wright PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 2458c85b835SJames Wright PetscCall(KSPSetType(projection->ksp, KSPPREONLY)); 2468c85b835SJames Wright } 2478c85b835SJames Wright PetscCall(KSPSetFromOptions_WithMatCeed(projection->ksp, mat_mass)); 2488c85b835SJames Wright PetscCall(MatDestroy(&mat_mass)); 2498c85b835SJames Wright } 2508c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 2518c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 2528c85b835SJames Wright } 2538c85b835SJames Wright 2548c85b835SJames Wright { // Create OperatorApplyContext to calculate divergence at quadrature points 2558c85b835SJames Wright CeedQFunction qf_calc_divergence; 2568c85b835SJames Wright CeedOperator op_calc_divergence; 257*0880fbb6SJames Wright CeedElemRestriction elem_restr_div_diff_flux = NULL; 258*0880fbb6SJames Wright PetscInt dim; 2598c85b835SJames Wright 260*0880fbb6SJames Wright PetscCall(DMGetDimension(projection->dm, &dim)); 261*0880fbb6SJames Wright PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_div_diff_flux, NULL, NULL, NULL)); 2628c85b835SJames Wright 2638c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux3D_4, ComputeDivDiffusiveFlux3D_4_loc, &qf_calc_divergence)); 2648c85b835SJames Wright 2658c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_calc_divergence, "Grad F_diff", projection->num_comp * dim, CEED_EVAL_GRAD)); 2668c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_calc_divergence, "qdata", q_data_size, CEED_EVAL_NONE)); 2678c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_calc_divergence, "Div F_diff", 4, CEED_EVAL_NONE)); 2688c85b835SJames Wright 2698c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_calc_divergence, NULL, NULL, &op_calc_divergence)); 2708c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "Grad F_diff", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE)); 2718c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 2728c85b835SJames Wright PetscCallCeed( 2738c85b835SJames Wright ceed, CeedOperatorSetField(op_calc_divergence, "Div F_diff", elem_restr_div_diff_flux, CEED_BASIS_NONE, diff_flux_proj->div_diff_flux_ceed)); 2748c85b835SJames Wright 27536038bbcSJames Wright PetscCall( 27636038bbcSJames Wright OperatorApplyContextCreate(projection->dm, NULL, ceed, op_calc_divergence, NULL, NULL, NULL, NULL, &diff_flux_proj->calc_div_diff_flux)); 2778c85b835SJames Wright 278*0880fbb6SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_div_diff_flux)); 2798c85b835SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_calc_divergence)); 2808c85b835SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_calc_divergence)); 2818c85b835SJames Wright } 2828c85b835SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 2838c85b835SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 2848c85b835SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 2858c85b835SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux)); 2868c85b835SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 2878c85b835SJames Wright } 2888c85b835SJames Wright 2898c85b835SJames Wright /** 2908c85b835SJames Wright @brief Setup projection of divergence of diffusive flux 2918c85b835SJames Wright 2928c85b835SJames Wright @param[in] ceed `Ceed` context 2938c85b835SJames Wright @param[in,out] user `User` context 2948c85b835SJames Wright @param[in] ceed_data `CeedData` context 2958c85b835SJames Wright @param[in] problem `ProblemData` context 2968c85b835SJames Wright **/ 2978c85b835SJames Wright PetscErrorCode DivDiffFluxProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem) { 2988c85b835SJames Wright PetscFunctionBeginUser; 2998c85b835SJames Wright switch (user->app_ctx->divFdiffproj_method) { 3008c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_DIRECT: 3018c85b835SJames Wright PetscCall(DivDiffFluxProjectionSetup_Direct(ceed, user, ceed_data, problem)); 3028c85b835SJames Wright break; 3038c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_INDIRECT: 3048c85b835SJames Wright PetscCall(DivDiffFluxProjectionSetup_Indirect(ceed, user, ceed_data, problem)); 3058c85b835SJames Wright break; 3068c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_NONE: 3078c85b835SJames Wright SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s", 3088c85b835SJames Wright DivDiffFluxProjectionMethods[user->app_ctx->divFdiffproj_method]); 3098c85b835SJames Wright break; 3108c85b835SJames Wright } 3118c85b835SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3128c85b835SJames Wright } 3138c85b835SJames Wright 3148c85b835SJames Wright /** 3158c85b835SJames Wright @brief Project the divergence of diffusive flux 3168c85b835SJames Wright 3178c85b835SJames Wright This implicitly sets the `CeedVector` input (`div_diff_flux_ceed`) to the divergence of diffusive flux. 3188c85b835SJames Wright 3198c85b835SJames Wright @param[in] diff_flux_proj `NodalProjectionData` for the projection 3208c85b835SJames Wright @param[in] Q_loc Localized solution vector 3218c85b835SJames Wright **/ 32236038bbcSJames Wright PetscErrorCode DivDiffFluxProjectionApply(DivDiffFluxProjectionData diff_flux_proj, Vec Q_loc) { 3238c85b835SJames Wright NodalProjectionData projection = diff_flux_proj->projection; 3248c85b835SJames Wright 3258c85b835SJames Wright PetscFunctionBeginUser; 3268c85b835SJames Wright PetscCall(PetscLogEventBegin(FLUIDS_DivDiffFluxProjection, Q_loc, 0, 0, 0)); 3278c85b835SJames Wright switch (diff_flux_proj->method) { 3288c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_DIRECT: { 3298c85b835SJames Wright Vec DivDiffFlux; 3308c85b835SJames Wright 3318c85b835SJames Wright PetscCall(DMGetGlobalVector(projection->dm, &DivDiffFlux)); 3328c85b835SJames Wright if (diff_flux_proj->ceed_vec_has_array) { 3338c85b835SJames Wright PetscCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->DivDiffFlux_loc)); 3348c85b835SJames Wright diff_flux_proj->ceed_vec_has_array = PETSC_FALSE; 3358c85b835SJames Wright } 3368c85b835SJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, DivDiffFlux, projection->l2_rhs_ctx)); 3378c85b835SJames Wright PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_rhs_view")); 3388c85b835SJames Wright 3398c85b835SJames Wright PetscCall(KSPSolve(projection->ksp, DivDiffFlux, DivDiffFlux)); 3408c85b835SJames Wright PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_view")); 3418c85b835SJames Wright 3428c85b835SJames Wright PetscCall(DMGlobalToLocal(projection->dm, DivDiffFlux, INSERT_VALUES, diff_flux_proj->DivDiffFlux_loc)); 3438c85b835SJames Wright PetscCall(VecReadPetscToCeed(diff_flux_proj->DivDiffFlux_loc, &diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->div_diff_flux_ceed)); 3448c85b835SJames Wright diff_flux_proj->ceed_vec_has_array = PETSC_TRUE; 3458c85b835SJames Wright 3468c85b835SJames Wright PetscCall(DMRestoreGlobalVector(projection->dm, &DivDiffFlux)); 3478c85b835SJames Wright break; 3488c85b835SJames Wright } 3498c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_INDIRECT: { 3508c85b835SJames Wright Vec DiffFlux; 3518c85b835SJames Wright 3528c85b835SJames Wright PetscCall(DMGetGlobalVector(projection->dm, &DiffFlux)); 3538c85b835SJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, DiffFlux, projection->l2_rhs_ctx)); 3548c85b835SJames Wright PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_rhs_view")); 3558c85b835SJames Wright 3568c85b835SJames Wright PetscCall(KSPSolve(projection->ksp, DiffFlux, DiffFlux)); 3578c85b835SJames Wright PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_view")); 3588c85b835SJames Wright 3598c85b835SJames Wright PetscCall(ApplyCeedOperatorGlobalToLocal(DiffFlux, NULL, diff_flux_proj->calc_div_diff_flux)); 3608c85b835SJames Wright PetscCall(DMRestoreGlobalVector(projection->dm, &DiffFlux)); 3618c85b835SJames Wright } break; 3628c85b835SJames Wright case DIV_DIFF_FLUX_PROJ_NONE: 3638c85b835SJames Wright SETERRQ(PetscObjectComm((PetscObject)projection->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s", 3648c85b835SJames Wright DivDiffFluxProjectionMethods[diff_flux_proj->method]); 3658c85b835SJames Wright break; 3668c85b835SJames Wright } 3678c85b835SJames Wright PetscCall(PetscLogEventEnd(FLUIDS_DivDiffFluxProjection, Q_loc, 0, 0, 0)); 3688c85b835SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3698c85b835SJames Wright } 3708c85b835SJames Wright 3718c85b835SJames Wright /** 3728c85b835SJames Wright @brief Destroy `DivDiffFluxProjectionData` object 3738c85b835SJames Wright 3748c85b835SJames Wright @param[in,out] diff_flux_proj Object to destroy 3758c85b835SJames Wright **/ 3768c85b835SJames Wright PetscErrorCode DivDiffFluxProjectionDataDestroy(DivDiffFluxProjectionData diff_flux_proj) { 3778c85b835SJames Wright PetscFunctionBeginUser; 3788c85b835SJames Wright if (diff_flux_proj == NULL) PetscFunctionReturn(PETSC_SUCCESS); 379*0880fbb6SJames Wright Ceed ceed = CeedVectorReturnCeed(diff_flux_proj->div_diff_flux_ceed); 380*0880fbb6SJames Wright 3818c85b835SJames Wright PetscCall(NodalProjectionDataDestroy(diff_flux_proj->projection)); 3828c85b835SJames Wright PetscCall(OperatorApplyContextDestroy(diff_flux_proj->calc_div_diff_flux)); 3838c85b835SJames Wright if (diff_flux_proj->ceed_vec_has_array) { 3848c85b835SJames Wright PetscCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->DivDiffFlux_loc)); 3858c85b835SJames Wright diff_flux_proj->ceed_vec_has_array = PETSC_FALSE; 3868c85b835SJames Wright } 387*0880fbb6SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&diff_flux_proj->div_diff_flux_ceed)); 388*0880fbb6SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&diff_flux_proj->elem_restr_div_diff_flux)); 389*0880fbb6SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&diff_flux_proj->basis_div_diff_flux)); 3908c85b835SJames Wright PetscCall(VecDestroy(&diff_flux_proj->DivDiffFlux_loc)); 3918c85b835SJames Wright PetscCall(PetscFree(diff_flux_proj)); 3928c85b835SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3938c85b835SJames Wright } 394