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 // --- Create Sub-Operator for freestream boundaries 228 for (CeedInt i=0; i < bc->num_freestream; i++) { 229 PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, 230 bc->freestreams[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, 231 ceed_data->qf_apply_freestream, ceed_data->qf_apply_freestream_jacobian, 232 op_apply, op_apply_ijacobian)); 233 } 234 } 235 236 // ----- Get Context Labels for Operator 237 CeedOperatorContextGetFieldLabel(*op_apply, "solution time", 238 &phys->solution_time_label); 239 CeedOperatorContextGetFieldLabel(*op_apply, "timestep size", 240 &phys->timestep_size_label); 241 242 PetscFunctionReturn(0); 243 } 244 245 PetscErrorCode SetupBCQFunctions(Ceed ceed, PetscInt dim_sur, 246 PetscInt num_comp_x, PetscInt num_comp_q, 247 PetscInt q_data_size_sur, PetscInt jac_data_size_sur, 248 ProblemQFunctionSpec apply_bc, ProblemQFunctionSpec apply_bc_jacobian, 249 CeedQFunction *qf_apply_bc, CeedQFunction *qf_apply_bc_jacobian) { 250 PetscFunctionBeginUser; 251 252 if (apply_bc.qfunction) { 253 // *INDENT-OFF* 254 CeedQFunctionCreateInterior(ceed, 1, apply_bc.qfunction, apply_bc.qfunction_loc, qf_apply_bc); 255 CeedQFunctionSetContext(*qf_apply_bc, apply_bc.qfunction_context); 256 CeedQFunctionContextDestroy(&apply_bc.qfunction_context); 257 CeedQFunctionAddInput(*qf_apply_bc, "q", num_comp_q, CEED_EVAL_INTERP); 258 CeedQFunctionAddInput(*qf_apply_bc, "Grad_q", num_comp_q*dim_sur, CEED_EVAL_GRAD); 259 CeedQFunctionAddInput(*qf_apply_bc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE); 260 CeedQFunctionAddInput(*qf_apply_bc, "x", num_comp_x, CEED_EVAL_INTERP); 261 CeedQFunctionAddOutput(*qf_apply_bc, "v", num_comp_q, CEED_EVAL_INTERP); 262 // *INDENT-ON* 263 if (jac_data_size_sur) 264 CeedQFunctionAddOutput(*qf_apply_bc, "surface jacobian data", jac_data_size_sur, 265 CEED_EVAL_NONE); 266 } 267 if (apply_bc_jacobian.qfunction) { 268 // *INDENT-OFF* 269 CeedQFunctionCreateInterior(ceed, 1, apply_bc_jacobian.qfunction, 270 apply_bc_jacobian.qfunction_loc, qf_apply_bc_jacobian); 271 CeedQFunctionSetContext(*qf_apply_bc_jacobian, apply_bc_jacobian.qfunction_context); 272 CeedQFunctionContextDestroy(&apply_bc_jacobian.qfunction_context); 273 CeedQFunctionAddInput(*qf_apply_bc_jacobian, "dq", num_comp_q, CEED_EVAL_INTERP); 274 CeedQFunctionAddInput(*qf_apply_bc_jacobian, "Grad_dq", num_comp_q*dim_sur, CEED_EVAL_GRAD); 275 CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface qdata", q_data_size_sur, CEED_EVAL_NONE); 276 CeedQFunctionAddInput(*qf_apply_bc_jacobian, "x", num_comp_x, CEED_EVAL_INTERP); 277 CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface jacobian data", 278 jac_data_size_sur, CEED_EVAL_NONE); 279 CeedQFunctionAddOutput(*qf_apply_bc_jacobian, "v", num_comp_q, CEED_EVAL_INTERP); 280 // *INDENT-ON* 281 } 282 PetscFunctionReturn(0); 283 } 284 285 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, 286 AppCtx app_ctx, ProblemData *problem, SimpleBC bc) { 287 PetscErrorCode ierr; 288 PetscFunctionBeginUser; 289 290 // ***************************************************************************** 291 // Set up CEED objects for the interior domain (volume) 292 // ***************************************************************************** 293 const PetscInt num_comp_q = 5; 294 const CeedInt dim = problem->dim, 295 num_comp_x = problem->dim, 296 q_data_size_vol = problem->q_data_size_vol, 297 jac_data_size_vol = num_comp_q + 6 + 3, 298 P = app_ctx->degree + 1, 299 Q = P + app_ctx->q_extra; 300 CeedElemRestriction elem_restr_jd_i; 301 CeedVector jac_data; 302 303 // ----------------------------------------------------------------------------- 304 // CEED Bases 305 // ----------------------------------------------------------------------------- 306 CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_q, P, Q, CEED_GAUSS, 307 &ceed_data->basis_q); 308 CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, CEED_GAUSS, 309 &ceed_data->basis_x); 310 CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, P, 311 CEED_GAUSS_LOBATTO, &ceed_data->basis_xc); 312 313 // ----------------------------------------------------------------------------- 314 // CEED Restrictions 315 // ----------------------------------------------------------------------------- 316 // -- Create restriction 317 ierr = GetRestrictionForDomain(ceed, dm, 0, 0, 0, Q, q_data_size_vol, 318 &ceed_data->elem_restr_q, &ceed_data->elem_restr_x, 319 &ceed_data->elem_restr_qd_i); CHKERRQ(ierr); 320 321 ierr = GetRestrictionForDomain(ceed, dm, 0, 0, 0, Q, jac_data_size_vol, 322 NULL, NULL, 323 &elem_restr_jd_i); CHKERRQ(ierr); 324 // -- Create E vectors 325 CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL); 326 CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed, 327 NULL); 328 CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL); 329 330 // ----------------------------------------------------------------------------- 331 // CEED QFunctions 332 // ----------------------------------------------------------------------------- 333 // -- Create QFunction for quadrature data 334 CeedQFunctionCreateInterior(ceed, 1, problem->setup_vol.qfunction, 335 problem->setup_vol.qfunction_loc, 336 &ceed_data->qf_setup_vol); 337 if (problem->setup_vol.qfunction_context) { 338 CeedQFunctionSetContext(ceed_data->qf_setup_vol, 339 problem->setup_vol.qfunction_context); 340 CeedQFunctionContextDestroy(&problem->setup_vol.qfunction_context); 341 } 342 CeedQFunctionAddInput(ceed_data->qf_setup_vol, "dx", num_comp_x*dim, 343 CEED_EVAL_GRAD); 344 CeedQFunctionAddInput(ceed_data->qf_setup_vol, "weight", 1, CEED_EVAL_WEIGHT); 345 CeedQFunctionAddOutput(ceed_data->qf_setup_vol, "qdata", q_data_size_vol, 346 CEED_EVAL_NONE); 347 348 // -- Create QFunction for ICs 349 CeedQFunctionCreateInterior(ceed, 1, problem->ics.qfunction, 350 problem->ics.qfunction_loc, 351 &ceed_data->qf_ics); 352 CeedQFunctionSetContext(ceed_data->qf_ics, problem->ics.qfunction_context); 353 CeedQFunctionContextDestroy(&problem->ics.qfunction_context); 354 CeedQFunctionAddInput(ceed_data->qf_ics, "x", num_comp_x, CEED_EVAL_INTERP); 355 CeedQFunctionAddOutput(ceed_data->qf_ics, "q0", num_comp_q, CEED_EVAL_NONE); 356 357 // -- Create QFunction for RHS 358 if (problem->apply_vol_rhs.qfunction) { 359 CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qfunction, 360 problem->apply_vol_rhs.qfunction_loc, &ceed_data->qf_rhs_vol); 361 CeedQFunctionSetContext(ceed_data->qf_rhs_vol, 362 problem->apply_vol_rhs.qfunction_context); 363 CeedQFunctionContextDestroy(&problem->apply_vol_rhs.qfunction_context); 364 CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP); 365 CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "Grad_q", num_comp_q*dim, 366 CEED_EVAL_GRAD); 367 CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "qdata", q_data_size_vol, 368 CEED_EVAL_NONE); 369 CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP); 370 CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "v", num_comp_q, 371 CEED_EVAL_INTERP); 372 CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "Grad_v", num_comp_q*dim, 373 CEED_EVAL_GRAD); 374 } 375 376 // -- Create QFunction for IFunction 377 if (problem->apply_vol_ifunction.qfunction) { 378 CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qfunction, 379 problem->apply_vol_ifunction.qfunction_loc, &ceed_data->qf_ifunction_vol); 380 CeedQFunctionSetContext(ceed_data->qf_ifunction_vol, 381 problem->apply_vol_ifunction.qfunction_context); 382 CeedQFunctionContextDestroy(&problem->apply_vol_ifunction.qfunction_context); 383 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q", num_comp_q, 384 CEED_EVAL_INTERP); 385 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "Grad_q", num_comp_q*dim, 386 CEED_EVAL_GRAD); 387 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q dot", num_comp_q, 388 CEED_EVAL_INTERP); 389 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "qdata", q_data_size_vol, 390 CEED_EVAL_NONE); 391 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "x", num_comp_x, 392 CEED_EVAL_INTERP); 393 CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "v", num_comp_q, 394 CEED_EVAL_INTERP); 395 CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "Grad_v", num_comp_q*dim, 396 CEED_EVAL_GRAD); 397 CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "jac_data", 398 jac_data_size_vol, CEED_EVAL_NONE); 399 } 400 401 CeedQFunction qf_ijacobian_vol = NULL; 402 if (problem->apply_vol_ijacobian.qfunction) { 403 CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qfunction, 404 problem->apply_vol_ijacobian.qfunction_loc, &qf_ijacobian_vol); 405 CeedQFunctionSetContext(qf_ijacobian_vol, 406 problem->apply_vol_ijacobian.qfunction_context); 407 CeedQFunctionContextDestroy(&problem->apply_vol_ijacobian.qfunction_context); 408 CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, 409 CEED_EVAL_INTERP); 410 CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q*dim, 411 CEED_EVAL_GRAD); 412 CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol, 413 CEED_EVAL_NONE); 414 CeedQFunctionAddInput(qf_ijacobian_vol, "x", num_comp_x, 415 CEED_EVAL_INTERP); 416 CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", 417 jac_data_size_vol, CEED_EVAL_NONE); 418 CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, 419 CEED_EVAL_INTERP); 420 CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q*dim, 421 CEED_EVAL_GRAD); 422 } 423 424 // --------------------------------------------------------------------------- 425 // Element coordinates 426 // --------------------------------------------------------------------------- 427 // -- Create CEED vector 428 CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord, 429 NULL); 430 431 // -- Copy PETSc vector in CEED vector 432 Vec X_loc; 433 const PetscScalar *X_loc_array; 434 { 435 DM cdm; 436 ierr = DMGetCellCoordinateDM(dm, &cdm); CHKERRQ(ierr); 437 if (cdm) {ierr = DMGetCellCoordinatesLocal(dm, &X_loc); CHKERRQ(ierr);} 438 else {ierr = DMGetCoordinatesLocal(dm, &X_loc); CHKERRQ(ierr);} 439 } 440 ierr = VecScale(X_loc, problem->dm_scale); CHKERRQ(ierr); 441 ierr = VecGetArrayRead(X_loc, &X_loc_array); CHKERRQ(ierr); 442 CeedVectorSetArray(ceed_data->x_coord, CEED_MEM_HOST, CEED_COPY_VALUES, 443 (PetscScalar *)X_loc_array); 444 ierr = VecRestoreArrayRead(X_loc, &X_loc_array); CHKERRQ(ierr); 445 446 // ----------------------------------------------------------------------------- 447 // CEED vectors 448 // ----------------------------------------------------------------------------- 449 // -- Create CEED vector for geometric data 450 CeedInt num_qpts_vol; 451 PetscInt loc_num_elem_vol; 452 CeedBasisGetNumQuadraturePoints(ceed_data->basis_q, &num_qpts_vol); 453 CeedElemRestrictionGetNumElements(ceed_data->elem_restr_q, &loc_num_elem_vol); 454 CeedVectorCreate(ceed, q_data_size_vol*loc_num_elem_vol*num_qpts_vol, 455 &ceed_data->q_data); 456 457 CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL); 458 // ----------------------------------------------------------------------------- 459 // CEED Operators 460 // ----------------------------------------------------------------------------- 461 // -- Create CEED operator for quadrature data 462 CeedOperatorCreate(ceed, ceed_data->qf_setup_vol, NULL, NULL, 463 &ceed_data->op_setup_vol); 464 CeedOperatorSetField(ceed_data->op_setup_vol, "dx", ceed_data->elem_restr_x, 465 ceed_data->basis_x, CEED_VECTOR_ACTIVE); 466 CeedOperatorSetField(ceed_data->op_setup_vol, "weight", 467 CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x, CEED_VECTOR_NONE); 468 CeedOperatorSetField(ceed_data->op_setup_vol, "qdata", 469 ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 470 471 // -- Create CEED operator for ICs 472 CeedOperatorCreate(ceed, ceed_data->qf_ics, NULL, NULL, &ceed_data->op_ics); 473 CeedOperatorSetField(ceed_data->op_ics, "x", ceed_data->elem_restr_x, 474 ceed_data->basis_xc, CEED_VECTOR_ACTIVE); 475 CeedOperatorSetField(ceed_data->op_ics, "q0", ceed_data->elem_restr_q, 476 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 477 CeedOperatorContextGetFieldLabel(ceed_data->op_ics, "evaluation time", 478 &user->phys->ics_time_label); 479 480 // Create CEED operator for RHS 481 if (ceed_data->qf_rhs_vol) { 482 CeedOperator op; 483 CeedOperatorCreate(ceed, ceed_data->qf_rhs_vol, NULL, NULL, &op); 484 CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, 485 CEED_VECTOR_ACTIVE); 486 CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, 487 CEED_VECTOR_ACTIVE); 488 CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, 489 CEED_BASIS_COLLOCATED, ceed_data->q_data); 490 CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, 491 ceed_data->x_coord); 492 CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, 493 CEED_VECTOR_ACTIVE); 494 CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, 495 CEED_VECTOR_ACTIVE); 496 user->op_rhs_vol = op; 497 } 498 499 // -- CEED operator for IFunction 500 if (ceed_data->qf_ifunction_vol) { 501 CeedOperator op; 502 CeedOperatorCreate(ceed, ceed_data->qf_ifunction_vol, NULL, NULL, &op); 503 CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, 504 CEED_VECTOR_ACTIVE); 505 CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, 506 CEED_VECTOR_ACTIVE); 507 CeedOperatorSetField(op, "q dot", ceed_data->elem_restr_q, ceed_data->basis_q, 508 user->q_dot_ceed); 509 CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, 510 CEED_BASIS_COLLOCATED, ceed_data->q_data); 511 CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, 512 ceed_data->x_coord); 513 CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, 514 CEED_VECTOR_ACTIVE); 515 CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, 516 CEED_VECTOR_ACTIVE); 517 CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, 518 CEED_BASIS_COLLOCATED, jac_data); 519 520 user->op_ifunction_vol = op; 521 } 522 523 CeedOperator op_ijacobian_vol = NULL; 524 if (qf_ijacobian_vol) { 525 CeedOperator op; 526 CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op); 527 CeedOperatorSetField(op, "dq", ceed_data->elem_restr_q, ceed_data->basis_q, 528 CEED_VECTOR_ACTIVE); 529 CeedOperatorSetField(op, "Grad_dq", ceed_data->elem_restr_q, ceed_data->basis_q, 530 CEED_VECTOR_ACTIVE); 531 CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, 532 CEED_BASIS_COLLOCATED, ceed_data->q_data); 533 CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, 534 ceed_data->x_coord); 535 CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, 536 CEED_BASIS_COLLOCATED, jac_data); 537 CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, 538 CEED_VECTOR_ACTIVE); 539 CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, 540 CEED_VECTOR_ACTIVE); 541 op_ijacobian_vol = op; 542 CeedQFunctionDestroy(&qf_ijacobian_vol); 543 } 544 545 // ***************************************************************************** 546 // Set up CEED objects for the exterior domain (surface) 547 // ***************************************************************************** 548 CeedInt height = 1, 549 dim_sur = dim - height, 550 P_sur = app_ctx->degree + 1, 551 Q_sur = P_sur + app_ctx->q_extra; 552 const CeedInt q_data_size_sur = problem->q_data_size_sur, 553 jac_data_size_sur = problem->jac_data_size_sur; 554 555 // ----------------------------------------------------------------------------- 556 // CEED Bases 557 // ----------------------------------------------------------------------------- 558 CeedBasisCreateTensorH1Lagrange(ceed, dim_sur, num_comp_q, P_sur, Q_sur, 559 CEED_GAUSS, &ceed_data->basis_q_sur); 560 CeedBasisCreateTensorH1Lagrange(ceed, dim_sur, num_comp_x, 2, Q_sur, CEED_GAUSS, 561 &ceed_data->basis_x_sur); 562 CeedBasisCreateTensorH1Lagrange(ceed, dim_sur, num_comp_x, 2, P_sur, 563 CEED_GAUSS_LOBATTO, &ceed_data->basis_xc_sur); 564 565 // ----------------------------------------------------------------------------- 566 // CEED QFunctions 567 // ----------------------------------------------------------------------------- 568 // -- Create QFunction for quadrature data 569 CeedQFunctionCreateInterior(ceed, 1, problem->setup_sur.qfunction, 570 problem->setup_sur.qfunction_loc, 571 &ceed_data->qf_setup_sur); 572 if (problem->setup_sur.qfunction_context) { 573 CeedQFunctionSetContext(ceed_data->qf_setup_sur, 574 problem->setup_sur.qfunction_context); 575 CeedQFunctionContextDestroy(&problem->setup_sur.qfunction_context); 576 } 577 CeedQFunctionAddInput(ceed_data->qf_setup_sur, "dx", num_comp_x*dim_sur, 578 CEED_EVAL_GRAD); 579 CeedQFunctionAddInput(ceed_data->qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT); 580 CeedQFunctionAddOutput(ceed_data->qf_setup_sur, "surface qdata", 581 q_data_size_sur, CEED_EVAL_NONE); 582 583 PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, 584 q_data_size_sur, jac_data_size_sur, 585 problem->apply_inflow, problem->apply_inflow_jacobian, 586 &ceed_data->qf_apply_inflow, &ceed_data->qf_apply_inflow_jacobian)); 587 588 PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, 589 q_data_size_sur, jac_data_size_sur, 590 problem->apply_outflow, problem->apply_outflow_jacobian, 591 &ceed_data->qf_apply_outflow, &ceed_data->qf_apply_outflow_jacobian)); 592 593 PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, 594 q_data_size_sur, jac_data_size_sur, 595 problem->apply_freestream, problem->apply_freestream_jacobian, 596 &ceed_data->qf_apply_freestream, &ceed_data->qf_apply_freestream_jacobian)); 597 598 // ***************************************************************************** 599 // CEED Operator Apply 600 // ***************************************************************************** 601 // -- Apply CEED Operator for the geometric data 602 CeedOperatorApply(ceed_data->op_setup_vol, ceed_data->x_coord, 603 ceed_data->q_data, CEED_REQUEST_IMMEDIATE); 604 605 // -- Create and apply CEED Composite Operator for the entire domain 606 if (!user->phys->implicit) { // RHS 607 ierr = CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys, 608 user->op_rhs_vol, NULL, height, P_sur, Q_sur, 609 q_data_size_sur, 0, 610 &user->op_rhs, NULL); CHKERRQ(ierr); 611 } else { // IFunction 612 ierr = CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys, 613 user->op_ifunction_vol, op_ijacobian_vol, 614 height, P_sur, Q_sur, 615 q_data_size_sur, jac_data_size_sur, 616 &user->op_ifunction, 617 op_ijacobian_vol ? &user->op_ijacobian : NULL); CHKERRQ(ierr); 618 if (user->op_ijacobian) { 619 CeedOperatorContextGetFieldLabel(user->op_ijacobian, "ijacobian time shift", 620 &user->phys->ijacobian_time_shift_label); 621 } 622 if (problem->use_dirichlet_ceed) { 623 PetscCall(SetupStrongBC_Ceed(ceed, ceed_data, dm, user, app_ctx, problem, bc, 624 Q_sur, q_data_size_sur)); 625 } 626 627 } 628 CeedElemRestrictionDestroy(&elem_restr_jd_i); 629 CeedOperatorDestroy(&op_ijacobian_vol); 630 CeedVectorDestroy(&jac_data); 631 PetscFunctionReturn(0); 632 } 633