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