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) { return i >= 0 ? i : -(i + 1); } 15 16 // Utility function to create local CEED restriction 17 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt value, CeedElemRestriction *elem_restr) { 18 PetscInt num_elem, elem_size, num_dof, num_comp, *elem_restr_offsets; 19 20 PetscFunctionBeginUser; 21 PetscCall(DMPlexGetLocalOffsets(dm, domain_label, value, height, 0, &num_elem, &elem_size, &num_comp, &num_dof, &elem_restr_offsets)); 22 23 CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp, 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES, elem_restr_offsets, elem_restr); 24 PetscCall(PetscFree(elem_restr_offsets)); 25 26 PetscFunctionReturn(0); 27 } 28 29 // Utility function to get Ceed Restriction for each domain 30 PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, PetscInt value, CeedInt Q, CeedInt q_data_size, 31 CeedElemRestriction *elem_restr_q, CeedElemRestriction *elem_restr_x, CeedElemRestriction *elem_restr_qd_i) { 32 DM dm_coord; 33 CeedInt dim, loc_num_elem; 34 CeedInt Q_dim; 35 CeedElemRestriction elem_restr_tmp; 36 PetscFunctionBeginUser; 37 38 PetscCall(DMGetDimension(dm, &dim)); 39 dim -= height; 40 Q_dim = CeedIntPow(Q, dim); 41 PetscCall(CreateRestrictionFromPlex(ceed, dm, height, domain_label, value, &elem_restr_tmp)); 42 if (elem_restr_q) *elem_restr_q = elem_restr_tmp; 43 if (elem_restr_x) { 44 PetscCall(DMGetCellCoordinateDM(dm, &dm_coord)); 45 if (!dm_coord) { 46 PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 47 } 48 PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL)); 49 PetscCall(CreateRestrictionFromPlex(ceed, dm_coord, height, domain_label, value, elem_restr_x)); 50 } 51 if (elem_restr_qd_i) { 52 CeedElemRestrictionGetNumElements(elem_restr_tmp, &loc_num_elem); 53 CeedElemRestrictionCreateStrided(ceed, loc_num_elem, Q_dim, q_data_size, q_data_size * loc_num_elem * Q_dim, CEED_STRIDES_BACKEND, 54 elem_restr_qd_i); 55 } 56 if (!elem_restr_q) CeedElemRestrictionDestroy(&elem_restr_tmp); 57 PetscFunctionReturn(0); 58 } 59 60 PetscErrorCode AddBCSubOperator(Ceed ceed, DM dm, CeedData ceed_data, DMLabel domain_label, PetscInt value, CeedInt height, CeedInt Q_sur, 61 CeedInt q_data_size_sur, CeedInt jac_data_size_sur, CeedQFunction qf_apply_bc, CeedQFunction qf_apply_bc_jacobian, 62 CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) { 63 CeedVector q_data_sur, jac_data_sur; 64 CeedOperator op_setup_sur, op_apply_bc, op_apply_bc_jacobian = NULL; 65 CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur, elem_restr_qd_i_sur, elem_restr_jd_i_sur; 66 CeedInt num_qpts_sur; 67 PetscFunctionBeginUser; 68 69 // --- Get number of quadrature points for the boundaries 70 CeedBasisGetNumQuadraturePoints(ceed_data->basis_q_sur, &num_qpts_sur); 71 72 // ---- CEED Restriction 73 PetscCall(GetRestrictionForDomain(ceed, dm, height, domain_label, value, Q_sur, q_data_size_sur, &elem_restr_q_sur, &elem_restr_x_sur, 74 &elem_restr_qd_i_sur)); 75 if (jac_data_size_sur > 0) { 76 // State-dependent data will be passed from residual to Jacobian. This will be collocated. 77 PetscCall(GetRestrictionForDomain(ceed, dm, height, domain_label, value, Q_sur, jac_data_size_sur, NULL, NULL, &elem_restr_jd_i_sur)); 78 CeedElemRestrictionCreateVector(elem_restr_jd_i_sur, &jac_data_sur, NULL); 79 } else { 80 elem_restr_jd_i_sur = NULL; 81 jac_data_sur = NULL; 82 } 83 84 // ---- CEED Vector 85 PetscInt loc_num_elem_sur; 86 CeedElemRestrictionGetNumElements(elem_restr_q_sur, &loc_num_elem_sur); 87 CeedVectorCreate(ceed, q_data_size_sur * loc_num_elem_sur * num_qpts_sur, &q_data_sur); 88 89 // ---- CEED Operator 90 // ----- CEED Operator for Setup (geometric factors) 91 CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur); 92 CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x_sur, ceed_data->basis_x_sur, CEED_VECTOR_ACTIVE); 93 CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x_sur, CEED_VECTOR_NONE); 94 CeedOperatorSetField(op_setup_sur, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 95 96 // ----- CEED Operator for Physics 97 CeedOperatorCreate(ceed, qf_apply_bc, NULL, NULL, &op_apply_bc); 98 CeedOperatorSetField(op_apply_bc, "q", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE); 99 CeedOperatorSetField(op_apply_bc, "Grad_q", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE); 100 CeedOperatorSetField(op_apply_bc, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, q_data_sur); 101 CeedOperatorSetField(op_apply_bc, "x", elem_restr_x_sur, ceed_data->basis_x_sur, ceed_data->x_coord); 102 CeedOperatorSetField(op_apply_bc, "v", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE); 103 if (elem_restr_jd_i_sur) CeedOperatorSetField(op_apply_bc, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_COLLOCATED, jac_data_sur); 104 105 if (qf_apply_bc_jacobian) { 106 CeedOperatorCreate(ceed, qf_apply_bc_jacobian, NULL, NULL, &op_apply_bc_jacobian); 107 CeedOperatorSetField(op_apply_bc_jacobian, "dq", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE); 108 CeedOperatorSetField(op_apply_bc_jacobian, "Grad_dq", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE); 109 CeedOperatorSetField(op_apply_bc_jacobian, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, q_data_sur); 110 CeedOperatorSetField(op_apply_bc_jacobian, "x", elem_restr_x_sur, ceed_data->basis_x_sur, ceed_data->x_coord); 111 CeedOperatorSetField(op_apply_bc_jacobian, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_COLLOCATED, jac_data_sur); 112 CeedOperatorSetField(op_apply_bc_jacobian, "v", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE); 113 } 114 115 // ----- Apply CEED operator for Setup 116 CeedOperatorApply(op_setup_sur, ceed_data->x_coord, q_data_sur, CEED_REQUEST_IMMEDIATE); 117 118 // ----- Apply Sub-Operator for Physics 119 CeedCompositeOperatorAddSub(*op_apply, op_apply_bc); 120 if (op_apply_bc_jacobian) CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_bc_jacobian); 121 122 // ----- Cleanup 123 CeedVectorDestroy(&q_data_sur); 124 CeedVectorDestroy(&jac_data_sur); 125 CeedElemRestrictionDestroy(&elem_restr_q_sur); 126 CeedElemRestrictionDestroy(&elem_restr_x_sur); 127 CeedElemRestrictionDestroy(&elem_restr_qd_i_sur); 128 CeedElemRestrictionDestroy(&elem_restr_jd_i_sur); 129 CeedOperatorDestroy(&op_setup_sur); 130 CeedOperatorDestroy(&op_apply_bc); 131 CeedOperatorDestroy(&op_apply_bc_jacobian); 132 133 PetscFunctionReturn(0); 134 } 135 136 // Utility function to create CEED Composite Operator for the entire domain 137 PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc, CeedData ceed_data, Physics phys, CeedOperator op_apply_vol, 138 CeedOperator op_apply_ijacobian_vol, CeedInt height, CeedInt P_sur, CeedInt Q_sur, CeedInt q_data_size_sur, 139 CeedInt jac_data_size_sur, CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) { 140 DMLabel domain_label; 141 PetscFunctionBeginUser; 142 143 // Create Composite Operaters 144 CeedCompositeOperatorCreate(ceed, op_apply); 145 if (op_apply_ijacobian) CeedCompositeOperatorCreate(ceed, op_apply_ijacobian); 146 147 // --Apply Sub-Operator for the volume 148 CeedCompositeOperatorAddSub(*op_apply, op_apply_vol); 149 if (op_apply_ijacobian) CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_ijacobian_vol); 150 151 // -- Create Sub-Operator for in/outflow BCs 152 if (phys->has_neumann || 1) { 153 // --- Setup 154 PetscCall(DMGetLabel(dm, "Face Sets", &domain_label)); 155 156 // --- Create Sub-Operator for inflow boundaries 157 for (CeedInt i = 0; i < bc->num_inflow; i++) { 158 PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->inflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, 159 ceed_data->qf_apply_inflow, ceed_data->qf_apply_inflow_jacobian, op_apply, op_apply_ijacobian)); 160 } 161 // --- Create Sub-Operator for outflow boundaries 162 for (CeedInt i = 0; i < bc->num_outflow; i++) { 163 PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->outflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, 164 ceed_data->qf_apply_outflow, ceed_data->qf_apply_outflow_jacobian, op_apply, op_apply_ijacobian)); 165 } 166 // --- Create Sub-Operator for freestream boundaries 167 for (CeedInt i = 0; i < bc->num_freestream; i++) { 168 PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->freestreams[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, 169 ceed_data->qf_apply_freestream, ceed_data->qf_apply_freestream_jacobian, op_apply, op_apply_ijacobian)); 170 } 171 } 172 173 // ----- Get Context Labels for Operator 174 CeedOperatorGetContextFieldLabel(*op_apply, "solution time", &phys->solution_time_label); 175 CeedOperatorGetContextFieldLabel(*op_apply, "timestep size", &phys->timestep_size_label); 176 177 PetscFunctionReturn(0); 178 } 179 180 PetscErrorCode SetupBCQFunctions(Ceed ceed, PetscInt dim_sur, PetscInt num_comp_x, PetscInt num_comp_q, PetscInt q_data_size_sur, 181 PetscInt jac_data_size_sur, ProblemQFunctionSpec apply_bc, ProblemQFunctionSpec apply_bc_jacobian, 182 CeedQFunction *qf_apply_bc, CeedQFunction *qf_apply_bc_jacobian) { 183 PetscFunctionBeginUser; 184 185 if (apply_bc.qfunction) { 186 CeedQFunctionCreateInterior(ceed, 1, apply_bc.qfunction, apply_bc.qfunction_loc, qf_apply_bc); 187 CeedQFunctionSetContext(*qf_apply_bc, apply_bc.qfunction_context); 188 CeedQFunctionAddInput(*qf_apply_bc, "q", num_comp_q, CEED_EVAL_INTERP); 189 CeedQFunctionAddInput(*qf_apply_bc, "Grad_q", num_comp_q * dim_sur, CEED_EVAL_GRAD); 190 CeedQFunctionAddInput(*qf_apply_bc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE); 191 CeedQFunctionAddInput(*qf_apply_bc, "x", num_comp_x, CEED_EVAL_INTERP); 192 CeedQFunctionAddOutput(*qf_apply_bc, "v", num_comp_q, CEED_EVAL_INTERP); 193 if (jac_data_size_sur) CeedQFunctionAddOutput(*qf_apply_bc, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE); 194 } 195 if (apply_bc_jacobian.qfunction) { 196 CeedQFunctionCreateInterior(ceed, 1, apply_bc_jacobian.qfunction, apply_bc_jacobian.qfunction_loc, qf_apply_bc_jacobian); 197 CeedQFunctionSetContext(*qf_apply_bc_jacobian, apply_bc_jacobian.qfunction_context); 198 CeedQFunctionAddInput(*qf_apply_bc_jacobian, "dq", num_comp_q, CEED_EVAL_INTERP); 199 CeedQFunctionAddInput(*qf_apply_bc_jacobian, "Grad_dq", num_comp_q * dim_sur, CEED_EVAL_GRAD); 200 CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface qdata", q_data_size_sur, CEED_EVAL_NONE); 201 CeedQFunctionAddInput(*qf_apply_bc_jacobian, "x", num_comp_x, CEED_EVAL_INTERP); 202 CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE); 203 CeedQFunctionAddOutput(*qf_apply_bc_jacobian, "v", num_comp_q, CEED_EVAL_INTERP); 204 } 205 PetscFunctionReturn(0); 206 } 207 208 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc) { 209 PetscFunctionBeginUser; 210 211 // ***************************************************************************** 212 // Set up CEED objects for the interior domain (volume) 213 // ***************************************************************************** 214 const PetscInt num_comp_q = 5; 215 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, 216 P = app_ctx->degree + 1, Q = P + app_ctx->q_extra; 217 CeedElemRestriction elem_restr_jd_i; 218 CeedVector jac_data; 219 220 // ----------------------------------------------------------------------------- 221 // CEED Bases 222 // ----------------------------------------------------------------------------- 223 CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_q, P, Q, CEED_GAUSS, &ceed_data->basis_q); 224 CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, CEED_GAUSS, &ceed_data->basis_x); 225 CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, P, CEED_GAUSS_LOBATTO, &ceed_data->basis_xc); 226 227 // ----------------------------------------------------------------------------- 228 // CEED Restrictions 229 // ----------------------------------------------------------------------------- 230 // -- Create restriction 231 PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, Q, q_data_size_vol, &ceed_data->elem_restr_q, &ceed_data->elem_restr_x, 232 &ceed_data->elem_restr_qd_i)); 233 234 PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, Q, jac_data_size_vol, NULL, NULL, &elem_restr_jd_i)); 235 // -- Create E vectors 236 CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL); 237 CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed, NULL); 238 CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL); 239 240 // ----------------------------------------------------------------------------- 241 // CEED QFunctions 242 // ----------------------------------------------------------------------------- 243 // -- Create QFunction for quadrature data 244 CeedQFunctionCreateInterior(ceed, 1, problem->setup_vol.qfunction, problem->setup_vol.qfunction_loc, &ceed_data->qf_setup_vol); 245 if (problem->setup_vol.qfunction_context) { 246 CeedQFunctionSetContext(ceed_data->qf_setup_vol, problem->setup_vol.qfunction_context); 247 } 248 CeedQFunctionAddInput(ceed_data->qf_setup_vol, "dx", num_comp_x * dim, CEED_EVAL_GRAD); 249 CeedQFunctionAddInput(ceed_data->qf_setup_vol, "weight", 1, CEED_EVAL_WEIGHT); 250 CeedQFunctionAddOutput(ceed_data->qf_setup_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE); 251 252 // -- Create QFunction for ICs 253 CeedQFunctionCreateInterior(ceed, 1, problem->ics.qfunction, problem->ics.qfunction_loc, &ceed_data->qf_ics); 254 CeedQFunctionSetContext(ceed_data->qf_ics, problem->ics.qfunction_context); 255 CeedQFunctionAddInput(ceed_data->qf_ics, "x", num_comp_x, CEED_EVAL_INTERP); 256 CeedQFunctionAddInput(ceed_data->qf_ics, "qdata", q_data_size_vol, CEED_EVAL_NONE); 257 CeedQFunctionAddOutput(ceed_data->qf_ics, "q0", num_comp_q, CEED_EVAL_NONE); 258 259 // -- Create QFunction for RHS 260 if (problem->apply_vol_rhs.qfunction) { 261 CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qfunction, problem->apply_vol_rhs.qfunction_loc, &ceed_data->qf_rhs_vol); 262 CeedQFunctionSetContext(ceed_data->qf_rhs_vol, problem->apply_vol_rhs.qfunction_context); 263 CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP); 264 CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD); 265 CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE); 266 CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP); 267 CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP); 268 CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD); 269 } 270 271 // -- Create QFunction for IFunction 272 if (problem->apply_vol_ifunction.qfunction) { 273 CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qfunction, problem->apply_vol_ifunction.qfunction_loc, 274 &ceed_data->qf_ifunction_vol); 275 CeedQFunctionSetContext(ceed_data->qf_ifunction_vol, problem->apply_vol_ifunction.qfunction_context); 276 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP); 277 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD); 278 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP); 279 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE); 280 CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP); 281 CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP); 282 CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD); 283 CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE); 284 } 285 286 CeedQFunction qf_ijacobian_vol = NULL; 287 if (problem->apply_vol_ijacobian.qfunction) { 288 CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qfunction, problem->apply_vol_ijacobian.qfunction_loc, &qf_ijacobian_vol); 289 CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfunction_context); 290 CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP); 291 CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD); 292 CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE); 293 CeedQFunctionAddInput(qf_ijacobian_vol, "x", num_comp_x, CEED_EVAL_INTERP); 294 CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE); 295 CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP); 296 CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD); 297 } 298 299 // --------------------------------------------------------------------------- 300 // Element coordinates 301 // --------------------------------------------------------------------------- 302 // -- Create CEED vector 303 CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord, NULL); 304 305 // -- Copy PETSc vector in CEED vector 306 Vec X_loc; 307 const PetscScalar *X_loc_array; 308 { 309 DM cdm; 310 PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 311 if (cdm) { 312 PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc)); 313 } else { 314 PetscCall(DMGetCoordinatesLocal(dm, &X_loc)); 315 } 316 } 317 PetscCall(VecScale(X_loc, problem->dm_scale)); 318 PetscCall(VecGetArrayRead(X_loc, &X_loc_array)); 319 CeedVectorSetArray(ceed_data->x_coord, CEED_MEM_HOST, CEED_COPY_VALUES, (PetscScalar *)X_loc_array); 320 PetscCall(VecRestoreArrayRead(X_loc, &X_loc_array)); 321 322 // ----------------------------------------------------------------------------- 323 // CEED vectors 324 // ----------------------------------------------------------------------------- 325 // -- Create CEED vector for geometric data 326 CeedElemRestrictionCreateVector(ceed_data->elem_restr_qd_i, &ceed_data->q_data, NULL); 327 CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL); 328 329 // ----------------------------------------------------------------------------- 330 // CEED Operators 331 // ----------------------------------------------------------------------------- 332 // -- Create CEED operator for quadrature data 333 CeedOperatorCreate(ceed, ceed_data->qf_setup_vol, NULL, NULL, &ceed_data->op_setup_vol); 334 CeedOperatorSetField(ceed_data->op_setup_vol, "dx", ceed_data->elem_restr_x, ceed_data->basis_x, CEED_VECTOR_ACTIVE); 335 CeedOperatorSetField(ceed_data->op_setup_vol, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x, CEED_VECTOR_NONE); 336 CeedOperatorSetField(ceed_data->op_setup_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 337 338 // -- Create CEED operator for ICs 339 CeedOperatorCreate(ceed, ceed_data->qf_ics, NULL, NULL, &ceed_data->op_ics); 340 CeedOperatorSetField(ceed_data->op_ics, "x", ceed_data->elem_restr_x, ceed_data->basis_xc, CEED_VECTOR_ACTIVE); 341 CeedOperatorSetField(ceed_data->op_ics, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data); 342 CeedOperatorSetField(ceed_data->op_ics, "q0", ceed_data->elem_restr_q, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 343 CeedOperatorGetContextFieldLabel(ceed_data->op_ics, "evaluation time", &user->phys->ics_time_label); 344 345 // Create CEED operator for RHS 346 if (ceed_data->qf_rhs_vol) { 347 CeedOperator op; 348 CeedOperatorCreate(ceed, ceed_data->qf_rhs_vol, NULL, NULL, &op); 349 CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 350 CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 351 CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data); 352 CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord); 353 CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 354 CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 355 user->op_rhs_vol = op; 356 } 357 358 // -- CEED operator for IFunction 359 if (ceed_data->qf_ifunction_vol) { 360 CeedOperator op; 361 CeedOperatorCreate(ceed, ceed_data->qf_ifunction_vol, NULL, NULL, &op); 362 CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 363 CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 364 CeedOperatorSetField(op, "q dot", ceed_data->elem_restr_q, ceed_data->basis_q, user->q_dot_ceed); 365 CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data); 366 CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord); 367 CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 368 CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 369 CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_COLLOCATED, jac_data); 370 371 user->op_ifunction_vol = op; 372 } 373 374 CeedOperator op_ijacobian_vol = NULL; 375 if (qf_ijacobian_vol) { 376 CeedOperator op; 377 CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op); 378 CeedOperatorSetField(op, "dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 379 CeedOperatorSetField(op, "Grad_dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 380 CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data); 381 CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord); 382 CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_COLLOCATED, jac_data); 383 CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 384 CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 385 op_ijacobian_vol = op; 386 CeedQFunctionDestroy(&qf_ijacobian_vol); 387 } 388 389 // ***************************************************************************** 390 // Set up CEED objects for the exterior domain (surface) 391 // ***************************************************************************** 392 CeedInt height = 1, dim_sur = dim - height, P_sur = app_ctx->degree + 1, Q_sur = P_sur + app_ctx->q_extra; 393 const CeedInt q_data_size_sur = problem->q_data_size_sur, jac_data_size_sur = problem->jac_data_size_sur; 394 395 // ----------------------------------------------------------------------------- 396 // CEED Bases 397 // ----------------------------------------------------------------------------- 398 CeedBasisCreateTensorH1Lagrange(ceed, dim_sur, num_comp_q, P_sur, Q_sur, CEED_GAUSS, &ceed_data->basis_q_sur); 399 CeedBasisCreateTensorH1Lagrange(ceed, dim_sur, num_comp_x, 2, Q_sur, CEED_GAUSS, &ceed_data->basis_x_sur); 400 CeedBasisCreateTensorH1Lagrange(ceed, dim_sur, num_comp_x, 2, P_sur, CEED_GAUSS_LOBATTO, &ceed_data->basis_xc_sur); 401 402 // ----------------------------------------------------------------------------- 403 // CEED QFunctions 404 // ----------------------------------------------------------------------------- 405 // -- Create QFunction for quadrature data 406 CeedQFunctionCreateInterior(ceed, 1, problem->setup_sur.qfunction, problem->setup_sur.qfunction_loc, &ceed_data->qf_setup_sur); 407 if (problem->setup_sur.qfunction_context) { 408 CeedQFunctionSetContext(ceed_data->qf_setup_sur, problem->setup_sur.qfunction_context); 409 } 410 CeedQFunctionAddInput(ceed_data->qf_setup_sur, "dx", num_comp_x * dim_sur, CEED_EVAL_GRAD); 411 CeedQFunctionAddInput(ceed_data->qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT); 412 CeedQFunctionAddOutput(ceed_data->qf_setup_sur, "surface qdata", q_data_size_sur, CEED_EVAL_NONE); 413 414 PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_inflow, 415 problem->apply_inflow_jacobian, &ceed_data->qf_apply_inflow, &ceed_data->qf_apply_inflow_jacobian)); 416 PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_outflow, 417 problem->apply_outflow_jacobian, &ceed_data->qf_apply_outflow, &ceed_data->qf_apply_outflow_jacobian)); 418 PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_freestream, 419 problem->apply_freestream_jacobian, &ceed_data->qf_apply_freestream, &ceed_data->qf_apply_freestream_jacobian)); 420 421 // ***************************************************************************** 422 // CEED Operator Apply 423 // ***************************************************************************** 424 // -- Apply CEED Operator for the geometric data 425 CeedOperatorApply(ceed_data->op_setup_vol, ceed_data->x_coord, ceed_data->q_data, CEED_REQUEST_IMMEDIATE); 426 427 // -- Create and apply CEED Composite Operator for the entire domain 428 if (!user->phys->implicit) { // RHS 429 PetscCall(CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys, user->op_rhs_vol, NULL, height, P_sur, Q_sur, q_data_size_sur, 0, 430 &user->op_rhs, NULL)); 431 } else { // IFunction 432 PetscCall(CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys, user->op_ifunction_vol, op_ijacobian_vol, height, P_sur, Q_sur, 433 q_data_size_sur, jac_data_size_sur, &user->op_ifunction, op_ijacobian_vol ? &user->op_ijacobian : NULL)); 434 if (user->op_ijacobian) { 435 CeedOperatorGetContextFieldLabel(user->op_ijacobian, "ijacobian time shift", &user->phys->ijacobian_time_shift_label); 436 } 437 if (problem->use_dirichlet_ceed) { 438 PetscCall(SetupStrongBC_Ceed(ceed, ceed_data, dm, user, app_ctx, problem, bc, Q_sur, q_data_size_sur)); 439 } 440 } 441 CeedElemRestrictionDestroy(&elem_restr_jd_i); 442 CeedOperatorDestroy(&op_ijacobian_vol); 443 CeedVectorDestroy(&jac_data); 444 PetscFunctionReturn(0); 445 } 446