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 CeedElemRestriction elem_restr_tmp; 48 PetscErrorCode ierr; 49 PetscFunctionBeginUser; 50 51 ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr); 52 dim -= height; 53 Q_dim = CeedIntPow(Q, dim); 54 ierr = CreateRestrictionFromPlex(ceed, dm, height, domain_label, value, 55 &elem_restr_tmp); 56 CHKERRQ(ierr); 57 if (elem_restr_q) *elem_restr_q = elem_restr_tmp; 58 if (elem_restr_x) { 59 ierr = DMGetCellCoordinateDM(dm, &dm_coord); CHKERRQ(ierr); 60 if (!dm_coord) { 61 ierr = DMGetCoordinateDM(dm, &dm_coord); CHKERRQ(ierr); 62 } 63 ierr = DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL); 64 CHKERRQ(ierr); 65 ierr = CreateRestrictionFromPlex(ceed, dm_coord, height, domain_label, value, 66 elem_restr_x); 67 CHKERRQ(ierr); 68 } 69 if (elem_restr_qd_i) { 70 CeedElemRestrictionGetNumElements(elem_restr_tmp, &loc_num_elem); 71 CeedElemRestrictionCreateStrided(ceed, loc_num_elem, Q_dim, 72 q_data_size, q_data_size*loc_num_elem*Q_dim, 73 CEED_STRIDES_BACKEND, elem_restr_qd_i); 74 } 75 if (!elem_restr_q) CeedElemRestrictionDestroy(&elem_restr_tmp); 76 PetscFunctionReturn(0); 77 } 78 79 PetscErrorCode AddBCSubOperator(Ceed ceed, DM dm, CeedData ceed_data, 80 DMLabel domain_label, PetscInt value, 81 CeedInt height, CeedInt Q_sur, 82 CeedInt q_data_size_sur, CeedInt jac_data_size_sur, 83 CeedQFunction qf_apply_bc, CeedQFunction qf_apply_bc_jacobian, 84 CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) { 85 CeedVector q_data_sur, jac_data_sur; 86 CeedOperator op_setup_sur, op_apply_bc, 87 op_apply_bc_jacobian = NULL; 88 CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur, elem_restr_qd_i_sur, 89 elem_restr_jd_i_sur; 90 CeedInt num_qpts_sur; 91 PetscFunctionBeginUser; 92 93 // --- Get number of quadrature points for the boundaries 94 CeedBasisGetNumQuadraturePoints(ceed_data->basis_q_sur, &num_qpts_sur); 95 96 97 // ---- CEED Restriction 98 PetscCall(GetRestrictionForDomain(ceed, dm, height, domain_label, value, Q_sur, 99 q_data_size_sur, &elem_restr_q_sur, &elem_restr_x_sur, &elem_restr_qd_i_sur)); 100 if (jac_data_size_sur > 0) { 101 // State-dependent data will be passed from residual to Jacobian. This will be collocated. 102 PetscCall(GetRestrictionForDomain(ceed, dm, height, domain_label, value, Q_sur, 103 jac_data_size_sur, NULL, NULL, &elem_restr_jd_i_sur)); 104 CeedElemRestrictionCreateVector(elem_restr_jd_i_sur, &jac_data_sur, NULL); 105 } else { 106 elem_restr_jd_i_sur = NULL; 107 jac_data_sur = NULL; 108 } 109 110 // ---- CEED Vector 111 PetscInt loc_num_elem_sur; 112 CeedElemRestrictionGetNumElements(elem_restr_q_sur, &loc_num_elem_sur); 113 CeedVectorCreate(ceed, q_data_size_sur*loc_num_elem_sur*num_qpts_sur, 114 &q_data_sur); 115 116 // ---- CEED Operator 117 // ----- CEED Operator for Setup (geometric factors) 118 CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur); 119 CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x_sur, 120 ceed_data->basis_x_sur, CEED_VECTOR_ACTIVE); 121 CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, 122 ceed_data->basis_x_sur, CEED_VECTOR_NONE); 123 CeedOperatorSetField(op_setup_sur, "surface qdata", elem_restr_qd_i_sur, 124 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 125 126 // ----- CEED Operator for Physics 127 CeedOperatorCreate(ceed, qf_apply_bc, NULL, NULL, &op_apply_bc); 128 CeedOperatorSetField(op_apply_bc, "q", elem_restr_q_sur, ceed_data->basis_q_sur, 129 CEED_VECTOR_ACTIVE); 130 CeedOperatorSetField(op_apply_bc, "Grad_q", elem_restr_q_sur, 131 ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE); 132 CeedOperatorSetField(op_apply_bc, "surface qdata", elem_restr_qd_i_sur, 133 CEED_BASIS_COLLOCATED, q_data_sur); 134 CeedOperatorSetField(op_apply_bc, "x", elem_restr_x_sur, ceed_data->basis_x_sur, 135 ceed_data->x_coord); 136 CeedOperatorSetField(op_apply_bc, "v", elem_restr_q_sur, ceed_data->basis_q_sur, 137 CEED_VECTOR_ACTIVE); 138 if (elem_restr_jd_i_sur) 139 CeedOperatorSetField(op_apply_bc, "surface jacobian data", 140 elem_restr_jd_i_sur, 141 CEED_BASIS_COLLOCATED, jac_data_sur); 142 143 if (qf_apply_bc_jacobian) { 144 CeedOperatorCreate(ceed, qf_apply_bc_jacobian, NULL, NULL, 145 &op_apply_bc_jacobian); 146 CeedOperatorSetField(op_apply_bc_jacobian, "dq", elem_restr_q_sur, 147 ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE); 148 CeedOperatorSetField(op_apply_bc_jacobian, "Grad_dq", elem_restr_q_sur, 149 ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE); 150 CeedOperatorSetField(op_apply_bc_jacobian, "surface qdata", elem_restr_qd_i_sur, 151 CEED_BASIS_COLLOCATED, q_data_sur); 152 CeedOperatorSetField(op_apply_bc_jacobian, "x", elem_restr_x_sur, 153 ceed_data->basis_x_sur, ceed_data->x_coord); 154 CeedOperatorSetField(op_apply_bc_jacobian, "surface jacobian data", 155 elem_restr_jd_i_sur, CEED_BASIS_COLLOCATED, jac_data_sur); 156 CeedOperatorSetField(op_apply_bc_jacobian, "v", elem_restr_q_sur, 157 ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE); 158 } 159 160 // ----- Apply CEED operator for Setup 161 CeedOperatorApply(op_setup_sur, ceed_data->x_coord, q_data_sur, 162 CEED_REQUEST_IMMEDIATE); 163 164 // ----- Apply Sub-Operator for Physics 165 CeedCompositeOperatorAddSub(*op_apply, op_apply_bc); 166 if (op_apply_bc_jacobian) 167 CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_bc_jacobian); 168 169 // ----- Cleanup 170 CeedVectorDestroy(&q_data_sur); 171 CeedVectorDestroy(&jac_data_sur); 172 CeedElemRestrictionDestroy(&elem_restr_q_sur); 173 CeedElemRestrictionDestroy(&elem_restr_x_sur); 174 CeedElemRestrictionDestroy(&elem_restr_qd_i_sur); 175 CeedElemRestrictionDestroy(&elem_restr_jd_i_sur); 176 CeedOperatorDestroy(&op_setup_sur); 177 CeedOperatorDestroy(&op_apply_bc); 178 CeedOperatorDestroy(&op_apply_bc_jacobian); 179 180 PetscFunctionReturn(0); 181 } 182 183 // Utility function to create CEED Composite Operator for the entire domain 184 PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc, 185 CeedData ceed_data, Physics phys, 186 CeedOperator op_apply_vol, 187 CeedOperator op_apply_ijacobian_vol, 188 CeedInt height, 189 CeedInt P_sur, CeedInt Q_sur, 190 CeedInt q_data_size_sur, CeedInt jac_data_size_sur, 191 CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) { 192 DMLabel domain_label; 193 PetscErrorCode ierr; 194 PetscFunctionBeginUser; 195 196 // Create Composite Operaters 197 CeedCompositeOperatorCreate(ceed, op_apply); 198 if (op_apply_ijacobian) 199 CeedCompositeOperatorCreate(ceed, op_apply_ijacobian); 200 201 // --Apply Sub-Operator for the volume 202 CeedCompositeOperatorAddSub(*op_apply, op_apply_vol); 203 if (op_apply_ijacobian) 204 CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_ijacobian_vol); 205 206 // -- Create Sub-Operator for in/outflow BCs 207 if (phys->has_neumann || 1) { 208 // --- Setup 209 ierr = DMGetLabel(dm, "Face Sets", &domain_label); CHKERRQ(ierr); 210 211 // --- Create Sub-Operator for inflow boundaries 212 for (CeedInt i=0; i < bc->num_inflow; i++) { 213 PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->inflows[i], 214 height, Q_sur, q_data_size_sur, jac_data_size_sur, 215 ceed_data->qf_apply_inflow, ceed_data->qf_apply_inflow_jacobian, 216 op_apply, op_apply_ijacobian)); 217 } 218 219 // --- Create Sub-Operator for outflow boundaries 220 for (CeedInt i=0; i < bc->num_outflow; i++) { 221 PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->outflows[i], 222 height, Q_sur, q_data_size_sur, jac_data_size_sur, 223 ceed_data->qf_apply_outflow, ceed_data->qf_apply_outflow_jacobian, 224 op_apply, op_apply_ijacobian)); 225 } 226 } 227 228 // ----- Get Context Labels for Operator 229 CeedOperatorContextGetFieldLabel(*op_apply, "solution time", 230 &phys->solution_time_label); 231 CeedOperatorContextGetFieldLabel(*op_apply, "timestep size", 232 &phys->timestep_size_label); 233 234 PetscFunctionReturn(0); 235 } 236 237 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, 238 AppCtx app_ctx, ProblemData *problem, SimpleBC bc) { 239 PetscErrorCode ierr; 240 PetscFunctionBeginUser; 241 242 // ***************************************************************************** 243 // Set up CEED objects for the interior domain (volume) 244 // ***************************************************************************** 245 const PetscInt num_comp_q = 5; 246 const CeedInt dim = problem->dim, 247 num_comp_x = problem->dim, 248 q_data_size_vol = problem->q_data_size_vol, 249 jac_data_size_vol = num_comp_q + 6 + 3, 250 P = app_ctx->degree + 1, 251 Q = P + app_ctx->q_extra; 252 CeedElemRestriction elem_restr_jd_i; 253 CeedVector jac_data; 254 255 // ----------------------------------------------------------------------------- 256 // CEED Bases 257 // ----------------------------------------------------------------------------- 258 CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_q, P, Q, CEED_GAUSS, 259 &ceed_data->basis_q); 260 CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, CEED_GAUSS, 261 &ceed_data->basis_x); 262 CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, P, 263 CEED_GAUSS_LOBATTO, &ceed_data->basis_xc); 264 265 // ----------------------------------------------------------------------------- 266 // CEED Restrictions 267 // ----------------------------------------------------------------------------- 268 // -- Create restriction 269 ierr = GetRestrictionForDomain(ceed, dm, 0, 0, 0, Q, q_data_size_vol, 270 &ceed_data->elem_restr_q, &ceed_data->elem_restr_x, 271 &ceed_data->elem_restr_qd_i); CHKERRQ(ierr); 272 273 ierr = GetRestrictionForDomain(ceed, dm, 0, 0, 0, Q, jac_data_size_vol, 274 NULL, NULL, 275 &elem_restr_jd_i); CHKERRQ(ierr); 276 // -- Create E vectors 277 CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL); 278 CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed, 279 NULL); 280 CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL); 281 282 // ----------------------------------------------------------------------------- 283 // CEED QFunctions 284 // ----------------------------------------------------------------------------- 285 // -- Create QFunction for quadrature data 286 CeedQFunctionCreateInterior(ceed, 1, problem->setup_vol.qfunction, 287 problem->setup_vol.qfunction_loc, 288 &ceed_data->qf_setup_vol); 289 if (problem->setup_vol.qfunction_context) { 290 CeedQFunctionSetContext(ceed_data->qf_setup_vol, 291 problem->setup_vol.qfunction_context); 292 CeedQFunctionContextDestroy(&problem->setup_vol.qfunction_context); 293 } 294 CeedQFunctionAddInput(ceed_data->qf_setup_vol, "dx", num_comp_x*dim, 295 CEED_EVAL_GRAD); 296 CeedQFunctionAddInput(ceed_data->qf_setup_vol, "weight", 1, CEED_EVAL_WEIGHT); 297 CeedQFunctionAddOutput(ceed_data->qf_setup_vol, "qdata", q_data_size_vol, 298 CEED_EVAL_NONE); 299 300 // -- Create QFunction for ICs 301 CeedQFunctionCreateInterior(ceed, 1, problem->ics.qfunction, 302 problem->ics.qfunction_loc, 303 &ceed_data->qf_ics); 304 CeedQFunctionSetContext(ceed_data->qf_ics, problem->ics.qfunction_context); 305 CeedQFunctionContextDestroy(&problem->ics.qfunction_context); 306 CeedQFunctionAddInput(ceed_data->qf_ics, "x", num_comp_x, CEED_EVAL_INTERP); 307 CeedQFunctionAddOutput(ceed_data->qf_ics, "q0", num_comp_q, CEED_EVAL_NONE); 308 309 // -- Create QFunction for RHS 310 if (problem->apply_vol_rhs.qfunction) { 311 CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qfunction, 312 problem->apply_vol_rhs.qfunction_loc, &ceed_data->qf_rhs_vol); 313 CeedQFunctionSetContext(ceed_data->qf_rhs_vol, 314 problem->apply_vol_rhs.qfunction_context); 315 CeedQFunctionContextDestroy(&problem->apply_vol_rhs.qfunction_context); 316 CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP); 317 CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "Grad_q", num_comp_q*dim, 318 CEED_EVAL_GRAD); 319 CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "qdata", q_data_size_vol, 320 CEED_EVAL_NONE); 321 CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP); 322 CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "v", num_comp_q, 323 CEED_EVAL_INTERP); 324 CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "Grad_v", num_comp_q*dim, 325 CEED_EVAL_GRAD); 326 } 327 328 // -- Create QFunction for IFunction 329 if (problem->apply_vol_ifunction.qfunction) { 330 CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qfunction, 331 problem->apply_vol_ifunction.qfunction_loc, &ceed_data->qf_ifunction_vol); 332 CeedQFunctionSetContext(ceed_data->qf_ifunction_vol, 333 problem->apply_vol_ifunction.qfunction_context); 334 CeedQFunctionContextDestroy(&problem->apply_vol_ifunction.qfunction_context); 335 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q", num_comp_q, 336 CEED_EVAL_INTERP); 337 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "Grad_q", num_comp_q*dim, 338 CEED_EVAL_GRAD); 339 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q dot", num_comp_q, 340 CEED_EVAL_INTERP); 341 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "qdata", q_data_size_vol, 342 CEED_EVAL_NONE); 343 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "x", num_comp_x, 344 CEED_EVAL_INTERP); 345 CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "v", num_comp_q, 346 CEED_EVAL_INTERP); 347 CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "Grad_v", num_comp_q*dim, 348 CEED_EVAL_GRAD); 349 CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "jac_data", 350 jac_data_size_vol, CEED_EVAL_NONE); 351 } 352 353 CeedQFunction qf_ijacobian_vol = NULL; 354 if (problem->apply_vol_ijacobian.qfunction) { 355 CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qfunction, 356 problem->apply_vol_ijacobian.qfunction_loc, &qf_ijacobian_vol); 357 CeedQFunctionSetContext(qf_ijacobian_vol, 358 problem->apply_vol_ijacobian.qfunction_context); 359 CeedQFunctionContextDestroy(&problem->apply_vol_ijacobian.qfunction_context); 360 CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, 361 CEED_EVAL_INTERP); 362 CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q*dim, 363 CEED_EVAL_GRAD); 364 CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol, 365 CEED_EVAL_NONE); 366 CeedQFunctionAddInput(qf_ijacobian_vol, "x", num_comp_x, 367 CEED_EVAL_INTERP); 368 CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", 369 jac_data_size_vol, CEED_EVAL_NONE); 370 CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, 371 CEED_EVAL_INTERP); 372 CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q*dim, 373 CEED_EVAL_GRAD); 374 } 375 376 // --------------------------------------------------------------------------- 377 // Element coordinates 378 // --------------------------------------------------------------------------- 379 // -- Create CEED vector 380 CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord, 381 NULL); 382 383 // -- Copy PETSc vector in CEED vector 384 Vec X_loc; 385 const PetscScalar *X_loc_array; 386 { 387 DM cdm; 388 ierr = DMGetCellCoordinateDM(dm, &cdm); CHKERRQ(ierr); 389 if (cdm) {ierr = DMGetCellCoordinatesLocal(dm, &X_loc); CHKERRQ(ierr);} 390 else {ierr = DMGetCoordinatesLocal(dm, &X_loc); CHKERRQ(ierr);} 391 } 392 ierr = VecScale(X_loc, problem->dm_scale); CHKERRQ(ierr); 393 ierr = VecGetArrayRead(X_loc, &X_loc_array); CHKERRQ(ierr); 394 CeedVectorSetArray(ceed_data->x_coord, CEED_MEM_HOST, CEED_COPY_VALUES, 395 (PetscScalar *)X_loc_array); 396 ierr = VecRestoreArrayRead(X_loc, &X_loc_array); CHKERRQ(ierr); 397 398 // ----------------------------------------------------------------------------- 399 // CEED vectors 400 // ----------------------------------------------------------------------------- 401 // -- Create CEED vector for geometric data 402 CeedInt num_qpts_vol; 403 PetscInt loc_num_elem_vol; 404 CeedBasisGetNumQuadraturePoints(ceed_data->basis_q, &num_qpts_vol); 405 CeedElemRestrictionGetNumElements(ceed_data->elem_restr_q, &loc_num_elem_vol); 406 CeedVectorCreate(ceed, q_data_size_vol*loc_num_elem_vol*num_qpts_vol, 407 &ceed_data->q_data); 408 409 CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL); 410 // ----------------------------------------------------------------------------- 411 // CEED Operators 412 // ----------------------------------------------------------------------------- 413 // -- Create CEED operator for quadrature data 414 CeedOperatorCreate(ceed, ceed_data->qf_setup_vol, NULL, NULL, 415 &ceed_data->op_setup_vol); 416 CeedOperatorSetField(ceed_data->op_setup_vol, "dx", ceed_data->elem_restr_x, 417 ceed_data->basis_x, CEED_VECTOR_ACTIVE); 418 CeedOperatorSetField(ceed_data->op_setup_vol, "weight", 419 CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x, CEED_VECTOR_NONE); 420 CeedOperatorSetField(ceed_data->op_setup_vol, "qdata", 421 ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 422 423 // -- Create CEED operator for ICs 424 CeedOperatorCreate(ceed, ceed_data->qf_ics, NULL, NULL, &ceed_data->op_ics); 425 CeedOperatorSetField(ceed_data->op_ics, "x", ceed_data->elem_restr_x, 426 ceed_data->basis_xc, CEED_VECTOR_ACTIVE); 427 CeedOperatorSetField(ceed_data->op_ics, "q0", ceed_data->elem_restr_q, 428 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 429 CeedOperatorContextGetFieldLabel(ceed_data->op_ics, "evaluation time", 430 &user->phys->ics_time_label); 431 432 // Create CEED operator for RHS 433 if (ceed_data->qf_rhs_vol) { 434 CeedOperator op; 435 CeedOperatorCreate(ceed, ceed_data->qf_rhs_vol, NULL, NULL, &op); 436 CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, 437 CEED_VECTOR_ACTIVE); 438 CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, 439 CEED_VECTOR_ACTIVE); 440 CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, 441 CEED_BASIS_COLLOCATED, ceed_data->q_data); 442 CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, 443 ceed_data->x_coord); 444 CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, 445 CEED_VECTOR_ACTIVE); 446 CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, 447 CEED_VECTOR_ACTIVE); 448 user->op_rhs_vol = op; 449 } 450 451 // -- CEED operator for IFunction 452 if (ceed_data->qf_ifunction_vol) { 453 CeedOperator op; 454 CeedOperatorCreate(ceed, ceed_data->qf_ifunction_vol, NULL, NULL, &op); 455 CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, 456 CEED_VECTOR_ACTIVE); 457 CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, 458 CEED_VECTOR_ACTIVE); 459 CeedOperatorSetField(op, "q dot", ceed_data->elem_restr_q, ceed_data->basis_q, 460 user->q_dot_ceed); 461 CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, 462 CEED_BASIS_COLLOCATED, ceed_data->q_data); 463 CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, 464 ceed_data->x_coord); 465 CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, 466 CEED_VECTOR_ACTIVE); 467 CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, 468 CEED_VECTOR_ACTIVE); 469 CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, 470 CEED_BASIS_COLLOCATED, jac_data); 471 472 user->op_ifunction_vol = op; 473 } 474 475 CeedOperator op_ijacobian_vol = NULL; 476 if (qf_ijacobian_vol) { 477 CeedOperator op; 478 CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op); 479 CeedOperatorSetField(op, "dq", ceed_data->elem_restr_q, ceed_data->basis_q, 480 CEED_VECTOR_ACTIVE); 481 CeedOperatorSetField(op, "Grad_dq", ceed_data->elem_restr_q, ceed_data->basis_q, 482 CEED_VECTOR_ACTIVE); 483 CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, 484 CEED_BASIS_COLLOCATED, ceed_data->q_data); 485 CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, 486 ceed_data->x_coord); 487 CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, 488 CEED_BASIS_COLLOCATED, jac_data); 489 CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, 490 CEED_VECTOR_ACTIVE); 491 CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, 492 CEED_VECTOR_ACTIVE); 493 op_ijacobian_vol = op; 494 CeedQFunctionDestroy(&qf_ijacobian_vol); 495 } 496 497 // ***************************************************************************** 498 // Set up CEED objects for the exterior domain (surface) 499 // ***************************************************************************** 500 CeedInt height = 1, 501 dim_sur = dim - height, 502 P_sur = app_ctx->degree + 1, 503 Q_sur = P_sur + app_ctx->q_extra; 504 const CeedInt q_data_size_sur = problem->q_data_size_sur, 505 jac_data_size_sur = problem->jac_data_size_sur; 506 507 // ----------------------------------------------------------------------------- 508 // CEED Bases 509 // ----------------------------------------------------------------------------- 510 CeedBasisCreateTensorH1Lagrange(ceed, dim_sur, num_comp_q, P_sur, Q_sur, 511 CEED_GAUSS, &ceed_data->basis_q_sur); 512 CeedBasisCreateTensorH1Lagrange(ceed, dim_sur, num_comp_x, 2, Q_sur, CEED_GAUSS, 513 &ceed_data->basis_x_sur); 514 CeedBasisCreateTensorH1Lagrange(ceed, dim_sur, num_comp_x, 2, P_sur, 515 CEED_GAUSS_LOBATTO, &ceed_data->basis_xc_sur); 516 517 // ----------------------------------------------------------------------------- 518 // CEED QFunctions 519 // ----------------------------------------------------------------------------- 520 // -- Create QFunction for quadrature data 521 CeedQFunctionCreateInterior(ceed, 1, problem->setup_sur.qfunction, 522 problem->setup_sur.qfunction_loc, 523 &ceed_data->qf_setup_sur); 524 if (problem->setup_sur.qfunction_context) { 525 CeedQFunctionSetContext(ceed_data->qf_setup_sur, 526 problem->setup_sur.qfunction_context); 527 CeedQFunctionContextDestroy(&problem->setup_sur.qfunction_context); 528 } 529 CeedQFunctionAddInput(ceed_data->qf_setup_sur, "dx", num_comp_x*dim_sur, 530 CEED_EVAL_GRAD); 531 CeedQFunctionAddInput(ceed_data->qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT); 532 CeedQFunctionAddOutput(ceed_data->qf_setup_sur, "surface qdata", 533 q_data_size_sur, CEED_EVAL_NONE); 534 535 // -- Creat QFunction for inflow boundaries 536 if (problem->apply_inflow.qfunction) { 537 CeedQFunctionCreateInterior(ceed, 1, problem->apply_inflow.qfunction, 538 problem->apply_inflow.qfunction_loc, &ceed_data->qf_apply_inflow); 539 CeedQFunctionSetContext(ceed_data->qf_apply_inflow, 540 problem->apply_inflow.qfunction_context); 541 CeedQFunctionContextDestroy(&problem->apply_inflow.qfunction_context); 542 CeedQFunctionAddInput(ceed_data->qf_apply_inflow, "q", num_comp_q, 543 CEED_EVAL_INTERP); 544 CeedQFunctionAddInput(ceed_data->qf_apply_inflow, "Grad_q", num_comp_q*(dim-1), 545 CEED_EVAL_GRAD); 546 CeedQFunctionAddInput(ceed_data->qf_apply_inflow, "surface qdata", 547 q_data_size_sur, CEED_EVAL_NONE); 548 CeedQFunctionAddInput(ceed_data->qf_apply_inflow, "x", num_comp_x, 549 CEED_EVAL_INTERP); 550 CeedQFunctionAddOutput(ceed_data->qf_apply_inflow, "v", num_comp_q, 551 CEED_EVAL_INTERP); 552 if (jac_data_size_sur) 553 CeedQFunctionAddOutput(ceed_data->qf_apply_inflow, "surface jacobian data", 554 jac_data_size_sur, 555 CEED_EVAL_NONE); 556 } 557 if (problem->apply_inflow_jacobian.qfunction) { 558 CeedQFunctionCreateInterior(ceed, 1, problem->apply_inflow_jacobian.qfunction, 559 problem->apply_inflow_jacobian.qfunction_loc, 560 &ceed_data->qf_apply_inflow_jacobian); 561 CeedQFunctionSetContext(ceed_data->qf_apply_inflow_jacobian, 562 problem->apply_inflow_jacobian.qfunction_context); 563 CeedQFunctionContextDestroy(&problem->apply_inflow_jacobian.qfunction_context); 564 CeedQFunctionAddInput(ceed_data->qf_apply_inflow_jacobian, "dq", num_comp_q, 565 CEED_EVAL_INTERP); 566 CeedQFunctionAddInput(ceed_data->qf_apply_inflow_jacobian, "Grad_dq", 567 num_comp_q*dim_sur, CEED_EVAL_GRAD); 568 CeedQFunctionAddInput(ceed_data->qf_apply_inflow_jacobian, "surface qdata", 569 q_data_size_sur, CEED_EVAL_NONE); 570 CeedQFunctionAddInput(ceed_data->qf_apply_inflow_jacobian, "x", num_comp_x, 571 CEED_EVAL_INTERP); 572 CeedQFunctionAddInput(ceed_data->qf_apply_inflow_jacobian, 573 "surface jacobian data", 574 jac_data_size_sur, CEED_EVAL_NONE); 575 CeedQFunctionAddOutput(ceed_data->qf_apply_inflow_jacobian, "v", num_comp_q, 576 CEED_EVAL_INTERP); 577 } 578 579 // -- Creat QFunction for outflow boundaries 580 if (problem->apply_outflow.qfunction) { 581 CeedQFunctionCreateInterior(ceed, 1, problem->apply_outflow.qfunction, 582 problem->apply_outflow.qfunction_loc, &ceed_data->qf_apply_outflow); 583 CeedQFunctionSetContext(ceed_data->qf_apply_outflow, 584 problem->apply_outflow.qfunction_context); 585 CeedQFunctionContextDestroy(&problem->apply_outflow.qfunction_context); 586 CeedQFunctionAddInput(ceed_data->qf_apply_outflow, "q", num_comp_q, 587 CEED_EVAL_INTERP); 588 CeedQFunctionAddInput(ceed_data->qf_apply_outflow, "Grad_q", num_comp_q*(dim-1), 589 CEED_EVAL_GRAD); 590 CeedQFunctionAddInput(ceed_data->qf_apply_outflow, "surface qdata", 591 q_data_size_sur, CEED_EVAL_NONE); 592 CeedQFunctionAddInput(ceed_data->qf_apply_outflow, "x", num_comp_x, 593 CEED_EVAL_INTERP); 594 CeedQFunctionAddOutput(ceed_data->qf_apply_outflow, "v", num_comp_q, 595 CEED_EVAL_INTERP); 596 if (jac_data_size_sur) 597 CeedQFunctionAddOutput(ceed_data->qf_apply_outflow, "surface jacobian data", 598 jac_data_size_sur, 599 CEED_EVAL_NONE); 600 } 601 if (problem->apply_outflow_jacobian.qfunction) { 602 CeedQFunctionCreateInterior(ceed, 1, problem->apply_outflow_jacobian.qfunction, 603 problem->apply_outflow_jacobian.qfunction_loc, 604 &ceed_data->qf_apply_outflow_jacobian); 605 CeedQFunctionSetContext(ceed_data->qf_apply_outflow_jacobian, 606 problem->apply_outflow_jacobian.qfunction_context); 607 CeedQFunctionContextDestroy(&problem->apply_outflow_jacobian.qfunction_context); 608 CeedQFunctionAddInput(ceed_data->qf_apply_outflow_jacobian, "dq", num_comp_q, 609 CEED_EVAL_INTERP); 610 CeedQFunctionAddInput(ceed_data->qf_apply_outflow_jacobian, "Grad_dq", 611 num_comp_q*dim_sur, CEED_EVAL_GRAD); 612 CeedQFunctionAddInput(ceed_data->qf_apply_outflow_jacobian, "surface qdata", 613 q_data_size_sur, CEED_EVAL_NONE); 614 CeedQFunctionAddInput(ceed_data->qf_apply_outflow_jacobian, "x", num_comp_x, 615 CEED_EVAL_INTERP); 616 CeedQFunctionAddInput(ceed_data->qf_apply_outflow_jacobian, 617 "surface jacobian data", 618 jac_data_size_sur, CEED_EVAL_NONE); 619 CeedQFunctionAddOutput(ceed_data->qf_apply_outflow_jacobian, "v", num_comp_q, 620 CEED_EVAL_INTERP); 621 } 622 623 // ***************************************************************************** 624 // CEED Operator Apply 625 // ***************************************************************************** 626 // -- Apply CEED Operator for the geometric data 627 CeedOperatorApply(ceed_data->op_setup_vol, ceed_data->x_coord, 628 ceed_data->q_data, CEED_REQUEST_IMMEDIATE); 629 630 // -- Create and apply CEED Composite Operator for the entire domain 631 if (!user->phys->implicit) { // RHS 632 ierr = CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys, 633 user->op_rhs_vol, NULL, height, P_sur, Q_sur, 634 q_data_size_sur, 0, 635 &user->op_rhs, NULL); CHKERRQ(ierr); 636 } else { // IFunction 637 ierr = CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys, 638 user->op_ifunction_vol, op_ijacobian_vol, 639 height, P_sur, Q_sur, 640 q_data_size_sur, jac_data_size_sur, 641 &user->op_ifunction, 642 op_ijacobian_vol ? &user->op_ijacobian : NULL); CHKERRQ(ierr); 643 if (user->op_ijacobian) { 644 CeedOperatorContextGetFieldLabel(user->op_ijacobian, "ijacobian time shift", 645 &user->phys->ijacobian_time_shift_label); 646 } 647 if (problem->use_dirichlet_ceed) { 648 PetscCall(SetupStrongBC_Ceed(ceed, ceed_data, dm, user, app_ctx, problem, bc, 649 Q_sur, q_data_size_sur)); 650 } 651 652 } 653 CeedElemRestrictionDestroy(&elem_restr_jd_i); 654 CeedOperatorDestroy(&op_ijacobian_vol); 655 CeedVectorDestroy(&jac_data); 656 PetscFunctionReturn(0); 657 } 658