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