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