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 <ceed.h> 12 #include <petscdmplex.h> 13 14 #include <petscds.h> 15 #include "../navierstokes.h" 16 17 // Utility function to create local CEED restriction 18 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt label_value, PetscInt dm_field, 19 CeedElemRestriction *elem_restr) { 20 PetscInt num_elem, elem_size, num_dof, num_comp, *elem_restr_offsets_petsc; 21 CeedInt *elem_restr_offsets_ceed; 22 23 PetscFunctionBeginUser; 24 PetscCall( 25 DMPlexGetLocalOffsets(dm, domain_label, label_value, height, dm_field, &num_elem, &elem_size, &num_comp, &num_dof, &elem_restr_offsets_petsc)); 26 27 PetscCall(IntArrayP2C(num_elem * elem_size, &elem_restr_offsets_petsc, &elem_restr_offsets_ceed)); 28 CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp, 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES, elem_restr_offsets_ceed, elem_restr); 29 PetscCall(PetscFree(elem_restr_offsets_ceed)); 30 31 PetscFunctionReturn(PETSC_SUCCESS); 32 } 33 34 // Utility function to get Ceed Restriction for each domain 35 PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, PetscInt label_value, PetscInt dm_field, CeedInt Q, 36 CeedInt q_data_size, CeedElemRestriction *elem_restr_q, CeedElemRestriction *elem_restr_x, 37 CeedElemRestriction *elem_restr_qd_i) { 38 DM dm_coord; 39 CeedInt loc_num_elem; 40 PetscInt dim; 41 CeedElemRestriction elem_restr_tmp; 42 PetscFunctionBeginUser; 43 44 PetscCall(DMGetDimension(dm, &dim)); 45 dim -= height; 46 PetscCall(CreateRestrictionFromPlex(ceed, dm, height, domain_label, label_value, dm_field, &elem_restr_tmp)); 47 if (elem_restr_q) *elem_restr_q = elem_restr_tmp; 48 if (elem_restr_x) { 49 PetscCall(DMGetCellCoordinateDM(dm, &dm_coord)); 50 if (!dm_coord) { 51 PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 52 } 53 PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL)); 54 PetscCall(CreateRestrictionFromPlex(ceed, dm_coord, height, domain_label, label_value, dm_field, elem_restr_x)); 55 } 56 if (elem_restr_qd_i) { 57 CeedElemRestrictionGetNumElements(elem_restr_tmp, &loc_num_elem); 58 CeedElemRestrictionCreateStrided(ceed, loc_num_elem, Q, q_data_size, q_data_size * loc_num_elem * Q, CEED_STRIDES_BACKEND, elem_restr_qd_i); 59 } 60 if (!elem_restr_q) CeedElemRestrictionDestroy(&elem_restr_tmp); 61 PetscFunctionReturn(PETSC_SUCCESS); 62 } 63 64 PetscErrorCode AddBCSubOperator(Ceed ceed, DM dm, CeedData ceed_data, DMLabel domain_label, PetscInt label_value, CeedInt height, CeedInt Q_sur, 65 CeedInt q_data_size_sur, CeedInt jac_data_size_sur, CeedQFunction qf_apply_bc, CeedQFunction qf_apply_bc_jacobian, 66 CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) { 67 CeedVector q_data_sur, jac_data_sur; 68 CeedOperator op_setup_sur, op_apply_bc, op_apply_bc_jacobian = NULL; 69 CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur, elem_restr_qd_i_sur, elem_restr_jd_i_sur; 70 CeedInt num_qpts_sur; 71 PetscFunctionBeginUser; 72 73 // --- Get number of quadrature points for the boundaries 74 CeedBasisGetNumQuadraturePoints(ceed_data->basis_q_sur, &num_qpts_sur); 75 76 // ---- CEED Restriction 77 PetscCall(GetRestrictionForDomain(ceed, dm, height, domain_label, label_value, 0, num_qpts_sur, q_data_size_sur, &elem_restr_q_sur, 78 &elem_restr_x_sur, &elem_restr_qd_i_sur)); 79 if (jac_data_size_sur > 0) { 80 // State-dependent data will be passed from residual to Jacobian. This will be collocated. 81 PetscCall( 82 GetRestrictionForDomain(ceed, dm, height, domain_label, label_value, 0, num_qpts_sur, jac_data_size_sur, NULL, NULL, &elem_restr_jd_i_sur)); 83 CeedElemRestrictionCreateVector(elem_restr_jd_i_sur, &jac_data_sur, NULL); 84 } else { 85 elem_restr_jd_i_sur = NULL; 86 jac_data_sur = NULL; 87 } 88 89 // ---- CEED Vector 90 CeedInt loc_num_elem_sur; 91 CeedElemRestrictionGetNumElements(elem_restr_q_sur, &loc_num_elem_sur); 92 CeedVectorCreate(ceed, q_data_size_sur * loc_num_elem_sur * num_qpts_sur, &q_data_sur); 93 94 // ---- CEED Operator 95 // ----- CEED Operator for Setup (geometric factors) 96 CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur); 97 CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x_sur, ceed_data->basis_x_sur, CEED_VECTOR_ACTIVE); 98 CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x_sur, CEED_VECTOR_NONE); 99 CeedOperatorSetField(op_setup_sur, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 100 101 // ----- CEED Operator for Physics 102 CeedOperatorCreate(ceed, qf_apply_bc, NULL, NULL, &op_apply_bc); 103 CeedOperatorSetField(op_apply_bc, "q", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE); 104 CeedOperatorSetField(op_apply_bc, "Grad_q", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE); 105 CeedOperatorSetField(op_apply_bc, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, q_data_sur); 106 CeedOperatorSetField(op_apply_bc, "x", elem_restr_x_sur, ceed_data->basis_x_sur, ceed_data->x_coord); 107 CeedOperatorSetField(op_apply_bc, "v", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE); 108 if (elem_restr_jd_i_sur) CeedOperatorSetField(op_apply_bc, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_COLLOCATED, jac_data_sur); 109 110 if (qf_apply_bc_jacobian) { 111 CeedOperatorCreate(ceed, qf_apply_bc_jacobian, NULL, NULL, &op_apply_bc_jacobian); 112 CeedOperatorSetField(op_apply_bc_jacobian, "dq", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE); 113 CeedOperatorSetField(op_apply_bc_jacobian, "Grad_dq", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE); 114 CeedOperatorSetField(op_apply_bc_jacobian, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, q_data_sur); 115 CeedOperatorSetField(op_apply_bc_jacobian, "x", elem_restr_x_sur, ceed_data->basis_x_sur, ceed_data->x_coord); 116 CeedOperatorSetField(op_apply_bc_jacobian, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_COLLOCATED, jac_data_sur); 117 CeedOperatorSetField(op_apply_bc_jacobian, "v", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE); 118 } 119 120 // ----- Apply CEED operator for Setup 121 CeedOperatorApply(op_setup_sur, ceed_data->x_coord, q_data_sur, CEED_REQUEST_IMMEDIATE); 122 123 // ----- Apply Sub-Operator for Physics 124 CeedCompositeOperatorAddSub(*op_apply, op_apply_bc); 125 if (op_apply_bc_jacobian) CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_bc_jacobian); 126 127 // ----- Cleanup 128 CeedVectorDestroy(&q_data_sur); 129 CeedVectorDestroy(&jac_data_sur); 130 CeedElemRestrictionDestroy(&elem_restr_q_sur); 131 CeedElemRestrictionDestroy(&elem_restr_x_sur); 132 CeedElemRestrictionDestroy(&elem_restr_qd_i_sur); 133 CeedElemRestrictionDestroy(&elem_restr_jd_i_sur); 134 CeedOperatorDestroy(&op_setup_sur); 135 CeedOperatorDestroy(&op_apply_bc); 136 CeedOperatorDestroy(&op_apply_bc_jacobian); 137 138 PetscFunctionReturn(PETSC_SUCCESS); 139 } 140 141 // Utility function to create CEED Composite Operator for the entire domain 142 PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc, CeedData ceed_data, Physics phys, CeedOperator op_apply_vol, 143 CeedOperator op_apply_ijacobian_vol, CeedInt height, CeedInt P_sur, CeedInt Q_sur, CeedInt q_data_size_sur, 144 CeedInt jac_data_size_sur, CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) { 145 DMLabel domain_label; 146 PetscFunctionBeginUser; 147 148 // Create Composite Operaters 149 CeedCompositeOperatorCreate(ceed, op_apply); 150 if (op_apply_ijacobian) CeedCompositeOperatorCreate(ceed, op_apply_ijacobian); 151 152 // --Apply Sub-Operator for the volume 153 CeedCompositeOperatorAddSub(*op_apply, op_apply_vol); 154 if (op_apply_ijacobian) CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_ijacobian_vol); 155 156 // -- Create Sub-Operator for in/outflow BCs 157 if (phys->has_neumann || 1) { 158 // --- Setup 159 PetscCall(DMGetLabel(dm, "Face Sets", &domain_label)); 160 161 // --- Create Sub-Operator for inflow boundaries 162 for (CeedInt i = 0; i < bc->num_inflow; i++) { 163 PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->inflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, 164 ceed_data->qf_apply_inflow, ceed_data->qf_apply_inflow_jacobian, op_apply, op_apply_ijacobian)); 165 } 166 // --- Create Sub-Operator for outflow boundaries 167 for (CeedInt i = 0; i < bc->num_outflow; i++) { 168 PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->outflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, 169 ceed_data->qf_apply_outflow, ceed_data->qf_apply_outflow_jacobian, op_apply, op_apply_ijacobian)); 170 } 171 // --- Create Sub-Operator for freestream boundaries 172 for (CeedInt i = 0; i < bc->num_freestream; i++) { 173 PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->freestreams[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, 174 ceed_data->qf_apply_freestream, ceed_data->qf_apply_freestream_jacobian, op_apply, op_apply_ijacobian)); 175 } 176 } 177 178 // ----- Get Context Labels for Operator 179 CeedOperatorGetContextFieldLabel(*op_apply, "solution time", &phys->solution_time_label); 180 CeedOperatorGetContextFieldLabel(*op_apply, "timestep size", &phys->timestep_size_label); 181 182 PetscFunctionReturn(PETSC_SUCCESS); 183 } 184 185 PetscErrorCode SetupBCQFunctions(Ceed ceed, PetscInt dim_sur, PetscInt num_comp_x, PetscInt num_comp_q, PetscInt q_data_size_sur, 186 PetscInt jac_data_size_sur, ProblemQFunctionSpec apply_bc, ProblemQFunctionSpec apply_bc_jacobian, 187 CeedQFunction *qf_apply_bc, CeedQFunction *qf_apply_bc_jacobian) { 188 PetscFunctionBeginUser; 189 190 if (apply_bc.qfunction) { 191 CeedQFunctionCreateInterior(ceed, 1, apply_bc.qfunction, apply_bc.qfunction_loc, qf_apply_bc); 192 CeedQFunctionSetContext(*qf_apply_bc, apply_bc.qfunction_context); 193 CeedQFunctionAddInput(*qf_apply_bc, "q", num_comp_q, CEED_EVAL_INTERP); 194 CeedQFunctionAddInput(*qf_apply_bc, "Grad_q", num_comp_q * dim_sur, CEED_EVAL_GRAD); 195 CeedQFunctionAddInput(*qf_apply_bc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE); 196 CeedQFunctionAddInput(*qf_apply_bc, "x", num_comp_x, CEED_EVAL_INTERP); 197 CeedQFunctionAddOutput(*qf_apply_bc, "v", num_comp_q, CEED_EVAL_INTERP); 198 if (jac_data_size_sur) CeedQFunctionAddOutput(*qf_apply_bc, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE); 199 } 200 if (apply_bc_jacobian.qfunction) { 201 CeedQFunctionCreateInterior(ceed, 1, apply_bc_jacobian.qfunction, apply_bc_jacobian.qfunction_loc, qf_apply_bc_jacobian); 202 CeedQFunctionSetContext(*qf_apply_bc_jacobian, apply_bc_jacobian.qfunction_context); 203 CeedQFunctionAddInput(*qf_apply_bc_jacobian, "dq", num_comp_q, CEED_EVAL_INTERP); 204 CeedQFunctionAddInput(*qf_apply_bc_jacobian, "Grad_dq", num_comp_q * dim_sur, CEED_EVAL_GRAD); 205 CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface qdata", q_data_size_sur, CEED_EVAL_NONE); 206 CeedQFunctionAddInput(*qf_apply_bc_jacobian, "x", num_comp_x, CEED_EVAL_INTERP); 207 CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE); 208 CeedQFunctionAddOutput(*qf_apply_bc_jacobian, "v", num_comp_q, CEED_EVAL_INTERP); 209 } 210 PetscFunctionReturn(PETSC_SUCCESS); 211 } 212 213 // ----------------------------------------------------------------------------- 214 // Convert DM field to DS field 215 // ----------------------------------------------------------------------------- 216 PetscErrorCode DMFieldToDSField(DM dm, DMLabel domain_label, PetscInt dm_field, PetscInt *ds_field) { 217 PetscDS ds; 218 IS field_is; 219 const PetscInt *fields; 220 PetscInt num_fields; 221 222 PetscFunctionBeginUser; 223 224 // Translate dm_field to ds_field 225 PetscCall(DMGetRegionDS(dm, domain_label, &field_is, &ds, NULL)); 226 PetscCall(ISGetIndices(field_is, &fields)); 227 PetscCall(ISGetSize(field_is, &num_fields)); 228 for (PetscInt i = 0; i < num_fields; i++) { 229 if (dm_field == fields[i]) { 230 *ds_field = i; 231 break; 232 } 233 } 234 PetscCall(ISRestoreIndices(field_is, &fields)); 235 236 if (*ds_field == -1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Could not find dm_field %" PetscInt_FMT " in DS", dm_field); 237 238 PetscFunctionReturn(PETSC_SUCCESS); 239 } 240 241 // ----------------------------------------------------------------------------- 242 // Utility function - convert from DMPolytopeType to CeedElemTopology 243 // ----------------------------------------------------------------------------- 244 CeedElemTopology ElemTopologyP2C(DMPolytopeType cell_type) { 245 switch (cell_type) { 246 case DM_POLYTOPE_TRIANGLE: 247 return CEED_TOPOLOGY_TRIANGLE; 248 case DM_POLYTOPE_QUADRILATERAL: 249 return CEED_TOPOLOGY_QUAD; 250 case DM_POLYTOPE_TETRAHEDRON: 251 return CEED_TOPOLOGY_TET; 252 case DM_POLYTOPE_HEXAHEDRON: 253 return CEED_TOPOLOGY_HEX; 254 default: 255 return 0; 256 } 257 } 258 259 // ----------------------------------------------------------------------------- 260 // Create libCEED Basis from PetscTabulation 261 // ----------------------------------------------------------------------------- 262 PetscErrorCode BasisCreateFromTabulation(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt face, PetscFE fe, 263 PetscTabulation basis_tabulation, PetscQuadrature quadrature, CeedBasis *basis) { 264 PetscInt first_point; 265 PetscInt ids[1] = {label_value}; 266 DMLabel depth_label; 267 DMPolytopeType cell_type; 268 CeedElemTopology elem_topo; 269 PetscScalar *q_points, *interp, *grad; 270 const PetscScalar *q_weights; 271 PetscDualSpace dual_space; 272 PetscInt num_dual_basis_vectors; 273 PetscInt dim, num_comp, P, Q; 274 275 PetscFunctionBeginUser; 276 277 // General basis information 278 PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 279 PetscCall(PetscFEGetNumComponents(fe, &num_comp)); 280 PetscCall(PetscFEGetDualSpace(fe, &dual_space)); 281 PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors)); 282 P = num_dual_basis_vectors / num_comp; 283 284 // Use depth label if no domain label present 285 if (!domain_label) { 286 PetscInt depth; 287 288 PetscCall(DMPlexGetDepth(dm, &depth)); 289 PetscCall(DMPlexGetDepthLabel(dm, &depth_label)); 290 ids[0] = depth - height; 291 } 292 293 // Get cell interp, grad, and quadrature data 294 PetscCall(DMGetFirstLabeledPoint(dm, dm, domain_label ? domain_label : depth_label, 1, ids, height, &first_point, NULL)); 295 PetscCall(DMPlexGetCellType(dm, first_point, &cell_type)); 296 elem_topo = ElemTopologyP2C(cell_type); 297 if (!elem_topo) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DMPlex topology not supported"); 298 { 299 size_t q_points_size; 300 const PetscScalar *q_points_petsc; 301 PetscInt q_dim; 302 303 PetscCall(PetscQuadratureGetData(quadrature, &q_dim, NULL, &Q, &q_points_petsc, &q_weights)); 304 q_points_size = Q * dim * sizeof(CeedScalar); 305 PetscCall(PetscCalloc(q_points_size, &q_points)); 306 for (PetscInt q = 0; q < Q; q++) { 307 for (PetscInt d = 0; d < q_dim; d++) q_points[q * dim + d] = q_points_petsc[q * q_dim + d]; 308 } 309 } 310 311 { // Convert to libCEED orientation 312 PetscBool is_simplex = PETSC_FALSE; 313 IS permutation = NULL; 314 const PetscInt *permutation_indices; 315 316 PetscCall(DMPlexIsSimplex(dm, &is_simplex)); 317 if (!is_simplex) { 318 PetscSection section; 319 320 // -- Get permutation 321 PetscCall(DMGetLocalSection(dm, §ion)); 322 PetscCall(PetscSectionGetClosurePermutation(section, (PetscObject)dm, dim, num_comp * P, &permutation)); 323 PetscCall(ISGetIndices(permutation, &permutation_indices)); 324 } 325 326 // -- Copy interp, grad matrices 327 PetscCall(PetscCalloc(P * Q * sizeof(CeedScalar), &interp)); 328 PetscCall(PetscCalloc(P * Q * dim * sizeof(CeedScalar), &grad)); 329 const CeedInt c = 0; 330 for (CeedInt q = 0; q < Q; q++) { 331 for (CeedInt p_ceed = 0; p_ceed < P; p_ceed++) { 332 CeedInt p_petsc = is_simplex ? (p_ceed * num_comp) : permutation_indices[p_ceed * num_comp]; 333 334 interp[q * P + p_ceed] = basis_tabulation->T[0][((face * Q + q) * P * num_comp + p_petsc) * num_comp + c]; 335 for (CeedInt d = 0; d < dim; d++) { 336 grad[(d * Q + q) * P + p_ceed] = basis_tabulation->T[1][(((face * Q + q) * P * num_comp + p_petsc) * num_comp + c) * dim + d]; 337 } 338 } 339 } 340 341 // -- Cleanup 342 if (permutation) PetscCall(ISRestoreIndices(permutation, &permutation_indices)); 343 PetscCall(ISDestroy(&permutation)); 344 } 345 346 // Finally, create libCEED basis 347 CeedBasisCreateH1(ceed, elem_topo, num_comp, P, Q, interp, grad, q_points, q_weights, basis); 348 PetscCall(PetscFree(q_points)); 349 PetscCall(PetscFree(interp)); 350 PetscCall(PetscFree(grad)); 351 352 PetscFunctionReturn(PETSC_SUCCESS); 353 } 354 355 // ----------------------------------------------------------------------------- 356 // Get CEED Basis from DMPlex 357 // ----------------------------------------------------------------------------- 358 PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, CeedBasis *basis) { 359 PetscDS ds; 360 PetscFE fe; 361 PetscQuadrature quadrature; 362 PetscBool is_simplex = PETSC_TRUE; 363 PetscInt ds_field = -1; 364 365 PetscFunctionBeginUser; 366 367 // Get element information 368 PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL)); 369 PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field)); 370 PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe)); 371 PetscCall(PetscFEGetHeightSubspace(fe, height, &fe)); 372 PetscCall(PetscFEGetQuadrature(fe, &quadrature)); 373 374 // Check if simplex or tensor-product mesh 375 PetscCall(DMPlexIsSimplex(dm, &is_simplex)); 376 377 // Build libCEED basis 378 if (is_simplex) { 379 PetscTabulation basis_tabulation; 380 PetscInt num_derivatives = 1, face = 0; 381 382 PetscCall(PetscFEGetCellTabulation(fe, num_derivatives, &basis_tabulation)); 383 PetscCall(BasisCreateFromTabulation(ceed, dm, domain_label, label_value, height, face, fe, basis_tabulation, quadrature, basis)); 384 } else { 385 PetscDualSpace dual_space; 386 PetscInt num_dual_basis_vectors; 387 PetscInt dim, num_comp, P, Q; 388 389 PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 390 PetscCall(PetscFEGetNumComponents(fe, &num_comp)); 391 PetscCall(PetscFEGetDualSpace(fe, &dual_space)); 392 PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors)); 393 P = num_dual_basis_vectors / num_comp; 394 PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &Q, NULL, NULL)); 395 396 CeedInt P_1d = (CeedInt)round(pow(P, 1.0 / dim)); 397 CeedInt Q_1d = (CeedInt)round(pow(Q, 1.0 / dim)); 398 399 CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp, P_1d, Q_1d, CEED_GAUSS, basis); 400 } 401 402 PetscFunctionReturn(PETSC_SUCCESS); 403 } 404 405 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc) { 406 PetscFunctionBeginUser; 407 408 // ***************************************************************************** 409 // Set up CEED objects for the interior domain (volume) 410 // ***************************************************************************** 411 const PetscInt num_comp_q = 5; 412 const CeedInt dim = problem->dim, num_comp_x = problem->dim, q_data_size_vol = problem->q_data_size_vol, jac_data_size_vol = num_comp_q + 6 + 3; 413 CeedElemRestriction elem_restr_jd_i; 414 CeedVector jac_data; 415 CeedInt num_qpts; 416 417 // ----------------------------------------------------------------------------- 418 // CEED Bases 419 // ----------------------------------------------------------------------------- 420 DM dm_coord; 421 PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 422 423 PetscCall(CreateBasisFromPlex(ceed, dm, 0, 0, 0, 0, &ceed_data->basis_q)); 424 PetscCall(CreateBasisFromPlex(ceed, dm_coord, 0, 0, 0, 0, &ceed_data->basis_x)); 425 PetscCall(CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &ceed_data->basis_xc)); 426 CeedBasisGetNumQuadraturePoints(ceed_data->basis_q, &num_qpts); 427 428 // ----------------------------------------------------------------------------- 429 // CEED Restrictions 430 // ----------------------------------------------------------------------------- 431 // -- Create restriction 432 PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, 0, num_qpts, q_data_size_vol, &ceed_data->elem_restr_q, &ceed_data->elem_restr_x, 433 &ceed_data->elem_restr_qd_i)); 434 435 PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, 0, num_qpts, jac_data_size_vol, NULL, NULL, &elem_restr_jd_i)); 436 // -- Create E vectors 437 CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL); 438 CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed, NULL); 439 CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL); 440 441 // ----------------------------------------------------------------------------- 442 // CEED QFunctions 443 // ----------------------------------------------------------------------------- 444 // -- Create QFunction for quadrature data 445 CeedQFunctionCreateInterior(ceed, 1, problem->setup_vol.qfunction, problem->setup_vol.qfunction_loc, &ceed_data->qf_setup_vol); 446 if (problem->setup_vol.qfunction_context) { 447 CeedQFunctionSetContext(ceed_data->qf_setup_vol, problem->setup_vol.qfunction_context); 448 } 449 CeedQFunctionAddInput(ceed_data->qf_setup_vol, "dx", num_comp_x * dim, CEED_EVAL_GRAD); 450 CeedQFunctionAddInput(ceed_data->qf_setup_vol, "weight", 1, CEED_EVAL_WEIGHT); 451 CeedQFunctionAddOutput(ceed_data->qf_setup_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE); 452 453 // -- Create QFunction for ICs 454 CeedQFunctionCreateInterior(ceed, 1, problem->ics.qfunction, problem->ics.qfunction_loc, &ceed_data->qf_ics); 455 CeedQFunctionSetContext(ceed_data->qf_ics, problem->ics.qfunction_context); 456 CeedQFunctionAddInput(ceed_data->qf_ics, "x", num_comp_x, CEED_EVAL_INTERP); 457 CeedQFunctionAddInput(ceed_data->qf_ics, "dx", num_comp_x * dim, CEED_EVAL_GRAD); 458 CeedQFunctionAddOutput(ceed_data->qf_ics, "q0", num_comp_q, CEED_EVAL_NONE); 459 460 // -- Create QFunction for RHS 461 if (problem->apply_vol_rhs.qfunction) { 462 CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qfunction, problem->apply_vol_rhs.qfunction_loc, &ceed_data->qf_rhs_vol); 463 CeedQFunctionSetContext(ceed_data->qf_rhs_vol, problem->apply_vol_rhs.qfunction_context); 464 CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP); 465 CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD); 466 CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE); 467 CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP); 468 CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP); 469 CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD); 470 } 471 472 // -- Create QFunction for IFunction 473 if (problem->apply_vol_ifunction.qfunction) { 474 CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qfunction, problem->apply_vol_ifunction.qfunction_loc, 475 &ceed_data->qf_ifunction_vol); 476 CeedQFunctionSetContext(ceed_data->qf_ifunction_vol, problem->apply_vol_ifunction.qfunction_context); 477 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP); 478 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD); 479 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP); 480 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE); 481 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP); 482 CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP); 483 CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD); 484 CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE); 485 } 486 487 CeedQFunction qf_ijacobian_vol = NULL; 488 if (problem->apply_vol_ijacobian.qfunction) { 489 CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qfunction, problem->apply_vol_ijacobian.qfunction_loc, &qf_ijacobian_vol); 490 CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfunction_context); 491 CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP); 492 CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD); 493 CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE); 494 CeedQFunctionAddInput(qf_ijacobian_vol, "x", num_comp_x, CEED_EVAL_INTERP); 495 CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE); 496 CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP); 497 CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD); 498 } 499 500 // --------------------------------------------------------------------------- 501 // Element coordinates 502 // --------------------------------------------------------------------------- 503 // -- Create CEED vector 504 CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord, NULL); 505 506 // -- Copy PETSc vector in CEED vector 507 Vec X_loc; 508 { 509 DM cdm; 510 PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 511 if (cdm) { 512 PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc)); 513 } else { 514 PetscCall(DMGetCoordinatesLocal(dm, &X_loc)); 515 } 516 } 517 PetscCall(VecScale(X_loc, problem->dm_scale)); 518 PetscCall(VecCopyP2C(X_loc, ceed_data->x_coord)); 519 520 // ----------------------------------------------------------------------------- 521 // CEED vectors 522 // ----------------------------------------------------------------------------- 523 // -- Create CEED vector for geometric data 524 CeedElemRestrictionCreateVector(ceed_data->elem_restr_qd_i, &ceed_data->q_data, NULL); 525 CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL); 526 527 // ----------------------------------------------------------------------------- 528 // CEED Operators 529 // ----------------------------------------------------------------------------- 530 // -- Create CEED operator for quadrature data 531 CeedOperatorCreate(ceed, ceed_data->qf_setup_vol, NULL, NULL, &ceed_data->op_setup_vol); 532 CeedOperatorSetField(ceed_data->op_setup_vol, "dx", ceed_data->elem_restr_x, ceed_data->basis_x, CEED_VECTOR_ACTIVE); 533 CeedOperatorSetField(ceed_data->op_setup_vol, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x, CEED_VECTOR_NONE); 534 CeedOperatorSetField(ceed_data->op_setup_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 535 536 // -- Create CEED operator for ICs 537 CeedOperator op_ics; 538 CeedOperatorCreate(ceed, ceed_data->qf_ics, NULL, NULL, &op_ics); 539 CeedOperatorSetField(op_ics, "x", ceed_data->elem_restr_x, ceed_data->basis_xc, CEED_VECTOR_ACTIVE); 540 CeedOperatorSetField(op_ics, "dx", ceed_data->elem_restr_x, ceed_data->basis_xc, CEED_VECTOR_ACTIVE); 541 CeedOperatorSetField(op_ics, "q0", ceed_data->elem_restr_q, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 542 CeedOperatorGetContextFieldLabel(op_ics, "evaluation time", &user->phys->ics_time_label); 543 PetscCall(OperatorApplyContextCreate(NULL, dm, user->ceed, op_ics, ceed_data->x_coord, NULL, NULL, user->Q_loc, &ceed_data->op_ics_ctx)); 544 CeedOperatorDestroy(&op_ics); 545 546 // Create CEED operator for RHS 547 if (ceed_data->qf_rhs_vol) { 548 CeedOperator op; 549 CeedOperatorCreate(ceed, ceed_data->qf_rhs_vol, NULL, NULL, &op); 550 CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 551 CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 552 CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data); 553 CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord); 554 CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 555 CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 556 user->op_rhs_vol = op; 557 } 558 559 // -- CEED operator for IFunction 560 if (ceed_data->qf_ifunction_vol) { 561 CeedOperator op; 562 CeedOperatorCreate(ceed, ceed_data->qf_ifunction_vol, NULL, NULL, &op); 563 CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 564 CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 565 CeedOperatorSetField(op, "q dot", ceed_data->elem_restr_q, ceed_data->basis_q, user->q_dot_ceed); 566 CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data); 567 CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord); 568 CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 569 CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 570 CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_COLLOCATED, jac_data); 571 572 user->op_ifunction_vol = op; 573 } 574 575 CeedOperator op_ijacobian_vol = NULL; 576 if (qf_ijacobian_vol) { 577 CeedOperator op; 578 CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op); 579 CeedOperatorSetField(op, "dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 580 CeedOperatorSetField(op, "Grad_dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 581 CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data); 582 CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord); 583 CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_COLLOCATED, jac_data); 584 CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 585 CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 586 op_ijacobian_vol = op; 587 CeedQFunctionDestroy(&qf_ijacobian_vol); 588 } 589 590 // ***************************************************************************** 591 // Set up CEED objects for the exterior domain (surface) 592 // ***************************************************************************** 593 CeedInt height = 1, dim_sur = dim - height, P_sur = app_ctx->degree + 1, Q_sur = P_sur + app_ctx->q_extra; 594 const CeedInt q_data_size_sur = problem->q_data_size_sur, jac_data_size_sur = problem->jac_data_size_sur; 595 596 // ----------------------------------------------------------------------------- 597 // CEED Bases 598 // ----------------------------------------------------------------------------- 599 600 DMLabel label = 0; 601 PetscInt face_id = 0; 602 PetscInt field = 0; // Still want the normal, default field 603 PetscCall(CreateBasisFromPlex(ceed, dm, label, face_id, height, field, &ceed_data->basis_q_sur)); 604 PetscCall(CreateBasisFromPlex(ceed, dm_coord, label, face_id, height, field, &ceed_data->basis_x_sur)); 605 PetscCall(CeedBasisCreateProjection(ceed_data->basis_x_sur, ceed_data->basis_q_sur, &ceed_data->basis_xc_sur)); 606 607 // ----------------------------------------------------------------------------- 608 // CEED QFunctions 609 // ----------------------------------------------------------------------------- 610 // -- Create QFunction for quadrature data 611 CeedQFunctionCreateInterior(ceed, 1, problem->setup_sur.qfunction, problem->setup_sur.qfunction_loc, &ceed_data->qf_setup_sur); 612 if (problem->setup_sur.qfunction_context) { 613 CeedQFunctionSetContext(ceed_data->qf_setup_sur, problem->setup_sur.qfunction_context); 614 } 615 CeedQFunctionAddInput(ceed_data->qf_setup_sur, "dx", num_comp_x * dim_sur, CEED_EVAL_GRAD); 616 CeedQFunctionAddInput(ceed_data->qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT); 617 CeedQFunctionAddOutput(ceed_data->qf_setup_sur, "surface qdata", q_data_size_sur, CEED_EVAL_NONE); 618 619 PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_inflow, 620 problem->apply_inflow_jacobian, &ceed_data->qf_apply_inflow, &ceed_data->qf_apply_inflow_jacobian)); 621 PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_outflow, 622 problem->apply_outflow_jacobian, &ceed_data->qf_apply_outflow, &ceed_data->qf_apply_outflow_jacobian)); 623 PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_freestream, 624 problem->apply_freestream_jacobian, &ceed_data->qf_apply_freestream, &ceed_data->qf_apply_freestream_jacobian)); 625 626 // ***************************************************************************** 627 // CEED Operator Apply 628 // ***************************************************************************** 629 // -- Apply CEED Operator for the geometric data 630 CeedOperatorApply(ceed_data->op_setup_vol, ceed_data->x_coord, ceed_data->q_data, CEED_REQUEST_IMMEDIATE); 631 632 // -- Create and apply CEED Composite Operator for the entire domain 633 if (!user->phys->implicit) { // RHS 634 CeedOperator op_rhs; 635 PetscCall(CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys, user->op_rhs_vol, NULL, height, P_sur, Q_sur, q_data_size_sur, 0, &op_rhs, 636 NULL)); 637 PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_rhs, user->q_ceed, user->g_ceed, user->Q_loc, NULL, &user->op_rhs_ctx)); 638 CeedOperatorDestroy(&op_rhs); 639 } else { // IFunction 640 PetscCall(CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys, user->op_ifunction_vol, op_ijacobian_vol, height, P_sur, Q_sur, 641 q_data_size_sur, jac_data_size_sur, &user->op_ifunction, op_ijacobian_vol ? &user->op_ijacobian : NULL)); 642 if (user->op_ijacobian) { 643 CeedOperatorGetContextFieldLabel(user->op_ijacobian, "ijacobian time shift", &user->phys->ijacobian_time_shift_label); 644 } 645 if (problem->use_strong_bc_ceed) PetscCall(SetupStrongBC_Ceed(ceed, ceed_data, dm, user, problem, bc, q_data_size_sur)); 646 if (app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) PetscCall(SGS_DD_ModelSetup(ceed, user, ceed_data, problem)); 647 } 648 649 if (app_ctx->turb_spanstats_enable) PetscCall(TurbulenceStatisticsSetup(ceed, user, ceed_data, problem)); 650 if (app_ctx->diff_filter_monitor) PetscCall(DifferentialFilterSetup(ceed, user, ceed_data, problem)); 651 652 CeedElemRestrictionDestroy(&elem_restr_jd_i); 653 CeedOperatorDestroy(&op_ijacobian_vol); 654 CeedVectorDestroy(&jac_data); 655 PetscFunctionReturn(PETSC_SUCCESS); 656 } 657