1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3a515125bSLeila Ghaffari 4a515125bSLeila Ghaffari /// @file 5ea615d4cSJames Wright /// Setup libCEED Operators for HONEE 6a515125bSLeila Ghaffari 7e419654dSJeremy L Thompson #include <ceed.h> 8e419654dSJeremy L Thompson #include <petscdmplex.h> 99ae013d6SJames Wright #include <smartsim.h> 10e419654dSJeremy L Thompson 11149fb536SJames Wright #include <navierstokes.h> 12a515125bSLeila Ghaffari 13a78efa86SJames Wright // @brief Create CeedOperator for unstabilized mass KSP for explicit timestepping 140c373b74SJames Wright static PetscErrorCode CreateKSPMassOperator_Unstabilized(Honee honee, CeedOperator *op_mass) { 150c373b74SJames Wright Ceed ceed = honee->ceed; 16a78efa86SJames Wright CeedInt num_comp_q, q_data_size; 17a78efa86SJames Wright CeedQFunction qf_mass; 1801e19bfaSJames Wright CeedElemRestriction elem_restr_q, elem_restr_qd; 19a78efa86SJames Wright CeedBasis basis_q; 20a78efa86SJames Wright CeedVector q_data; 21a78efa86SJames Wright 22a78efa86SJames Wright PetscFunctionBeginUser; 23a78efa86SJames Wright { // Get restriction and basis from the RHS function 24a78efa86SJames Wright CeedOperator *sub_ops; 2501e19bfaSJames Wright CeedOperatorField op_field; 26a78efa86SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 27a78efa86SJames Wright 28da7f3a07SJames Wright PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_rhs_ctx->op, &sub_ops)); 2901e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &op_field)); 3001e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_q, &basis_q, NULL)); 3101e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &op_field)); 3201e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_qd, NULL, &q_data)); 33a78efa86SJames Wright } 34a78efa86SJames Wright 35a78efa86SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q)); 3601e19bfaSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd, &q_data_size)); 37a78efa86SJames Wright 3864dd23feSJames Wright PetscCall(HoneeMassQFunctionCreate(ceed, num_comp_q, q_data_size, &qf_mass)); 39a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass)); 4071f2ed29SJames Wright PetscCallCeed(ceed, CeedOperatorSetName(*op_mass, "RHS Mass Operator")); 41a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "u", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 4201e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 43a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 44a78efa86SJames Wright 45fff85bd3SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 46fff85bd3SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 47fff85bd3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 4801e19bfaSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 49a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 50a78efa86SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 51a78efa86SJames Wright } 52a78efa86SJames Wright 53a78efa86SJames Wright // @brief Create KSP to solve the inverse mass operator for explicit time stepping schemes 540c373b74SJames Wright static PetscErrorCode CreateKSPMass(Honee honee, ProblemData problem) { 550c373b74SJames Wright Ceed ceed = honee->ceed; 560c373b74SJames Wright DM dm = honee->dm; 57a78efa86SJames Wright CeedOperator op_mass; 58a78efa86SJames Wright 59a78efa86SJames Wright PetscFunctionBeginUser; 600c373b74SJames Wright if (problem->create_mass_operator) PetscCall(problem->create_mass_operator(honee, &op_mass)); 610c373b74SJames Wright else PetscCall(CreateKSPMassOperator_Unstabilized(honee, &op_mass)); 62a78efa86SJames Wright 63a78efa86SJames Wright { // -- Setup KSP for mass operator 64a78efa86SJames Wright Mat mat_mass; 65a78efa86SJames Wright Vec Zeros_loc; 66a78efa86SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)dm); 67a78efa86SJames Wright 68a78efa86SJames Wright PetscCall(DMCreateLocalVector(dm, &Zeros_loc)); 69a78efa86SJames Wright PetscCall(VecZeroEntries(Zeros_loc)); 70000d2032SJeremy L Thompson PetscCall(MatCreateCeed(dm, dm, op_mass, NULL, &mat_mass)); 71a78efa86SJames Wright PetscCall(MatCeedSetLocalVectors(mat_mass, Zeros_loc, NULL)); 72a78efa86SJames Wright 730c373b74SJames Wright PetscCall(KSPCreate(comm, &honee->mass_ksp)); 740c373b74SJames Wright PetscCall(KSPSetOptionsPrefix(honee->mass_ksp, "mass_")); 7571f2ed29SJames Wright PetscCall(PetscObjectSetName((PetscObject)honee->mass_ksp, "Explicit Mass")); 76a78efa86SJames Wright { // lumped by default 77a78efa86SJames Wright PC pc; 780c373b74SJames Wright PetscCall(KSPGetPC(honee->mass_ksp, &pc)); 79a78efa86SJames Wright PetscCall(PCSetType(pc, PCJACOBI)); 80a78efa86SJames Wright PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 810c373b74SJames Wright PetscCall(KSPSetType(honee->mass_ksp, KSPPREONLY)); 82a78efa86SJames Wright } 830c373b74SJames Wright PetscCall(KSPSetFromOptions_WithMatCeed(honee->mass_ksp, mat_mass)); 84a78efa86SJames Wright PetscCall(VecDestroy(&Zeros_loc)); 85a78efa86SJames Wright PetscCall(MatDestroy(&mat_mass)); 86a78efa86SJames Wright } 87a78efa86SJames Wright 88a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 89a78efa86SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 90a78efa86SJames Wright } 91a78efa86SJames Wright 92d3c60affSJames Wright PetscErrorCode SetupLibceed(Ceed ceed, DM dm, Honee honee, AppCtx app_ctx, ProblemData problem) { 93f27c2204SJames Wright const PetscInt num_comp_q = problem->num_components; 94cf8f12c8SJames Wright PetscInt dim, num_comp_x; 95cf8f12c8SJames Wright CeedInt num_comps_jac_data = problem->num_comps_jac_data, q_data_size_vol; 96cf8f12c8SJames Wright CeedElemRestriction elem_restr_jd_i = NULL, elem_restr_qd, elem_restr_q; 97cf8f12c8SJames Wright CeedBasis basis_q; 98be29160dSJames Wright CeedVector jac_data = NULL, q_data; 99bbc90103SJames Wright CeedOperator op_ifunction_vol = NULL, op_rhs_vol = NULL, op_ijacobian_vol = NULL; 100bbc90103SJames Wright 101bbc90103SJames Wright PetscFunctionBeginUser; 10266d54740SJames Wright PetscCall(DMGetDimension(dm, &dim)); 103cf8f12c8SJames Wright PetscCall(DMGetCoordinateNumComps(dm, &num_comp_x)); 10494a7b3d2SKenneth E. Jansen 1058c85b835SJames Wright CeedElemRestriction elem_restr_diff_flux = NULL; 1060880fbb6SJames Wright CeedVector div_diff_flux_ceed = NULL; 1078c85b835SJames Wright CeedBasis basis_diff_flux = NULL; 1088c85b835SJames Wright CeedEvalMode eval_mode_diff_flux = -1; 109bbc90103SJames Wright { // Create bases and element restrictions 110e3db12f8SJames Wright PetscInt height = 0, dm_field = 0; 11167263decSKenneth E. Jansen DM dm_coord; 11267263decSKenneth E. Jansen 113bbc90103SJames Wright PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 114cf8f12c8SJames Wright PetscCall(DMPlexCeedBasisCreate(ceed, dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &basis_q)); 115a515125bSLeila Ghaffari 116cf8f12c8SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, 0, &elem_restr_q)); 1171abc2837SJames Wright if (num_comps_jac_data) { 118e3db12f8SJames Wright PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, num_comps_jac_data, &elem_restr_jd_i)); 11928160fc2SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL)); 12028160fc2SJames Wright } 121ce8bebb6SJames Wright 122cf8f12c8SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_q, &honee->q_ceed, NULL)); 123cf8f12c8SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_q, &honee->q_dot_ceed, NULL)); 124cf8f12c8SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_q, &honee->g_ceed, NULL)); 125a515125bSLeila Ghaffari 126*77efe802SJames Wright PetscCall(DMPlexCeedCoordinateGetField(ceed, dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, &honee->elem_restr_x, &honee->basis_x, 127*77efe802SJames Wright &honee->x_coord)); 128a515125bSLeila Ghaffari 129e3db12f8SJames Wright PetscCall(QDataGet(ceed, dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, &q_data, 130be29160dSJames Wright &q_data_size_vol)); 131ce8bebb6SJames Wright } 132ce8bebb6SJames Wright 1337cea0f50SJames Wright if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) { 1340c373b74SJames Wright PetscCheck(honee->diff_flux_proj, honee->comm, PETSC_ERR_ARG_WRONGSTATE, 1357cea0f50SJames Wright "Divergence of diffusive flux projection requested but object not created"); 1360c373b74SJames Wright PetscCall(DivDiffFluxProjectionGetOperatorFieldData(honee->diff_flux_proj, &elem_restr_diff_flux, &basis_diff_flux, &div_diff_flux_ceed, 1370880fbb6SJames Wright &eval_mode_diff_flux)); 1387cea0f50SJames Wright } 1398c85b835SJames Wright 140ce8bebb6SJames Wright { // -- Create QFunction for ICs 141866f9b4aSJames Wright CeedBasis basis_xc; 142ce8bebb6SJames Wright CeedQFunction qf_ics; 1438f18bb8bSJames Wright CeedOperator op_ics; 144ce8bebb6SJames Wright 145cf8f12c8SJames Wright PetscCallCeed(ceed, CeedBasisCreateProjection(honee->basis_x, basis_q, &basis_xc)); 146e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->ics.qf_func_ptr, problem->ics.qf_loc, &qf_ics)); 147e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ics, problem->ics.qfctx)); 148ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ics, 0)); 149ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "x", num_comp_x, CEED_EVAL_INTERP)); 150ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "dx", num_comp_x * dim, CEED_EVAL_GRAD)); 151ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ics, "q0", num_comp_q, CEED_EVAL_NONE)); 152ce8bebb6SJames Wright 153ce8bebb6SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics)); 154e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "x", honee->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE)); 155e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "dx", honee->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE)); 156cf8f12c8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "q0", elem_restr_q, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 1570c373b74SJames Wright PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ics, "evaluation time", &honee->phys->ics_time_label)); 158e3663b90SJames Wright PetscCall(OperatorApplyContextCreate(NULL, dm, honee->ceed, op_ics, honee->x_coord, NULL, NULL, honee->Q_loc, &honee->op_ics_ctx)); 159a515125bSLeila Ghaffari 160866f9b4aSJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_xc)); 161ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ics)); 162ce8bebb6SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_ics)); 163ce8bebb6SJames Wright } 164ce8bebb6SJames Wright 165e07531f7SJames Wright if (problem->apply_vol_rhs.qf_func_ptr) { 166ce8bebb6SJames Wright CeedQFunction qf_rhs_vol; 167ce8bebb6SJames Wright 168e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qf_func_ptr, problem->apply_vol_rhs.qf_loc, &qf_rhs_vol)); 169e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_vol, problem->apply_vol_rhs.qfctx)); 170ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_rhs_vol, 0)); 171ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP)); 172ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 17389de3142SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE)); 174ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP)); 1755f952e8dSJames Wright if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) 1760c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "div F_diff", honee->diff_flux_proj->num_diff_flux_comps, eval_mode_diff_flux)); 177ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP)); 178ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD)); 179ce8bebb6SJames Wright 180bbc90103SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_vol, NULL, NULL, &op_rhs_vol)); 181cf8f12c8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 182cf8f12c8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 183be29160dSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 184e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "x", honee->elem_restr_x, honee->basis_x, honee->x_coord)); 1855f952e8dSJames Wright if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) 1865f952e8dSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "div F_diff", elem_restr_diff_flux, basis_diff_flux, div_diff_flux_ceed)); 187cf8f12c8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 188cf8f12c8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 189ce8bebb6SJames Wright 190ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_vol)); 191a515125bSLeila Ghaffari } 192a515125bSLeila Ghaffari 193e07531f7SJames Wright if (problem->apply_vol_ifunction.qf_func_ptr) { 194ce8bebb6SJames Wright CeedQFunction qf_ifunction_vol; 195ce8bebb6SJames Wright 1969eadbee4SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qf_func_ptr, problem->apply_vol_ifunction.qf_loc, 1979eadbee4SJames Wright &qf_ifunction_vol)); 198e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ifunction_vol, problem->apply_vol_ifunction.qfctx)); 199ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ifunction_vol, 0)); 200ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP)); 201ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 202ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP)); 20389de3142SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE)); 204ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP)); 2050880fbb6SJames Wright if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) 2060c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "div F_diff", honee->diff_flux_proj->num_diff_flux_comps, eval_mode_diff_flux)); 207ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP)); 208ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD)); 2091abc2837SJames Wright if (num_comps_jac_data) PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "jac_data", num_comps_jac_data, CEED_EVAL_NONE)); 210ce8bebb6SJames Wright 211bbc90103SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ifunction_vol, NULL, NULL, &op_ifunction_vol)); 212cf8f12c8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 213cf8f12c8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 214cf8f12c8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "q dot", elem_restr_q, basis_q, honee->q_dot_ceed)); 215be29160dSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 216e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "x", honee->elem_restr_x, honee->basis_x, honee->x_coord)); 2178c85b835SJames Wright if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) 2180880fbb6SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "div F_diff", elem_restr_diff_flux, basis_diff_flux, div_diff_flux_ceed)); 219cf8f12c8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 220cf8f12c8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 2211abc2837SJames Wright if (num_comps_jac_data) PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data)); 222752f40e3SJed Brown 223ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ifunction_vol)); 224a515125bSLeila Ghaffari } 225a515125bSLeila Ghaffari 226e07531f7SJames Wright if (problem->apply_vol_ijacobian.qf_func_ptr) { 227ce8bebb6SJames Wright CeedQFunction qf_ijacobian_vol; 228ce8bebb6SJames Wright 2299eadbee4SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qf_func_ptr, problem->apply_vol_ijacobian.qf_loc, 2309eadbee4SJames Wright &qf_ijacobian_vol)); 231e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfctx)); 232ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ijacobian_vol, 0)); 233ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP)); 234ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD)); 23589de3142SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE)); 2361abc2837SJames Wright if (num_comps_jac_data) PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", num_comps_jac_data, CEED_EVAL_NONE)); 237ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP)); 238ce8bebb6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD)); 239ce8bebb6SJames Wright 240bbc90103SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op_ijacobian_vol)); 241cf8f12c8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "dq", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 242cf8f12c8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "Grad_dq", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 243be29160dSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 2441abc2837SJames Wright if (num_comps_jac_data) PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data)); 245cf8f12c8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 246cf8f12c8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 247ce8bebb6SJames Wright 248b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ijacobian_vol)); 249f0b65372SJed Brown } 250f0b65372SJed Brown 251a515125bSLeila Ghaffari // -- Create and apply CEED Composite Operator for the entire domain 2520c373b74SJames Wright if (!honee->phys->implicit) { // RHS 253da5fe0e4SJames Wright CeedOperator op_rhs; 254866f9b4aSJames Wright 255da7f3a07SJames Wright PetscCallCeed(ceed, CeedOperatorCreateComposite(ceed, &op_rhs)); 256da7f3a07SJames Wright PetscCallCeed(ceed, CeedOperatorCompositeAddSub(op_rhs, op_rhs_vol)); 257d6cac220SJames Wright for (PetscCount b = 0; b < problem->num_bc_defs; b++) PetscCall(BCDefinitionAddOperators(problem->bc_defs[b], op_rhs, NULL)); 258866f9b4aSJames Wright 2590c373b74SJames Wright PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_rhs, honee->q_ceed, honee->g_ceed, honee->Q_loc, NULL, &honee->op_rhs_ctx)); 260866f9b4aSJames Wright 261866f9b4aSJames Wright // ----- Get Context Labels for Operator 2620c373b74SJames Wright PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "solution time", &honee->phys->solution_time_label)); 2630c373b74SJames Wright PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "timestep size", &honee->phys->timestep_size_label)); 264866f9b4aSJames Wright 265b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs)); 2660c373b74SJames Wright PetscCall(CreateKSPMass(honee, problem)); 2670c373b74SJames Wright PetscCheck(app_ctx->sgs_model_type == SGS_MODEL_NONE, honee->comm, PETSC_ERR_SUP, "SGS modeling not implemented for explicit timestepping"); 268a515125bSLeila Ghaffari } else { // IFunction 269ebfabadfSJames Wright CeedOperator op_ijacobian = NULL; 270ebfabadfSJames Wright 271866f9b4aSJames Wright // Create Composite Operaters 272da7f3a07SJames Wright PetscCallCeed(ceed, CeedOperatorCreateComposite(ceed, &honee->op_ifunction)); 273da7f3a07SJames Wright PetscCallCeed(ceed, CeedOperatorCompositeAddSub(honee->op_ifunction, op_ifunction_vol)); 274866f9b4aSJames Wright if (op_ijacobian_vol) { 275da7f3a07SJames Wright PetscCallCeed(ceed, CeedOperatorCreateComposite(ceed, &op_ijacobian)); 276da7f3a07SJames Wright PetscCallCeed(ceed, CeedOperatorCompositeAddSub(op_ijacobian, op_ijacobian_vol)); 277866f9b4aSJames Wright } 278d6cac220SJames Wright for (PetscCount b = 0; b < problem->num_bc_defs; b++) PetscCall(BCDefinitionAddOperators(problem->bc_defs[b], honee->op_ifunction, op_ijacobian)); 279866f9b4aSJames Wright 280866f9b4aSJames Wright // ----- Get Context Labels for Operator 2810c373b74SJames Wright PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(honee->op_ifunction, "solution time", &honee->phys->solution_time_label)); 2820c373b74SJames Wright PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(honee->op_ifunction, "timestep size", &honee->phys->timestep_size_label)); 283866f9b4aSJames Wright 284ebfabadfSJames Wright if (op_ijacobian) { 2850c373b74SJames Wright PetscCall(MatCreateCeed(honee->dm, honee->dm, op_ijacobian, NULL, &honee->mat_ijacobian)); 2860c373b74SJames Wright PetscCall(MatCeedSetLocalVectors(honee->mat_ijacobian, honee->Q_dot_loc, NULL)); 287ebfabadfSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian)); 288f0b65372SJed Brown } 289e3663b90SJames Wright if (app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) PetscCall(SgsDDSetup(ceed, honee, problem)); 2906d0190e2SJames Wright } 29191933550SJames Wright 292d3c60affSJames Wright if (problem->use_strong_bc_ceed) PetscCall(SetupStrongBC_Ceed(ceed, dm, honee, problem)); 2939ae013d6SJames Wright if (app_ctx->sgs_train_enable) PetscCall(SGS_DD_TrainingSetup(ceed, honee)); 294e3663b90SJames Wright if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) PetscCall(DivDiffFluxProjectionSetup(honee, honee->diff_flux_proj)); 29591933550SJames Wright 296be29160dSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 2978c85b835SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&jac_data)); 298cf8f12c8SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 299be29160dSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 3008c85b835SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i)); 3018c85b835SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux)); 3028c85b835SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 303cf8f12c8SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 3040880fbb6SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&div_diff_flux_ceed)); 305b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian_vol)); 306bbc90103SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_ifunction_vol)); 307bbc90103SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_vol)); 308d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 309a515125bSLeila Ghaffari } 310