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