1 // Copyright (c) 2017-2022, 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 // Utility function - essential BC dofs are encoded in closure indices as -(i+1). 17 PetscInt Involute(PetscInt i) { return i >= 0 ? i : -(i + 1); } 18 19 // Utility function to create local CEED restriction 20 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt label_value, PetscInt dm_field, 21 CeedElemRestriction *elem_restr) { 22 PetscInt num_elem, elem_size, num_dof, num_comp, *elem_restr_offsets_petsc; 23 CeedInt *elem_restr_offsets_ceed; 24 25 PetscFunctionBeginUser; 26 PetscCall( 27 DMPlexGetLocalOffsets(dm, domain_label, label_value, height, dm_field, &num_elem, &elem_size, &num_comp, &num_dof, &elem_restr_offsets_petsc)); 28 29 PetscCall(IntArrayP2C(num_elem * elem_size, &elem_restr_offsets_petsc, &elem_restr_offsets_ceed)); 30 CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp, 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES, elem_restr_offsets_ceed, elem_restr); 31 PetscCall(PetscFree(elem_restr_offsets_ceed)); 32 33 PetscFunctionReturn(PETSC_SUCCESS); 34 } 35 36 // Utility function to get Ceed Restriction for each domain 37 PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, PetscInt label_value, PetscInt dm_field, CeedInt Q, 38 CeedInt q_data_size, CeedElemRestriction *elem_restr_q, CeedElemRestriction *elem_restr_x, 39 CeedElemRestriction *elem_restr_qd_i) { 40 DM dm_coord; 41 CeedInt Q_dim, loc_num_elem; 42 PetscInt dim; 43 CeedElemRestriction elem_restr_tmp; 44 PetscFunctionBeginUser; 45 46 PetscCall(DMGetDimension(dm, &dim)); 47 dim -= height; 48 Q_dim = CeedIntPow(Q, dim); 49 PetscCall(CreateRestrictionFromPlex(ceed, dm, height, domain_label, label_value, dm_field, &elem_restr_tmp)); 50 if (elem_restr_q) *elem_restr_q = elem_restr_tmp; 51 if (elem_restr_x) { 52 PetscCall(DMGetCellCoordinateDM(dm, &dm_coord)); 53 if (!dm_coord) { 54 PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 55 } 56 PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL)); 57 PetscCall(CreateRestrictionFromPlex(ceed, dm_coord, height, domain_label, label_value, dm_field, elem_restr_x)); 58 } 59 if (elem_restr_qd_i) { 60 CeedElemRestrictionGetNumElements(elem_restr_tmp, &loc_num_elem); 61 CeedElemRestrictionCreateStrided(ceed, loc_num_elem, Q_dim, q_data_size, q_data_size * loc_num_elem * Q_dim, CEED_STRIDES_BACKEND, 62 elem_restr_qd_i); 63 } 64 if (!elem_restr_q) CeedElemRestrictionDestroy(&elem_restr_tmp); 65 PetscFunctionReturn(PETSC_SUCCESS); 66 } 67 68 PetscErrorCode AddBCSubOperator(Ceed ceed, DM dm, CeedData ceed_data, DMLabel domain_label, PetscInt label_value, CeedInt height, CeedInt Q_sur, 69 CeedInt q_data_size_sur, CeedInt jac_data_size_sur, CeedQFunction qf_apply_bc, CeedQFunction qf_apply_bc_jacobian, 70 CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) { 71 CeedVector q_data_sur, jac_data_sur; 72 CeedOperator op_setup_sur, op_apply_bc, op_apply_bc_jacobian = NULL; 73 CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur, elem_restr_qd_i_sur, elem_restr_jd_i_sur; 74 CeedInt num_qpts_sur; 75 PetscFunctionBeginUser; 76 77 // --- Get number of quadrature points for the boundaries 78 CeedBasisGetNumQuadraturePoints(ceed_data->basis_q_sur, &num_qpts_sur); 79 80 // ---- CEED Restriction 81 PetscCall(GetRestrictionForDomain(ceed, dm, height, domain_label, label_value, 0, Q_sur, q_data_size_sur, &elem_restr_q_sur, &elem_restr_x_sur, 82 &elem_restr_qd_i_sur)); 83 if (jac_data_size_sur > 0) { 84 // State-dependent data will be passed from residual to Jacobian. This will be collocated. 85 PetscCall(GetRestrictionForDomain(ceed, dm, height, domain_label, label_value, 0, Q_sur, jac_data_size_sur, NULL, NULL, &elem_restr_jd_i_sur)); 86 CeedElemRestrictionCreateVector(elem_restr_jd_i_sur, &jac_data_sur, NULL); 87 } else { 88 elem_restr_jd_i_sur = NULL; 89 jac_data_sur = NULL; 90 } 91 92 // ---- CEED Vector 93 CeedInt loc_num_elem_sur; 94 CeedElemRestrictionGetNumElements(elem_restr_q_sur, &loc_num_elem_sur); 95 CeedVectorCreate(ceed, q_data_size_sur * loc_num_elem_sur * num_qpts_sur, &q_data_sur); 96 97 // ---- CEED Operator 98 // ----- CEED Operator for Setup (geometric factors) 99 CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur); 100 CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x_sur, ceed_data->basis_x_sur, CEED_VECTOR_ACTIVE); 101 CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x_sur, CEED_VECTOR_NONE); 102 CeedOperatorSetField(op_setup_sur, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 103 104 // ----- CEED Operator for Physics 105 CeedOperatorCreate(ceed, qf_apply_bc, NULL, NULL, &op_apply_bc); 106 CeedOperatorSetField(op_apply_bc, "q", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE); 107 CeedOperatorSetField(op_apply_bc, "Grad_q", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE); 108 CeedOperatorSetField(op_apply_bc, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, q_data_sur); 109 CeedOperatorSetField(op_apply_bc, "x", elem_restr_x_sur, ceed_data->basis_x_sur, ceed_data->x_coord); 110 CeedOperatorSetField(op_apply_bc, "v", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE); 111 if (elem_restr_jd_i_sur) CeedOperatorSetField(op_apply_bc, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_COLLOCATED, jac_data_sur); 112 113 if (qf_apply_bc_jacobian) { 114 CeedOperatorCreate(ceed, qf_apply_bc_jacobian, NULL, NULL, &op_apply_bc_jacobian); 115 CeedOperatorSetField(op_apply_bc_jacobian, "dq", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE); 116 CeedOperatorSetField(op_apply_bc_jacobian, "Grad_dq", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE); 117 CeedOperatorSetField(op_apply_bc_jacobian, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, q_data_sur); 118 CeedOperatorSetField(op_apply_bc_jacobian, "x", elem_restr_x_sur, ceed_data->basis_x_sur, ceed_data->x_coord); 119 CeedOperatorSetField(op_apply_bc_jacobian, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_COLLOCATED, jac_data_sur); 120 CeedOperatorSetField(op_apply_bc_jacobian, "v", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE); 121 } 122 123 // ----- Apply CEED operator for Setup 124 CeedOperatorApply(op_setup_sur, ceed_data->x_coord, q_data_sur, CEED_REQUEST_IMMEDIATE); 125 126 // ----- Apply Sub-Operator for Physics 127 CeedCompositeOperatorAddSub(*op_apply, op_apply_bc); 128 if (op_apply_bc_jacobian) CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_bc_jacobian); 129 130 // ----- Cleanup 131 CeedVectorDestroy(&q_data_sur); 132 CeedVectorDestroy(&jac_data_sur); 133 CeedElemRestrictionDestroy(&elem_restr_q_sur); 134 CeedElemRestrictionDestroy(&elem_restr_x_sur); 135 CeedElemRestrictionDestroy(&elem_restr_qd_i_sur); 136 CeedElemRestrictionDestroy(&elem_restr_jd_i_sur); 137 CeedOperatorDestroy(&op_setup_sur); 138 CeedOperatorDestroy(&op_apply_bc); 139 CeedOperatorDestroy(&op_apply_bc_jacobian); 140 141 PetscFunctionReturn(PETSC_SUCCESS); 142 } 143 144 // Utility function to create CEED Composite Operator for the entire domain 145 PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc, CeedData ceed_data, Physics phys, CeedOperator op_apply_vol, 146 CeedOperator op_apply_ijacobian_vol, CeedInt height, CeedInt P_sur, CeedInt Q_sur, CeedInt q_data_size_sur, 147 CeedInt jac_data_size_sur, CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) { 148 DMLabel domain_label; 149 PetscFunctionBeginUser; 150 151 // Create Composite Operaters 152 CeedCompositeOperatorCreate(ceed, op_apply); 153 if (op_apply_ijacobian) CeedCompositeOperatorCreate(ceed, op_apply_ijacobian); 154 155 // --Apply Sub-Operator for the volume 156 CeedCompositeOperatorAddSub(*op_apply, op_apply_vol); 157 if (op_apply_ijacobian) CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_ijacobian_vol); 158 159 // -- Create Sub-Operator for in/outflow BCs 160 if (phys->has_neumann || 1) { 161 // --- Setup 162 PetscCall(DMGetLabel(dm, "Face Sets", &domain_label)); 163 164 // --- Create Sub-Operator for inflow boundaries 165 for (CeedInt i = 0; i < bc->num_inflow; i++) { 166 PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->inflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, 167 ceed_data->qf_apply_inflow, ceed_data->qf_apply_inflow_jacobian, op_apply, op_apply_ijacobian)); 168 } 169 // --- Create Sub-Operator for outflow boundaries 170 for (CeedInt i = 0; i < bc->num_outflow; i++) { 171 PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->outflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, 172 ceed_data->qf_apply_outflow, ceed_data->qf_apply_outflow_jacobian, op_apply, op_apply_ijacobian)); 173 } 174 // --- Create Sub-Operator for freestream boundaries 175 for (CeedInt i = 0; i < bc->num_freestream; i++) { 176 PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->freestreams[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, 177 ceed_data->qf_apply_freestream, ceed_data->qf_apply_freestream_jacobian, op_apply, op_apply_ijacobian)); 178 } 179 } 180 181 // ----- Get Context Labels for Operator 182 CeedOperatorGetContextFieldLabel(*op_apply, "solution time", &phys->solution_time_label); 183 CeedOperatorGetContextFieldLabel(*op_apply, "timestep size", &phys->timestep_size_label); 184 185 PetscFunctionReturn(PETSC_SUCCESS); 186 } 187 188 PetscErrorCode SetupBCQFunctions(Ceed ceed, PetscInt dim_sur, PetscInt num_comp_x, PetscInt num_comp_q, PetscInt q_data_size_sur, 189 PetscInt jac_data_size_sur, ProblemQFunctionSpec apply_bc, ProblemQFunctionSpec apply_bc_jacobian, 190 CeedQFunction *qf_apply_bc, CeedQFunction *qf_apply_bc_jacobian) { 191 PetscFunctionBeginUser; 192 193 if (apply_bc.qfunction) { 194 CeedQFunctionCreateInterior(ceed, 1, apply_bc.qfunction, apply_bc.qfunction_loc, qf_apply_bc); 195 CeedQFunctionSetContext(*qf_apply_bc, apply_bc.qfunction_context); 196 CeedQFunctionAddInput(*qf_apply_bc, "q", num_comp_q, CEED_EVAL_INTERP); 197 CeedQFunctionAddInput(*qf_apply_bc, "Grad_q", num_comp_q * dim_sur, CEED_EVAL_GRAD); 198 CeedQFunctionAddInput(*qf_apply_bc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE); 199 CeedQFunctionAddInput(*qf_apply_bc, "x", num_comp_x, CEED_EVAL_INTERP); 200 CeedQFunctionAddOutput(*qf_apply_bc, "v", num_comp_q, CEED_EVAL_INTERP); 201 if (jac_data_size_sur) CeedQFunctionAddOutput(*qf_apply_bc, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE); 202 } 203 if (apply_bc_jacobian.qfunction) { 204 CeedQFunctionCreateInterior(ceed, 1, apply_bc_jacobian.qfunction, apply_bc_jacobian.qfunction_loc, qf_apply_bc_jacobian); 205 CeedQFunctionSetContext(*qf_apply_bc_jacobian, apply_bc_jacobian.qfunction_context); 206 CeedQFunctionAddInput(*qf_apply_bc_jacobian, "dq", num_comp_q, CEED_EVAL_INTERP); 207 CeedQFunctionAddInput(*qf_apply_bc_jacobian, "Grad_dq", num_comp_q * dim_sur, CEED_EVAL_GRAD); 208 CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface qdata", q_data_size_sur, CEED_EVAL_NONE); 209 CeedQFunctionAddInput(*qf_apply_bc_jacobian, "x", num_comp_x, CEED_EVAL_INTERP); 210 CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE); 211 CeedQFunctionAddOutput(*qf_apply_bc_jacobian, "v", num_comp_q, CEED_EVAL_INTERP); 212 } 213 PetscFunctionReturn(PETSC_SUCCESS); 214 } 215 216 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc) { 217 PetscFunctionBeginUser; 218 219 // ***************************************************************************** 220 // Set up CEED objects for the interior domain (volume) 221 // ***************************************************************************** 222 const PetscInt num_comp_q = 5; 223 const CeedInt dim = problem->dim, num_comp_x = problem->dim, q_data_size_vol = problem->q_data_size_vol, jac_data_size_vol = num_comp_q + 6 + 3, 224 P = app_ctx->degree + 1, Q = P + app_ctx->q_extra; 225 CeedElemRestriction elem_restr_jd_i; 226 CeedVector jac_data; 227 228 // ----------------------------------------------------------------------------- 229 // CEED Bases 230 // ----------------------------------------------------------------------------- 231 CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_q, P, Q, CEED_GAUSS, &ceed_data->basis_q); 232 CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, CEED_GAUSS, &ceed_data->basis_x); 233 CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, P, CEED_GAUSS_LOBATTO, &ceed_data->basis_xc); 234 235 // ----------------------------------------------------------------------------- 236 // CEED Restrictions 237 // ----------------------------------------------------------------------------- 238 // -- Create restriction 239 PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, 0, Q, q_data_size_vol, &ceed_data->elem_restr_q, &ceed_data->elem_restr_x, 240 &ceed_data->elem_restr_qd_i)); 241 242 PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, 0, Q, jac_data_size_vol, NULL, NULL, &elem_restr_jd_i)); 243 // -- Create E vectors 244 CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL); 245 CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed, NULL); 246 CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL); 247 248 // ----------------------------------------------------------------------------- 249 // CEED QFunctions 250 // ----------------------------------------------------------------------------- 251 // -- Create QFunction for quadrature data 252 CeedQFunctionCreateInterior(ceed, 1, problem->setup_vol.qfunction, problem->setup_vol.qfunction_loc, &ceed_data->qf_setup_vol); 253 if (problem->setup_vol.qfunction_context) { 254 CeedQFunctionSetContext(ceed_data->qf_setup_vol, problem->setup_vol.qfunction_context); 255 } 256 CeedQFunctionAddInput(ceed_data->qf_setup_vol, "dx", num_comp_x * dim, CEED_EVAL_GRAD); 257 CeedQFunctionAddInput(ceed_data->qf_setup_vol, "weight", 1, CEED_EVAL_WEIGHT); 258 CeedQFunctionAddOutput(ceed_data->qf_setup_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE); 259 260 // -- Create QFunction for ICs 261 CeedQFunctionCreateInterior(ceed, 1, problem->ics.qfunction, problem->ics.qfunction_loc, &ceed_data->qf_ics); 262 CeedQFunctionSetContext(ceed_data->qf_ics, problem->ics.qfunction_context); 263 CeedQFunctionAddInput(ceed_data->qf_ics, "x", num_comp_x, CEED_EVAL_INTERP); 264 CeedQFunctionAddInput(ceed_data->qf_ics, "dx", num_comp_x * dim, CEED_EVAL_GRAD); 265 CeedQFunctionAddOutput(ceed_data->qf_ics, "q0", num_comp_q, CEED_EVAL_NONE); 266 267 // -- Create QFunction for RHS 268 if (problem->apply_vol_rhs.qfunction) { 269 CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qfunction, problem->apply_vol_rhs.qfunction_loc, &ceed_data->qf_rhs_vol); 270 CeedQFunctionSetContext(ceed_data->qf_rhs_vol, problem->apply_vol_rhs.qfunction_context); 271 CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP); 272 CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD); 273 CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE); 274 CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP); 275 CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP); 276 CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD); 277 } 278 279 // -- Create QFunction for IFunction 280 if (problem->apply_vol_ifunction.qfunction) { 281 CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qfunction, problem->apply_vol_ifunction.qfunction_loc, 282 &ceed_data->qf_ifunction_vol); 283 CeedQFunctionSetContext(ceed_data->qf_ifunction_vol, problem->apply_vol_ifunction.qfunction_context); 284 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP); 285 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD); 286 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP); 287 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE); 288 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP); 289 CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP); 290 CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD); 291 CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE); 292 } 293 294 CeedQFunction qf_ijacobian_vol = NULL; 295 if (problem->apply_vol_ijacobian.qfunction) { 296 CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qfunction, problem->apply_vol_ijacobian.qfunction_loc, &qf_ijacobian_vol); 297 CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfunction_context); 298 CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP); 299 CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD); 300 CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE); 301 CeedQFunctionAddInput(qf_ijacobian_vol, "x", num_comp_x, CEED_EVAL_INTERP); 302 CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE); 303 CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP); 304 CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD); 305 } 306 307 // --------------------------------------------------------------------------- 308 // Element coordinates 309 // --------------------------------------------------------------------------- 310 // -- Create CEED vector 311 CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord, NULL); 312 313 // -- Copy PETSc vector in CEED vector 314 Vec X_loc; 315 { 316 DM cdm; 317 PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 318 if (cdm) { 319 PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc)); 320 } else { 321 PetscCall(DMGetCoordinatesLocal(dm, &X_loc)); 322 } 323 } 324 PetscCall(VecScale(X_loc, problem->dm_scale)); 325 PetscCall(VecCopyP2C(X_loc, ceed_data->x_coord)); 326 327 // ----------------------------------------------------------------------------- 328 // CEED vectors 329 // ----------------------------------------------------------------------------- 330 // -- Create CEED vector for geometric data 331 CeedElemRestrictionCreateVector(ceed_data->elem_restr_qd_i, &ceed_data->q_data, NULL); 332 CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL); 333 334 // ----------------------------------------------------------------------------- 335 // CEED Operators 336 // ----------------------------------------------------------------------------- 337 // -- Create CEED operator for quadrature data 338 CeedOperatorCreate(ceed, ceed_data->qf_setup_vol, NULL, NULL, &ceed_data->op_setup_vol); 339 CeedOperatorSetField(ceed_data->op_setup_vol, "dx", ceed_data->elem_restr_x, ceed_data->basis_x, CEED_VECTOR_ACTIVE); 340 CeedOperatorSetField(ceed_data->op_setup_vol, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x, CEED_VECTOR_NONE); 341 CeedOperatorSetField(ceed_data->op_setup_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 342 343 // -- Create CEED operator for ICs 344 CeedOperator op_ics; 345 CeedOperatorCreate(ceed, ceed_data->qf_ics, NULL, NULL, &op_ics); 346 CeedOperatorSetField(op_ics, "x", ceed_data->elem_restr_x, ceed_data->basis_xc, CEED_VECTOR_ACTIVE); 347 CeedOperatorSetField(op_ics, "dx", ceed_data->elem_restr_x, ceed_data->basis_xc, CEED_VECTOR_ACTIVE); 348 CeedOperatorSetField(op_ics, "q0", ceed_data->elem_restr_q, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 349 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 CeedOperatorDestroy(&op_ics); 352 353 // Create CEED operator for RHS 354 if (ceed_data->qf_rhs_vol) { 355 CeedOperator op; 356 CeedOperatorCreate(ceed, ceed_data->qf_rhs_vol, NULL, NULL, &op); 357 CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 358 CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 359 CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data); 360 CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord); 361 CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 362 CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 363 user->op_rhs_vol = op; 364 } 365 366 // -- CEED operator for IFunction 367 if (ceed_data->qf_ifunction_vol) { 368 CeedOperator op; 369 CeedOperatorCreate(ceed, ceed_data->qf_ifunction_vol, NULL, NULL, &op); 370 CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 371 CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 372 CeedOperatorSetField(op, "q dot", ceed_data->elem_restr_q, ceed_data->basis_q, user->q_dot_ceed); 373 CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data); 374 CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord); 375 CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 376 CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 377 CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_COLLOCATED, jac_data); 378 379 user->op_ifunction_vol = op; 380 } 381 382 CeedOperator op_ijacobian_vol = NULL; 383 if (qf_ijacobian_vol) { 384 CeedOperator op; 385 CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op); 386 CeedOperatorSetField(op, "dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 387 CeedOperatorSetField(op, "Grad_dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 388 CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data); 389 CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord); 390 CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_COLLOCATED, jac_data); 391 CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 392 CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 393 op_ijacobian_vol = op; 394 CeedQFunctionDestroy(&qf_ijacobian_vol); 395 } 396 397 // ***************************************************************************** 398 // Set up CEED objects for the exterior domain (surface) 399 // ***************************************************************************** 400 CeedInt height = 1, dim_sur = dim - height, P_sur = app_ctx->degree + 1, Q_sur = P_sur + app_ctx->q_extra; 401 const CeedInt q_data_size_sur = problem->q_data_size_sur, jac_data_size_sur = problem->jac_data_size_sur; 402 403 // ----------------------------------------------------------------------------- 404 // CEED Bases 405 // ----------------------------------------------------------------------------- 406 CeedBasisCreateTensorH1Lagrange(ceed, dim_sur, num_comp_q, P_sur, Q_sur, CEED_GAUSS, &ceed_data->basis_q_sur); 407 CeedBasisCreateTensorH1Lagrange(ceed, dim_sur, num_comp_x, 2, Q_sur, CEED_GAUSS, &ceed_data->basis_x_sur); 408 CeedBasisCreateTensorH1Lagrange(ceed, dim_sur, num_comp_x, 2, P_sur, CEED_GAUSS_LOBATTO, &ceed_data->basis_xc_sur); 409 410 // ----------------------------------------------------------------------------- 411 // CEED QFunctions 412 // ----------------------------------------------------------------------------- 413 // -- Create QFunction for quadrature data 414 CeedQFunctionCreateInterior(ceed, 1, problem->setup_sur.qfunction, problem->setup_sur.qfunction_loc, &ceed_data->qf_setup_sur); 415 if (problem->setup_sur.qfunction_context) { 416 CeedQFunctionSetContext(ceed_data->qf_setup_sur, problem->setup_sur.qfunction_context); 417 } 418 CeedQFunctionAddInput(ceed_data->qf_setup_sur, "dx", num_comp_x * dim_sur, CEED_EVAL_GRAD); 419 CeedQFunctionAddInput(ceed_data->qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT); 420 CeedQFunctionAddOutput(ceed_data->qf_setup_sur, "surface qdata", q_data_size_sur, CEED_EVAL_NONE); 421 422 PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_inflow, 423 problem->apply_inflow_jacobian, &ceed_data->qf_apply_inflow, &ceed_data->qf_apply_inflow_jacobian)); 424 PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_outflow, 425 problem->apply_outflow_jacobian, &ceed_data->qf_apply_outflow, &ceed_data->qf_apply_outflow_jacobian)); 426 PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_freestream, 427 problem->apply_freestream_jacobian, &ceed_data->qf_apply_freestream, &ceed_data->qf_apply_freestream_jacobian)); 428 429 // ***************************************************************************** 430 // CEED Operator Apply 431 // ***************************************************************************** 432 // -- Apply CEED Operator for the geometric data 433 CeedOperatorApply(ceed_data->op_setup_vol, ceed_data->x_coord, ceed_data->q_data, CEED_REQUEST_IMMEDIATE); 434 435 // -- Create and apply CEED Composite Operator for the entire domain 436 if (!user->phys->implicit) { // RHS 437 CeedOperator op_rhs; 438 PetscCall(CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys, user->op_rhs_vol, NULL, height, P_sur, Q_sur, q_data_size_sur, 0, &op_rhs, 439 NULL)); 440 PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_rhs, user->q_ceed, user->g_ceed, user->Q_loc, NULL, &user->op_rhs_ctx)); 441 CeedOperatorDestroy(&op_rhs); 442 } else { // IFunction 443 PetscCall(CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys, user->op_ifunction_vol, op_ijacobian_vol, height, P_sur, Q_sur, 444 q_data_size_sur, jac_data_size_sur, &user->op_ifunction, op_ijacobian_vol ? &user->op_ijacobian : NULL)); 445 if (user->op_ijacobian) { 446 CeedOperatorGetContextFieldLabel(user->op_ijacobian, "ijacobian time shift", &user->phys->ijacobian_time_shift_label); 447 } 448 if (problem->use_strong_bc_ceed) PetscCall(SetupStrongBC_Ceed(ceed, ceed_data, dm, user, problem, bc, Q_sur, q_data_size_sur)); 449 if (app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) PetscCall(SGS_DD_ModelSetup(ceed, user, ceed_data, problem)); 450 } 451 452 if (app_ctx->turb_spanstats_enable) PetscCall(TurbulenceStatisticsSetup(ceed, user, ceed_data, problem)); 453 if (app_ctx->diff_filter_monitor) PetscCall(DifferentialFilterSetup(ceed, user, ceed_data, problem)); 454 455 CeedElemRestrictionDestroy(&elem_restr_jd_i); 456 CeedOperatorDestroy(&op_ijacobian_vol); 457 CeedVectorDestroy(&jac_data); 458 PetscFunctionReturn(PETSC_SUCCESS); 459 } 460