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