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