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