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