1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3 4 /// @file 5 /// Setup libCEED Operators for HONEE 6 7 #include <ceed.h> 8 #include <petscdmplex.h> 9 #include <smartsim.h> 10 11 #include <navierstokes.h> 12 13 // @brief Create CeedOperator for unstabilized mass KSP for explicit timestepping 14 static PetscErrorCode CreateKSPMassOperator_Unstabilized(Honee honee, CeedOperator *op_mass) { 15 Ceed ceed = honee->ceed; 16 CeedInt num_comp_q, q_data_size; 17 CeedQFunction qf_mass; 18 CeedElemRestriction elem_restr_q, elem_restr_qd; 19 CeedBasis basis_q; 20 CeedVector q_data; 21 22 PetscFunctionBeginUser; 23 { // Get restriction and basis from the RHS function 24 CeedOperator *sub_ops; 25 CeedOperatorField op_field; 26 PetscInt sub_op_index = 0; // will be 0 for the volume op 27 28 PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_rhs_ctx->op, &sub_ops)); 29 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &op_field)); 30 PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_q, &basis_q, NULL)); 31 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &op_field)); 32 PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_qd, NULL, &q_data)); 33 } 34 35 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q)); 36 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd, &q_data_size)); 37 38 PetscCall(HoneeMassQFunctionCreate(ceed, num_comp_q, q_data_size, &qf_mass)); 39 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass)); 40 PetscCallCeed(ceed, CeedOperatorSetName(*op_mass, "RHS Mass Operator")); 41 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "u", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 42 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 43 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 44 45 PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 46 PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 47 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 48 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 49 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 50 PetscFunctionReturn(PETSC_SUCCESS); 51 } 52 53 // @brief Create KSP to solve the inverse mass operator for explicit time stepping schemes 54 static PetscErrorCode CreateKSPMass(Honee honee, ProblemData problem) { 55 Ceed ceed = honee->ceed; 56 DM dm = honee->dm; 57 CeedOperator op_mass; 58 59 PetscFunctionBeginUser; 60 if (problem->create_mass_operator) PetscCall(problem->create_mass_operator(honee, &op_mass)); 61 else PetscCall(CreateKSPMassOperator_Unstabilized(honee, &op_mass)); 62 63 { // -- Setup KSP for mass operator 64 Mat mat_mass; 65 Vec Zeros_loc; 66 MPI_Comm comm = PetscObjectComm((PetscObject)dm); 67 68 PetscCall(DMCreateLocalVector(dm, &Zeros_loc)); 69 PetscCall(VecZeroEntries(Zeros_loc)); 70 PetscCall(MatCreateCeed(dm, dm, op_mass, NULL, &mat_mass)); 71 PetscCall(MatCeedSetLocalVectors(mat_mass, Zeros_loc, NULL)); 72 73 PetscCall(KSPCreate(comm, &honee->mass_ksp)); 74 PetscCall(KSPSetOptionsPrefix(honee->mass_ksp, "mass_")); 75 PetscCall(PetscObjectSetName((PetscObject)honee->mass_ksp, "Explicit Mass")); 76 { // lumped by default 77 PC pc; 78 PetscCall(KSPGetPC(honee->mass_ksp, &pc)); 79 PetscCall(PCSetType(pc, PCJACOBI)); 80 PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 81 PetscCall(KSPSetType(honee->mass_ksp, KSPPREONLY)); 82 } 83 PetscCall(KSPSetFromOptions_WithMatCeed(honee->mass_ksp, mat_mass)); 84 PetscCall(VecDestroy(&Zeros_loc)); 85 PetscCall(MatDestroy(&mat_mass)); 86 } 87 88 PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 89 PetscFunctionReturn(PETSC_SUCCESS); 90 } 91 92 PetscErrorCode SetupLibceed(Ceed ceed, DM dm, Honee honee, AppCtx app_ctx, ProblemData problem) { 93 const PetscInt num_comp_q = problem->num_components; 94 PetscInt dim; 95 CeedInt num_comps_jac_data = problem->num_comps_jac_data, num_comp_x, q_data_size_vol; 96 CeedElemRestriction elem_restr_jd_i = NULL, elem_restr_qd; 97 CeedVector jac_data = NULL, q_data; 98 CeedOperator op_ifunction_vol = NULL, op_rhs_vol = NULL, op_ijacobian_vol = NULL; 99 100 PetscFunctionBeginUser; 101 PetscCall(DMGetDimension(dm, &dim)); 102 num_comp_x = dim; 103 104 CeedElemRestriction elem_restr_diff_flux = NULL; 105 CeedVector div_diff_flux_ceed = NULL; 106 CeedBasis basis_diff_flux = NULL; 107 CeedEvalMode eval_mode_diff_flux = -1; 108 { // Create bases and element restrictions 109 DMLabel domain_label = NULL; 110 PetscInt label_value = 0, height = 0, dm_field = 0; 111 DM dm_coord; 112 113 PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 114 PetscCall(CreateBasisFromPlex(ceed, dm, domain_label, label_value, height, dm_field, &honee->basis_q)); 115 PetscCall(CreateBasisFromPlex(ceed, dm_coord, domain_label, label_value, height, dm_field, &honee->basis_x)); 116 117 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, 0, &honee->elem_restr_q)); 118 PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &honee->elem_restr_x)); 119 if (num_comps_jac_data) { 120 PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, num_comps_jac_data, &elem_restr_jd_i)); 121 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL)); 122 } 123 124 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(honee->elem_restr_q, &honee->q_ceed, NULL)); 125 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(honee->elem_restr_q, &honee->q_dot_ceed, NULL)); 126 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(honee->elem_restr_q, &honee->g_ceed, NULL)); 127 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(honee->elem_restr_x, &honee->x_coord, NULL)); 128 129 { // -- Copy PETSc coordinate vector into CEED vector 130 Vec X_loc; 131 DM cdm; 132 133 PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 134 if (cdm) { 135 PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc)); 136 } else { 137 PetscCall(DMGetCoordinatesLocal(dm, &X_loc)); 138 } 139 PetscCall(VecScale(X_loc, honee->units->meter)); 140 PetscCall(VecCopyPetscToCeed(X_loc, honee->x_coord)); 141 } 142 143 PetscCall(QDataGet(ceed, dm, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, &q_data, 144 &q_data_size_vol)); 145 } 146 147 if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) { 148 PetscCheck(honee->diff_flux_proj, honee->comm, PETSC_ERR_ARG_WRONGSTATE, 149 "Divergence of diffusive flux projection requested but object not created"); 150 PetscCall(DivDiffFluxProjectionGetOperatorFieldData(honee->diff_flux_proj, &elem_restr_diff_flux, &basis_diff_flux, &div_diff_flux_ceed, 151 &eval_mode_diff_flux)); 152 } 153 154 { // -- Create QFunction for ICs 155 CeedBasis basis_xc; 156 CeedQFunction qf_ics; 157 CeedOperator op_ics; 158 159 PetscCallCeed(ceed, CeedBasisCreateProjection(honee->basis_x, honee->basis_q, &basis_xc)); 160 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->ics.qf_func_ptr, problem->ics.qf_loc, &qf_ics)); 161 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ics, problem->ics.qfctx)); 162 PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ics, 0)); 163 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "x", num_comp_x, CEED_EVAL_INTERP)); 164 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "dx", num_comp_x * dim, CEED_EVAL_GRAD)); 165 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ics, "q0", num_comp_q, CEED_EVAL_NONE)); 166 167 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics)); 168 PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "x", honee->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE)); 169 PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "dx", honee->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE)); 170 PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "q0", honee->elem_restr_q, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 171 PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ics, "evaluation time", &honee->phys->ics_time_label)); 172 PetscCall(OperatorApplyContextCreate(NULL, dm, honee->ceed, op_ics, honee->x_coord, NULL, NULL, honee->Q_loc, &honee->op_ics_ctx)); 173 174 PetscCallCeed(ceed, CeedBasisDestroy(&basis_xc)); 175 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ics)); 176 PetscCallCeed(ceed, CeedOperatorDestroy(&op_ics)); 177 } 178 179 if (problem->apply_vol_rhs.qf_func_ptr) { 180 CeedQFunction qf_rhs_vol; 181 182 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qf_func_ptr, problem->apply_vol_rhs.qf_loc, &qf_rhs_vol)); 183 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_vol, problem->apply_vol_rhs.qfctx)); 184 PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_rhs_vol, 0)); 185 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP)); 186 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 187 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE)); 188 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP)); 189 if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) 190 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "div F_diff", honee->diff_flux_proj->num_diff_flux_comps, eval_mode_diff_flux)); 191 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP)); 192 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD)); 193 194 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_vol, NULL, NULL, &op_rhs_vol)); 195 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 196 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 197 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 198 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "x", honee->elem_restr_x, honee->basis_x, honee->x_coord)); 199 if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) 200 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "div F_diff", elem_restr_diff_flux, basis_diff_flux, div_diff_flux_ceed)); 201 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 202 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "Grad_v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 203 204 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_vol)); 205 } 206 207 if (problem->apply_vol_ifunction.qf_func_ptr) { 208 CeedQFunction qf_ifunction_vol; 209 210 PetscCallCeed( 211 ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qf_func_ptr, problem->apply_vol_ifunction.qf_loc, &qf_ifunction_vol)); 212 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ifunction_vol, problem->apply_vol_ifunction.qfctx)); 213 PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ifunction_vol, 0)); 214 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP)); 215 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 216 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP)); 217 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE)); 218 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP)); 219 if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) 220 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "div F_diff", honee->diff_flux_proj->num_diff_flux_comps, eval_mode_diff_flux)); 221 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP)); 222 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD)); 223 if (num_comps_jac_data) PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "jac_data", num_comps_jac_data, CEED_EVAL_NONE)); 224 225 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ifunction_vol, NULL, NULL, &op_ifunction_vol)); 226 PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 227 PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 228 PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "q dot", honee->elem_restr_q, honee->basis_q, honee->q_dot_ceed)); 229 PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 230 PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "x", honee->elem_restr_x, honee->basis_x, honee->x_coord)); 231 if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) 232 PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "div F_diff", elem_restr_diff_flux, basis_diff_flux, div_diff_flux_ceed)); 233 PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 234 PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "Grad_v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 235 if (num_comps_jac_data) PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data)); 236 237 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ifunction_vol)); 238 } 239 240 if (problem->apply_vol_ijacobian.qf_func_ptr) { 241 CeedQFunction qf_ijacobian_vol; 242 243 PetscCallCeed( 244 ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qf_func_ptr, problem->apply_vol_ijacobian.qf_loc, &qf_ijacobian_vol)); 245 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfctx)); 246 PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ijacobian_vol, 0)); 247 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP)); 248 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD)); 249 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE)); 250 if (num_comps_jac_data) PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", num_comps_jac_data, CEED_EVAL_NONE)); 251 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP)); 252 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD)); 253 254 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op_ijacobian_vol)); 255 PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "dq", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 256 PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "Grad_dq", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 257 PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 258 if (num_comps_jac_data) PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data)); 259 PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 260 PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "Grad_v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 261 262 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ijacobian_vol)); 263 } 264 265 // -- Create and apply CEED Composite Operator for the entire domain 266 if (!honee->phys->implicit) { // RHS 267 CeedOperator op_rhs; 268 269 PetscCallCeed(ceed, CeedOperatorCreateComposite(ceed, &op_rhs)); 270 PetscCallCeed(ceed, CeedOperatorCompositeAddSub(op_rhs, op_rhs_vol)); 271 for (PetscCount b = 0; b < problem->num_bc_defs; b++) PetscCall(BCDefinitionAddOperators(problem->bc_defs[b], op_rhs, NULL)); 272 273 PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_rhs, honee->q_ceed, honee->g_ceed, honee->Q_loc, NULL, &honee->op_rhs_ctx)); 274 275 // ----- Get Context Labels for Operator 276 PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "solution time", &honee->phys->solution_time_label)); 277 PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "timestep size", &honee->phys->timestep_size_label)); 278 279 PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs)); 280 PetscCall(CreateKSPMass(honee, problem)); 281 PetscCheck(app_ctx->sgs_model_type == SGS_MODEL_NONE, honee->comm, PETSC_ERR_SUP, "SGS modeling not implemented for explicit timestepping"); 282 } else { // IFunction 283 CeedOperator op_ijacobian = NULL; 284 285 // Create Composite Operaters 286 PetscCallCeed(ceed, CeedOperatorCreateComposite(ceed, &honee->op_ifunction)); 287 PetscCallCeed(ceed, CeedOperatorCompositeAddSub(honee->op_ifunction, op_ifunction_vol)); 288 if (op_ijacobian_vol) { 289 PetscCallCeed(ceed, CeedOperatorCreateComposite(ceed, &op_ijacobian)); 290 PetscCallCeed(ceed, CeedOperatorCompositeAddSub(op_ijacobian, op_ijacobian_vol)); 291 } 292 for (PetscCount b = 0; b < problem->num_bc_defs; b++) PetscCall(BCDefinitionAddOperators(problem->bc_defs[b], honee->op_ifunction, op_ijacobian)); 293 294 // ----- Get Context Labels for Operator 295 PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(honee->op_ifunction, "solution time", &honee->phys->solution_time_label)); 296 PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(honee->op_ifunction, "timestep size", &honee->phys->timestep_size_label)); 297 298 if (op_ijacobian) { 299 PetscCall(MatCreateCeed(honee->dm, honee->dm, op_ijacobian, NULL, &honee->mat_ijacobian)); 300 PetscCall(MatCeedSetLocalVectors(honee->mat_ijacobian, honee->Q_dot_loc, NULL)); 301 PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian)); 302 } 303 if (app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) PetscCall(SgsDDSetup(ceed, honee, problem)); 304 } 305 306 if (problem->use_strong_bc_ceed) PetscCall(SetupStrongBC_Ceed(ceed, dm, honee, problem)); 307 if (app_ctx->sgs_train_enable) PetscCall(SGS_DD_TrainingSetup(ceed, honee)); 308 if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) PetscCall(DivDiffFluxProjectionSetup(honee, honee->diff_flux_proj)); 309 310 PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 311 PetscCallCeed(ceed, CeedVectorDestroy(&jac_data)); 312 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 313 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i)); 314 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux)); 315 PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 316 PetscCallCeed(ceed, CeedVectorDestroy(&div_diff_flux_ceed)); 317 PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian_vol)); 318 PetscCallCeed(ceed, CeedOperatorDestroy(&op_ifunction_vol)); 319 PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_vol)); 320 PetscFunctionReturn(PETSC_SUCCESS); 321 } 322