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