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 10 #include <navierstokes.h> 11 12 // @brief Create CeedOperator for unstabilized mass KSP for explicit timestepping 13 static PetscErrorCode CreateKSPMassOperator_Unstabilized(Honee honee, CeedOperator *op_mass) { 14 Ceed ceed = honee->ceed; 15 CeedInt num_comp_q, q_data_size; 16 CeedQFunction qf_mass; 17 CeedElemRestriction elem_restr_q, elem_restr_qd; 18 CeedBasis basis_q; 19 CeedVector q_data; 20 21 PetscFunctionBeginUser; 22 { // Get restriction and basis from the RHS function 23 CeedOperator *sub_ops; 24 CeedOperatorField op_field; 25 PetscInt sub_op_index = 0; // will be 0 for the volume op 26 27 PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_rhs_ctx->op, &sub_ops)); 28 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &op_field)); 29 PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_q, &basis_q, NULL)); 30 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &op_field)); 31 PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_qd, NULL, &q_data)); 32 } 33 34 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q)); 35 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd, &q_data_size)); 36 37 PetscCall(HoneeMassQFunctionCreate(ceed, num_comp_q, q_data_size, &qf_mass)); 38 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass)); 39 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "u", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 40 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 41 PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 42 43 PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 44 PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 45 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 46 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 47 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 48 PetscFunctionReturn(PETSC_SUCCESS); 49 } 50 51 // @brief Create KSP to solve the inverse mass operator for explicit time stepping schemes 52 static PetscErrorCode CreateKSPMass(Honee honee, ProblemData problem) { 53 Ceed ceed = honee->ceed; 54 DM dm = honee->dm; 55 CeedOperator op_mass; 56 57 PetscFunctionBeginUser; 58 if (problem->create_mass_operator) PetscCall(problem->create_mass_operator(honee, &op_mass)); 59 else PetscCall(CreateKSPMassOperator_Unstabilized(honee, &op_mass)); 60 61 { // -- Setup KSP for mass operator 62 Mat mat_mass; 63 Vec Zeros_loc; 64 MPI_Comm comm = PetscObjectComm((PetscObject)dm); 65 66 PetscCall(DMCreateLocalVector(dm, &Zeros_loc)); 67 PetscCall(VecZeroEntries(Zeros_loc)); 68 PetscCall(MatCreateCeed(dm, dm, op_mass, NULL, &mat_mass)); 69 PetscCall(MatCeedSetLocalVectors(mat_mass, Zeros_loc, NULL)); 70 71 PetscCall(KSPCreate(comm, &honee->mass_ksp)); 72 PetscCall(KSPSetOptionsPrefix(honee->mass_ksp, "mass_")); 73 { // lumped by default 74 PC pc; 75 PetscCall(KSPGetPC(honee->mass_ksp, &pc)); 76 PetscCall(PCSetType(pc, PCJACOBI)); 77 PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 78 PetscCall(KSPSetType(honee->mass_ksp, KSPPREONLY)); 79 } 80 PetscCall(KSPSetFromOptions_WithMatCeed(honee->mass_ksp, mat_mass)); 81 PetscCall(VecDestroy(&Zeros_loc)); 82 PetscCall(MatDestroy(&mat_mass)); 83 } 84 85 PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 86 PetscFunctionReturn(PETSC_SUCCESS); 87 } 88 89 static PetscErrorCode AddBCSubOperator(Ceed ceed, DM dm, Honee honee, DMLabel domain_label, PetscInt label_value, CeedInt height, CeedInt Q_sur, 90 CeedInt q_data_size_sur, CeedInt jac_data_size_sur, CeedBasis basis_q_sur, CeedBasis basis_x_sur, 91 CeedQFunction qf_apply_bc, CeedQFunction qf_apply_bc_jacobian, CeedOperator op_apply, 92 CeedOperator op_apply_ijacobian) { 93 CeedVector q_data_sur, jac_data_sur = NULL; 94 CeedOperator op_apply_bc, op_apply_bc_jacobian = NULL; 95 CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur, elem_restr_qd_i_sur, elem_restr_jd_i_sur = NULL; 96 PetscInt dm_field = 0; 97 98 PetscFunctionBeginUser; 99 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, dm_field, &elem_restr_q_sur)); 100 PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &elem_restr_x_sur)); 101 if (jac_data_size_sur > 0) { 102 // State-dependent data will be passed from residual to Jacobian. This will be collocated. 103 PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_sur, &elem_restr_jd_i_sur)); 104 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i_sur, &jac_data_sur, NULL)); 105 } 106 107 PetscCall(QDataBoundaryGet(ceed, dm, domain_label, label_value, elem_restr_x_sur, basis_x_sur, honee->x_coord, &elem_restr_qd_i_sur, &q_data_sur, 108 &q_data_size_sur)); 109 110 // CEED Operator for Physics 111 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc, NULL, NULL, &op_apply_bc)); 112 PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "q", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE)); 113 PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "Grad_q", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE)); 114 PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur)); 115 PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "x", elem_restr_x_sur, basis_x_sur, honee->x_coord)); 116 PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "v", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE)); 117 if (elem_restr_jd_i_sur) 118 PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur)); 119 120 if (qf_apply_bc_jacobian && elem_restr_jd_i_sur) { 121 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc_jacobian, NULL, NULL, &op_apply_bc_jacobian)); 122 PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "dq", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE)); 123 PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "Grad_dq", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE)); 124 PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur)); 125 PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "x", elem_restr_x_sur, basis_x_sur, honee->x_coord)); 126 PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur)); 127 PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "v", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE)); 128 } 129 130 // Apply Sub-Operator for Physics 131 PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_apply, op_apply_bc)); 132 if (op_apply_bc_jacobian) PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_apply_ijacobian, op_apply_bc_jacobian)); 133 134 PetscCallCeed(ceed, CeedVectorDestroy(&q_data_sur)); 135 PetscCallCeed(ceed, CeedVectorDestroy(&jac_data_sur)); 136 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q_sur)); 137 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_sur)); 138 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_i_sur)); 139 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i_sur)); 140 PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc)); 141 PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc_jacobian)); 142 PetscFunctionReturn(PETSC_SUCCESS); 143 } 144 145 static PetscErrorCode SetupBCQFunctions(Ceed ceed, PetscInt dim_sur, PetscInt num_comp_x, PetscInt num_comp_q, PetscInt q_data_size_sur, 146 PetscInt jac_data_size_sur, ProblemQFunctionSpec apply_bc, ProblemQFunctionSpec apply_bc_jacobian, 147 CeedQFunction *qf_apply_bc, CeedQFunction *qf_apply_bc_jacobian) { 148 PetscFunctionBeginUser; 149 if (apply_bc.qf_func_ptr) { 150 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc.qf_func_ptr, apply_bc.qf_loc, qf_apply_bc)); 151 PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc, apply_bc.qfctx)); 152 PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc, 0)); 153 PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "q", num_comp_q, CEED_EVAL_INTERP)); 154 PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "Grad_q", num_comp_q * dim_sur, CEED_EVAL_GRAD)); 155 PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE)); 156 PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "x", num_comp_x, CEED_EVAL_INTERP)); 157 PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "v", num_comp_q, CEED_EVAL_INTERP)); 158 if (jac_data_size_sur) PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE)); 159 } 160 if (apply_bc_jacobian.qf_func_ptr) { 161 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc_jacobian.qf_func_ptr, apply_bc_jacobian.qf_loc, qf_apply_bc_jacobian)); 162 PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc_jacobian, apply_bc_jacobian.qfctx)); 163 PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc_jacobian, 0)); 164 PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "dq", num_comp_q, CEED_EVAL_INTERP)); 165 PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "Grad_dq", num_comp_q * dim_sur, CEED_EVAL_GRAD)); 166 PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface qdata", q_data_size_sur, CEED_EVAL_NONE)); 167 PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "x", num_comp_x, CEED_EVAL_INTERP)); 168 PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE)); 169 PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc_jacobian, "v", num_comp_q, CEED_EVAL_INTERP)); 170 } 171 PetscFunctionReturn(PETSC_SUCCESS); 172 } 173 174 // Utility function to add boundary operators to the composite operator 175 static PetscErrorCode AddBCSubOperators(Honee honee, Ceed ceed, DM dm, SimpleBC bc, ProblemData problem, CeedOperator op_apply, 176 CeedOperator op_apply_ijacobian) { 177 CeedInt height = 1, num_comp_q, num_comp_x; 178 CeedInt P_sur = honee->app_ctx->degree + 1, Q_sur = P_sur + honee->app_ctx->q_extra, dim_sur, q_data_size_sur; 179 const CeedInt jac_data_size_sur = honee->phys->implicit ? problem->jac_data_size_sur : 0; 180 PetscInt dim; 181 DMLabel face_sets_label; 182 CeedBasis basis_q_sur, basis_x_sur; 183 184 PetscFunctionBeginUser; 185 PetscCall(DMGetDimension(dm, &dim)); 186 PetscCall(QDataBoundaryGetNumComponents(dm, &q_data_size_sur)); 187 dim_sur = dim - height; 188 { // Get number of components and coordinate dimension from op_apply 189 CeedOperator *sub_ops; 190 CeedOperatorField field; 191 PetscInt sub_op_index = 0; // will be 0 for the volume op 192 CeedElemRestriction elem_restr_q, elem_restr_x; 193 194 PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(op_apply, &sub_ops)); 195 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field)); 196 PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q)); 197 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q)); 198 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 199 200 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "x", &field)); 201 PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_x)); 202 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x)); 203 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x)); 204 } 205 206 { // Get bases 207 DM dm_coord; 208 209 PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 210 DMLabel label = NULL; 211 PetscInt label_value = 0; 212 PetscInt field = 0; 213 PetscCall(CreateBasisFromPlex(ceed, dm, label, label_value, height, field, &basis_q_sur)); 214 PetscCall(CreateBasisFromPlex(ceed, dm_coord, label, label_value, height, field, &basis_x_sur)); 215 } 216 217 PetscCall(DMGetLabel(dm, "Face Sets", &face_sets_label)); 218 219 { // --- Create Sub-Operator for inflow boundaries 220 CeedQFunction qf_apply_inflow = NULL, qf_apply_inflow_jacobian = NULL; 221 222 PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_inflow, 223 problem->apply_inflow_jacobian, &qf_apply_inflow, &qf_apply_inflow_jacobian)); 224 for (CeedInt i = 0; i < bc->num_inflow; i++) { 225 PetscCall(AddBCSubOperator(ceed, dm, honee, face_sets_label, bc->inflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, basis_q_sur, 226 basis_x_sur, qf_apply_inflow, qf_apply_inflow_jacobian, op_apply, op_apply_ijacobian)); 227 } 228 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_inflow)); 229 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_inflow_jacobian)); 230 } 231 232 { // --- Create Sub-Operator for outflow boundaries 233 CeedQFunction qf_apply_outflow = NULL, qf_apply_outflow_jacobian = NULL; 234 235 PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_outflow, 236 problem->apply_outflow_jacobian, &qf_apply_outflow, &qf_apply_outflow_jacobian)); 237 for (CeedInt i = 0; i < bc->num_outflow; i++) { 238 PetscCall(AddBCSubOperator(ceed, dm, honee, face_sets_label, bc->outflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, basis_q_sur, 239 basis_x_sur, qf_apply_outflow, qf_apply_outflow_jacobian, op_apply, op_apply_ijacobian)); 240 } 241 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_outflow)); 242 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_outflow_jacobian)); 243 } 244 245 { // --- Create Sub-Operator for freestream boundaries 246 CeedQFunction qf_apply_freestream = NULL, qf_apply_freestream_jacobian = NULL; 247 248 PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_freestream, 249 problem->apply_freestream_jacobian, &qf_apply_freestream, &qf_apply_freestream_jacobian)); 250 for (CeedInt i = 0; i < bc->num_freestream; i++) { 251 PetscCall(AddBCSubOperator(ceed, dm, honee, face_sets_label, bc->freestreams[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, basis_q_sur, 252 basis_x_sur, qf_apply_freestream, qf_apply_freestream_jacobian, op_apply, op_apply_ijacobian)); 253 } 254 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_freestream)); 255 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_freestream_jacobian)); 256 } 257 258 { // --- Create Sub-Operator for slip boundaries 259 CeedQFunction qf_apply_slip = NULL, qf_apply_slip_jacobian = NULL; 260 261 PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_slip, 262 problem->apply_slip_jacobian, &qf_apply_slip, &qf_apply_slip_jacobian)); 263 for (CeedInt i = 0; i < bc->num_slip; i++) { 264 PetscCall(AddBCSubOperator(ceed, dm, honee, face_sets_label, bc->slips[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, basis_q_sur, 265 basis_x_sur, qf_apply_slip, qf_apply_slip_jacobian, op_apply, op_apply_ijacobian)); 266 } 267 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_slip)); 268 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_slip_jacobian)); 269 } 270 271 PetscCallCeed(ceed, CeedBasisDestroy(&basis_q_sur)); 272 PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_sur)); 273 PetscFunctionReturn(PETSC_SUCCESS); 274 } 275 276 PetscErrorCode SetupLibceed(Ceed ceed, DM dm, Honee honee, AppCtx app_ctx, ProblemData problem, SimpleBC bc) { 277 const PetscInt num_comp_q = 5; 278 PetscInt dim; 279 CeedInt jac_data_size_vol = problem->jac_data_size_vol, num_comp_x, q_data_size_vol; 280 CeedElemRestriction elem_restr_jd_i = NULL, elem_restr_qd; 281 CeedVector jac_data = NULL, q_data; 282 CeedOperator op_ifunction_vol = NULL, op_rhs_vol = NULL, op_ijacobian_vol = NULL; 283 284 PetscFunctionBeginUser; 285 PetscCall(DMGetDimension(dm, &dim)); 286 num_comp_x = dim; 287 288 CeedElemRestriction elem_restr_diff_flux = NULL; 289 CeedVector div_diff_flux_ceed = NULL; 290 CeedBasis basis_diff_flux = NULL; 291 CeedEvalMode eval_mode_diff_flux = -1; 292 { // Create bases and element restrictions 293 DMLabel domain_label = NULL; 294 PetscInt label_value = 0, height = 0, dm_field = 0; 295 DM dm_coord; 296 297 PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 298 PetscCall(CreateBasisFromPlex(ceed, dm, domain_label, label_value, height, dm_field, &honee->basis_q)); 299 PetscCall(CreateBasisFromPlex(ceed, dm_coord, domain_label, label_value, height, dm_field, &honee->basis_x)); 300 301 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, 0, &honee->elem_restr_q)); 302 PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &honee->elem_restr_x)); 303 if (jac_data_size_vol) { 304 PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_vol, &elem_restr_jd_i)); 305 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL)); 306 } 307 308 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(honee->elem_restr_q, &honee->q_ceed, NULL)); 309 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(honee->elem_restr_q, &honee->q_dot_ceed, NULL)); 310 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(honee->elem_restr_q, &honee->g_ceed, NULL)); 311 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(honee->elem_restr_x, &honee->x_coord, NULL)); 312 313 { // -- Copy PETSc coordinate vector into CEED vector 314 Vec X_loc; 315 DM cdm; 316 317 PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 318 if (cdm) { 319 PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc)); 320 } else { 321 PetscCall(DMGetCoordinatesLocal(dm, &X_loc)); 322 } 323 PetscCall(VecScale(X_loc, honee->units->meter)); 324 PetscCall(VecCopyPetscToCeed(X_loc, honee->x_coord)); 325 } 326 327 PetscCall(QDataGet(ceed, dm, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, &q_data, 328 &q_data_size_vol)); 329 } 330 331 if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) { 332 PetscCheck(honee->diff_flux_proj, honee->comm, PETSC_ERR_ARG_WRONGSTATE, 333 "Divergence of diffusive flux projection requested but object not created"); 334 PetscCall(DivDiffFluxProjectionGetOperatorFieldData(honee->diff_flux_proj, &elem_restr_diff_flux, &basis_diff_flux, &div_diff_flux_ceed, 335 &eval_mode_diff_flux)); 336 } 337 338 { // -- Create QFunction for ICs 339 CeedBasis basis_xc; 340 CeedQFunction qf_ics; 341 CeedOperator op_ics; 342 343 PetscCallCeed(ceed, CeedBasisCreateProjection(honee->basis_x, honee->basis_q, &basis_xc)); 344 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->ics.qf_func_ptr, problem->ics.qf_loc, &qf_ics)); 345 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ics, problem->ics.qfctx)); 346 PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ics, 0)); 347 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "x", num_comp_x, CEED_EVAL_INTERP)); 348 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "dx", num_comp_x * dim, CEED_EVAL_GRAD)); 349 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ics, "q0", num_comp_q, CEED_EVAL_NONE)); 350 351 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics)); 352 PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "x", honee->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE)); 353 PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "dx", honee->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE)); 354 PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "q0", honee->elem_restr_q, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 355 PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ics, "evaluation time", &honee->phys->ics_time_label)); 356 PetscCall(OperatorApplyContextCreate(NULL, dm, honee->ceed, op_ics, honee->x_coord, NULL, NULL, honee->Q_loc, &honee->op_ics_ctx)); 357 358 PetscCallCeed(ceed, CeedBasisDestroy(&basis_xc)); 359 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ics)); 360 PetscCallCeed(ceed, CeedOperatorDestroy(&op_ics)); 361 } 362 363 if (problem->apply_vol_rhs.qf_func_ptr) { 364 CeedQFunction qf_rhs_vol; 365 366 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qf_func_ptr, problem->apply_vol_rhs.qf_loc, &qf_rhs_vol)); 367 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_vol, problem->apply_vol_rhs.qfctx)); 368 PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_rhs_vol, 0)); 369 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP)); 370 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 371 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE)); 372 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP)); 373 if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) 374 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "div F_diff", honee->diff_flux_proj->num_diff_flux_comps, eval_mode_diff_flux)); 375 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP)); 376 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD)); 377 378 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_vol, NULL, NULL, &op_rhs_vol)); 379 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 380 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 381 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 382 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "x", honee->elem_restr_x, honee->basis_x, honee->x_coord)); 383 if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) 384 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "div F_diff", elem_restr_diff_flux, basis_diff_flux, div_diff_flux_ceed)); 385 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 386 PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "Grad_v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 387 388 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_vol)); 389 } 390 391 if (problem->apply_vol_ifunction.qf_func_ptr) { 392 CeedQFunction qf_ifunction_vol; 393 394 PetscCallCeed( 395 ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qf_func_ptr, problem->apply_vol_ifunction.qf_loc, &qf_ifunction_vol)); 396 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ifunction_vol, problem->apply_vol_ifunction.qfctx)); 397 PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ifunction_vol, 0)); 398 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP)); 399 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 400 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP)); 401 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE)); 402 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP)); 403 if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) 404 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "div F_diff", honee->diff_flux_proj->num_diff_flux_comps, eval_mode_diff_flux)); 405 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP)); 406 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD)); 407 if (jac_data_size_vol) PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE)); 408 409 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ifunction_vol, NULL, NULL, &op_ifunction_vol)); 410 PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 411 PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 412 PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "q dot", honee->elem_restr_q, honee->basis_q, honee->q_dot_ceed)); 413 PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 414 PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "x", honee->elem_restr_x, honee->basis_x, honee->x_coord)); 415 if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) 416 PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "div F_diff", elem_restr_diff_flux, basis_diff_flux, div_diff_flux_ceed)); 417 PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 418 PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "Grad_v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 419 if (jac_data_size_vol) PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data)); 420 421 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ifunction_vol)); 422 } 423 424 if (problem->apply_vol_ijacobian.qf_func_ptr) { 425 CeedQFunction qf_ijacobian_vol; 426 427 PetscCallCeed( 428 ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qf_func_ptr, problem->apply_vol_ijacobian.qf_loc, &qf_ijacobian_vol)); 429 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfctx)); 430 PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ijacobian_vol, 0)); 431 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP)); 432 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD)); 433 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE)); 434 if (jac_data_size_vol) PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE)); 435 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP)); 436 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD)); 437 438 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op_ijacobian_vol)); 439 PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "dq", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 440 PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "Grad_dq", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 441 PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 442 if (jac_data_size_vol) PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data)); 443 PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 444 PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "Grad_v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 445 446 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ijacobian_vol)); 447 } 448 449 // -- Create and apply CEED Composite Operator for the entire domain 450 if (!honee->phys->implicit) { // RHS 451 CeedOperator op_rhs; 452 453 PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_rhs)); 454 PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_rhs, op_rhs_vol)); 455 PetscCall(AddBCSubOperators(honee, ceed, dm, bc, problem, op_rhs, NULL)); 456 457 PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_rhs, honee->q_ceed, honee->g_ceed, honee->Q_loc, NULL, &honee->op_rhs_ctx)); 458 459 // ----- Get Context Labels for Operator 460 PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "solution time", &honee->phys->solution_time_label)); 461 PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "timestep size", &honee->phys->timestep_size_label)); 462 463 PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs)); 464 PetscCall(CreateKSPMass(honee, problem)); 465 PetscCheck(app_ctx->sgs_model_type == SGS_MODEL_NONE, honee->comm, PETSC_ERR_SUP, "SGS modeling not implemented for explicit timestepping"); 466 } else { // IFunction 467 CeedOperator op_ijacobian = NULL; 468 469 // Create Composite Operaters 470 PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &honee->op_ifunction)); 471 PetscCallCeed(ceed, CeedCompositeOperatorAddSub(honee->op_ifunction, op_ifunction_vol)); 472 if (op_ijacobian_vol) { 473 PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_ijacobian)); 474 PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_ijacobian, op_ijacobian_vol)); 475 } 476 PetscCall(AddBCSubOperators(honee, ceed, dm, bc, problem, honee->op_ifunction, op_ijacobian)); 477 478 // ----- Get Context Labels for Operator 479 PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(honee->op_ifunction, "solution time", &honee->phys->solution_time_label)); 480 PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(honee->op_ifunction, "timestep size", &honee->phys->timestep_size_label)); 481 482 if (op_ijacobian) { 483 PetscCall(MatCreateCeed(honee->dm, honee->dm, op_ijacobian, NULL, &honee->mat_ijacobian)); 484 PetscCall(MatCeedSetLocalVectors(honee->mat_ijacobian, honee->Q_dot_loc, NULL)); 485 PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian)); 486 } 487 if (app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) PetscCall(SgsDDSetup(ceed, honee, problem)); 488 } 489 490 if (problem->use_strong_bc_ceed) PetscCall(SetupStrongBC_Ceed(ceed, dm, honee, problem, bc)); 491 if (app_ctx->turb_spanstats_enable) PetscCall(TurbulenceStatisticsSetup(ceed, honee, problem)); 492 if (app_ctx->diff_filter_monitor && !honee->diff_filter) PetscCall(DifferentialFilterSetup(ceed, honee, problem)); 493 if (app_ctx->sgs_train_enable) PetscCall(SGS_DD_TrainingSetup(ceed, honee, problem)); 494 if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) PetscCall(DivDiffFluxProjectionSetup(honee, honee->diff_flux_proj)); 495 496 PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 497 PetscCallCeed(ceed, CeedVectorDestroy(&jac_data)); 498 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 499 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i)); 500 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux)); 501 PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 502 PetscCallCeed(ceed, CeedVectorDestroy(&div_diff_flux_ceed)); 503 PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian_vol)); 504 PetscCallCeed(ceed, CeedOperatorDestroy(&op_ifunction_vol)); 505 PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_vol)); 506 PetscFunctionReturn(PETSC_SUCCESS); 507 } 508