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