1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3a515125bSLeila Ghaffari 4a515125bSLeila Ghaffari /// @file 5a515125bSLeila Ghaffari /// Utility functions for setting up ADVECTION 6a515125bSLeila Ghaffari 72b916ea7SJeremy L Thompson #include "../qfunctions/advection.h" 82b916ea7SJeremy L Thompson 9e419654dSJeremy L Thompson #include <ceed.h> 10e419654dSJeremy L Thompson #include <petscdm.h> 11e419654dSJeremy L Thompson 12149fb536SJames Wright #include <navierstokes.h> 13a515125bSLeila Ghaffari 14a78efa86SJames Wright // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping 15a78efa86SJames Wright // 16a78efa86SJames Wright // Only used for SUPG stabilization 170c373b74SJames Wright PetscErrorCode CreateKSPMassOperator_AdvectionStabilized(Honee honee, CeedOperator *op_mass) { 180c373b74SJames Wright Ceed ceed = honee->ceed; 19a78efa86SJames Wright CeedInt num_comp_q, q_data_size; 2087d3884fSJeremy L Thompson CeedQFunction qf_mass = NULL; 2101e19bfaSJames Wright CeedElemRestriction elem_restr_q, elem_restr_qd; 22a78efa86SJames Wright CeedBasis basis_q; 23a78efa86SJames Wright CeedVector q_data; 24236a8ba6SJames Wright CeedQFunctionContext qfctx = NULL; 25a78efa86SJames Wright PetscInt dim; 26a78efa86SJames Wright 27a78efa86SJames Wright PetscFunctionBeginUser; 280c373b74SJames Wright PetscCall(DMGetDimension(honee->dm, &dim)); 29a78efa86SJames Wright { // Get restriction and basis from the RHS function 30a78efa86SJames Wright CeedOperator *sub_ops; 3101e19bfaSJames Wright CeedOperatorField op_field; 32a78efa86SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 33a78efa86SJames Wright 340c373b74SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_rhs_ctx->op, &sub_ops)); 3501e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &op_field)); 3601e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_q, &basis_q, NULL)); 3701e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &op_field)); 3801e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_qd, NULL, &q_data)); 39a78efa86SJames Wright 40236a8ba6SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qfctx)); 41a78efa86SJames Wright } 42a78efa86SJames Wright 43a78efa86SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q)); 4401e19bfaSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd, &q_data_size)); 45a78efa86SJames Wright 46a78efa86SJames Wright switch (dim) { 47a78efa86SJames Wright case 2: 48a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection2D, MassFunction_Advection2D_loc, &qf_mass)); 49a78efa86SJames Wright break; 50a78efa86SJames Wright case 3: 51a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection, MassFunction_Advection_loc, &qf_mass)); 52a78efa86SJames Wright break; 53a78efa86SJames Wright } 54a78efa86SJames Wright 55236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qfctx)); 56a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0)); 57a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP)); 58a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP)); 59a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE)); 60a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP)); 61a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD)); 62a78efa86SJames Wright 63a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass)); 64a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 650c373b74SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, honee->q_ceed)); 6601e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 67a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 68a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 69a78efa86SJames Wright 70fff85bd3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 7101e19bfaSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 72fff85bd3SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 73fff85bd3SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 74236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&qfctx)); 75a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 76a78efa86SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 77a78efa86SJames Wright } 78a78efa86SJames Wright 79c2d90829SJames Wright /** 80c2d90829SJames Wright @brief Create RHS CeedOperator for direct projection of divergence of diffusive flux 81c2d90829SJames Wright 820c373b74SJames Wright @param[in] honee `Honee` context 83c2d90829SJames Wright @param[in] diff_flux_proj `DivDiffFluxProjectionData` object 84c2d90829SJames Wright @param[out] op_rhs Operator to calculate the RHS of the L^2 projection 85c2d90829SJames Wright **/ 86e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Direct_AdvDif(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) { 870c373b74SJames Wright Ceed ceed = honee->ceed; 88c2d90829SJames Wright NodalProjectionData projection = diff_flux_proj->projection; 89c2d90829SJames Wright CeedInt num_comp_q; 90c2d90829SJames Wright PetscInt dim, label_value = 0; 91c2d90829SJames Wright DMLabel domain_label = NULL; 92c2d90829SJames Wright CeedQFunctionContext advection_qfctx = NULL; 93c2d90829SJames Wright 94c2d90829SJames Wright PetscFunctionBeginUser; 95c2d90829SJames Wright // -- Get Pre-requisite things 96c2d90829SJames Wright PetscCall(DMGetDimension(projection->dm, &dim)); 97e3663b90SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q)); 98c2d90829SJames Wright 9940b78511SJames Wright { // Get advection-diffusion QF context 100c2d90829SJames Wright CeedOperator *sub_ops; 101c2d90829SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 102c2d90829SJames Wright 1030c373b74SJames Wright if (honee->op_ifunction) PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_ifunction, &sub_ops)); 1040c373b74SJames Wright else PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_rhs_ctx->op, &sub_ops)); 105c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &advection_qfctx)); 106c2d90829SJames Wright } 107c2d90829SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, op_rhs)); 108c2d90829SJames Wright { // Add the volume integral CeedOperator 109c2d90829SJames Wright CeedQFunction qf_rhs_volume = NULL; 110c2d90829SJames Wright CeedOperator op_rhs_volume; 111c2d90829SJames Wright CeedVector q_data; 112c2d90829SJames Wright CeedElemRestriction elem_restr_qd, elem_restr_diff_flux_volume = NULL; 113c2d90829SJames Wright CeedBasis basis_diff_flux = NULL; 114c2d90829SJames Wright CeedInt q_data_size; 115c2d90829SJames Wright 116c2d90829SJames Wright PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_diff_flux_volume, &basis_diff_flux, NULL, NULL)); 117e3663b90SJames Wright PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, &q_data, 118e3663b90SJames Wright &q_data_size)); 119c2d90829SJames Wright switch (dim) { 120c2d90829SJames Wright case 2: 121c2d90829SJames Wright PetscCallCeed( 122c2d90829SJames Wright ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_AdvDif_2D, DivDiffusiveFluxVolumeRHS_AdvDif_2D_loc, &qf_rhs_volume)); 123c2d90829SJames Wright break; 124c2d90829SJames Wright case 3: 125c2d90829SJames Wright PetscCallCeed( 126c2d90829SJames Wright ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_AdvDif_3D, DivDiffusiveFluxVolumeRHS_AdvDif_3D_loc, &qf_rhs_volume)); 127c2d90829SJames Wright break; 128c2d90829SJames Wright } 1290c373b74SJames Wright PetscCheck(qf_rhs_volume, honee->comm, PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim); 130c2d90829SJames Wright 131c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_volume, advection_qfctx)); 132c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 133c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "qdata", q_data_size, CEED_EVAL_NONE)); 134c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_volume, "diffusive flux RHS", projection->num_comp * dim, CEED_EVAL_GRAD)); 135c2d90829SJames Wright 136c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_volume, NULL, NULL, &op_rhs_volume)); 137e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 138c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 139c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "diffusive flux RHS", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE)); 140c2d90829SJames Wright 141c2d90829SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_volume)); 142c2d90829SJames Wright 143c2d90829SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 144c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 145c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_volume)); 146c2d90829SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 147c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_volume)); 148c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_volume)); 149c2d90829SJames Wright } 150c2d90829SJames Wright 151c2d90829SJames Wright { // Add the boundary integral CeedOperator 152c2d90829SJames Wright CeedQFunction qf_rhs_boundary; 153c2d90829SJames Wright DMLabel face_sets_label; 154c2d90829SJames Wright PetscInt num_face_set_values, *face_set_values; 155c2d90829SJames Wright CeedInt q_data_size; 156c2d90829SJames Wright 157c2d90829SJames Wright // -- Build RHS operator 158c2d90829SJames Wright switch (dim) { 159c2d90829SJames Wright case 2: 160c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_AdvDif_2D, DivDiffusiveFluxBoundaryRHS_AdvDif_2D_loc, 161c2d90829SJames Wright &qf_rhs_boundary)); 162c2d90829SJames Wright break; 163c2d90829SJames Wright case 3: 164c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_AdvDif_3D, DivDiffusiveFluxBoundaryRHS_AdvDif_3D_loc, 165c2d90829SJames Wright &qf_rhs_boundary)); 166c2d90829SJames Wright break; 167c2d90829SJames Wright } 168c2d90829SJames Wright 1690c373b74SJames Wright PetscCall(QDataBoundaryGradientGetNumComponents(honee->dm, &q_data_size)); 170c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_boundary, advection_qfctx)); 171c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 172c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "qdata", q_data_size, CEED_EVAL_NONE)); 173c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_boundary, "diffusive flux RHS", projection->num_comp, CEED_EVAL_INTERP)); 174c2d90829SJames Wright 175c2d90829SJames Wright PetscCall(DMGetLabel(projection->dm, "Face Sets", &face_sets_label)); 176c2d90829SJames Wright PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_sets_label, &num_face_set_values, &face_set_values)); 177c2d90829SJames Wright for (PetscInt f = 0; f < num_face_set_values; f++) { 178c2d90829SJames Wright DMLabel face_orientation_label; 179c2d90829SJames Wright PetscInt num_orientations_values, *orientation_values; 180c2d90829SJames Wright 181c2d90829SJames Wright { 182c2d90829SJames Wright char *face_orientation_label_name; 183c2d90829SJames Wright 184c2d90829SJames Wright PetscCall(DMPlexCreateFaceLabel(projection->dm, face_set_values[f], &face_orientation_label_name)); 185c2d90829SJames Wright PetscCall(DMGetLabel(projection->dm, face_orientation_label_name, &face_orientation_label)); 186c2d90829SJames Wright PetscCall(PetscFree(face_orientation_label_name)); 187c2d90829SJames Wright } 188c2d90829SJames Wright PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_orientation_label, &num_orientations_values, &orientation_values)); 189c2d90829SJames Wright for (PetscInt o = 0; o < num_orientations_values; o++) { 190c2d90829SJames Wright CeedOperator op_rhs_boundary; 191c2d90829SJames Wright CeedBasis basis_q, basis_diff_flux_boundary; 192c2d90829SJames Wright CeedElemRestriction elem_restr_qdata, elem_restr_q, elem_restr_diff_flux_boundary; 193c2d90829SJames Wright CeedVector q_data; 194c2d90829SJames Wright CeedInt q_data_size; 195c2d90829SJames Wright PetscInt orientation = orientation_values[o], dm_field_q = 0, height_cell = 0, height_face = 1; 196c2d90829SJames Wright 1970c373b74SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, face_orientation_label, orientation, height_cell, dm_field_q, &elem_restr_q)); 1980c373b74SJames Wright PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, honee->dm, face_orientation_label, orientation, orientation, dm_field_q, &basis_q)); 199c2d90829SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, 200c2d90829SJames Wright &elem_restr_diff_flux_boundary)); 201c2d90829SJames Wright PetscCall(CreateBasisFromPlex(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, &basis_diff_flux_boundary)); 202e3663b90SJames Wright PetscCall( 203e3663b90SJames Wright QDataBoundaryGradientGet(ceed, honee->dm, face_orientation_label, orientation, honee->x_coord, &elem_restr_qdata, &q_data, &q_data_size)); 204c2d90829SJames Wright 205c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_boundary, NULL, NULL, &op_rhs_boundary)); 206c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 207c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "qdata", elem_restr_qdata, CEED_BASIS_NONE, q_data)); 208c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "diffusive flux RHS", elem_restr_diff_flux_boundary, basis_diff_flux_boundary, 209c2d90829SJames Wright CEED_VECTOR_ACTIVE)); 210c2d90829SJames Wright 211c2d90829SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_boundary)); 212c2d90829SJames Wright 213c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_boundary)); 214c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qdata)); 215c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 216c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_boundary)); 217c2d90829SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 218c2d90829SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux_boundary)); 219c2d90829SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 220c2d90829SJames Wright } 221c2d90829SJames Wright PetscCall(PetscFree(orientation_values)); 222c2d90829SJames Wright } 223c2d90829SJames Wright PetscCall(PetscFree(face_set_values)); 224c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_boundary)); 225c2d90829SJames Wright } 226c2d90829SJames Wright 227c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&advection_qfctx)); 228c2d90829SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 229c2d90829SJames Wright } 230c2d90829SJames Wright 23140b78511SJames Wright /** 23240b78511SJames Wright @brief Create RHS CeedOperator for indirect projection of divergence of diffusive flux 23340b78511SJames Wright 2340c373b74SJames Wright @param[in] honee `Honee` context 23540b78511SJames Wright @param[in] diff_flux_proj `DivDiffFluxProjectionData` object 23640b78511SJames Wright @param[out] op_rhs Operator to calculate the RHS of the L^2 projection 23740b78511SJames Wright **/ 238e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Indirect_AdvDif(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) { 2390c373b74SJames Wright Ceed ceed = honee->ceed; 24040b78511SJames Wright NodalProjectionData projection = diff_flux_proj->projection; 24140b78511SJames Wright CeedBasis basis_diff_flux; 24240b78511SJames Wright CeedElemRestriction elem_restr_diff_flux, elem_restr_qd; 24340b78511SJames Wright CeedVector q_data; 24440b78511SJames Wright CeedInt num_comp_q, q_data_size; 24540b78511SJames Wright PetscInt dim; 24640b78511SJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 24740b78511SJames Wright DMLabel domain_label = NULL; 24840b78511SJames Wright CeedQFunction qf_rhs = NULL; 24940b78511SJames Wright CeedQFunctionContext advection_qfctx = NULL; 25040b78511SJames Wright 25140b78511SJames Wright PetscFunctionBeginUser; 25240b78511SJames Wright PetscCall(DMGetDimension(projection->dm, &dim)); 253e3663b90SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q)); 25440b78511SJames Wright 25540b78511SJames Wright { // Get advection-diffusion QF context 25640b78511SJames Wright CeedOperator *sub_ops; 25740b78511SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 25840b78511SJames Wright 2590c373b74SJames Wright if (honee->op_ifunction) PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_ifunction, &sub_ops)); 2600c373b74SJames Wright else PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_rhs_ctx->op, &sub_ops)); 26140b78511SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &advection_qfctx)); 26240b78511SJames Wright } 26340b78511SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, domain_label, label_value, height, dm_field, &elem_restr_diff_flux)); 26440b78511SJames Wright PetscCall(CreateBasisFromPlex(ceed, projection->dm, domain_label, label_value, height, dm_field, &basis_diff_flux)); 265e3663b90SJames Wright PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, &q_data, 266e3663b90SJames Wright &q_data_size)); 26740b78511SJames Wright 26840b78511SJames Wright switch (dim) { 26940b78511SJames Wright case 2: 27040b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_AdvDif_2D, DiffusiveFluxRHS_AdvDif_2D_loc, &qf_rhs)); 27140b78511SJames Wright break; 27240b78511SJames Wright case 3: 27340b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_AdvDif_3D, DiffusiveFluxRHS_AdvDif_3D_loc, &qf_rhs)); 27440b78511SJames Wright break; 27540b78511SJames Wright } 2760c373b74SJames Wright PetscCheck(qf_rhs, honee->comm, PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim); 27740b78511SJames Wright 27840b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs, advection_qfctx)); 27940b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 28040b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "qdata", q_data_size, CEED_EVAL_NONE)); 28140b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs, "F_diff RHS", projection->num_comp, CEED_EVAL_INTERP)); 28240b78511SJames Wright 28340b78511SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, op_rhs)); 284e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 28540b78511SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 28640b78511SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "F_diff RHS", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE)); 28740b78511SJames Wright 28840b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs)); 28940b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&advection_qfctx)); 29040b78511SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 29140b78511SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 29240b78511SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 29340b78511SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux)); 29440b78511SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 29540b78511SJames Wright } 29640b78511SJames Wright 297991aef52SJames Wright PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 2985f636aeaSJames Wright AdvDifWindType wind_type; 2995f636aeaSJames Wright AdvDifICType advectionic_type; 300108c6689SJames Wright AdvDifBubbleContinuityType bubble_continuity_type = -1; 301a515125bSLeila Ghaffari StabilizationType stab; 30257272ee0SJames Wright StabilizationTauType stab_tau; 3033636f6a4SJames Wright SetupContextAdv setup_context; 3040c373b74SJames Wright Honee honee = *(Honee *)ctx; 3050c373b74SJames Wright MPI_Comm comm = honee->comm; 3060c373b74SJames Wright Ceed ceed = honee->ceed; 307a515125bSLeila Ghaffari PetscBool implicit; 30815a3537eSJed Brown AdvectionContext advection_ctx; 309e07531f7SJames Wright CeedQFunctionContext advection_qfctx; 3109529d636SJames Wright PetscInt dim; 311a515125bSLeila Ghaffari 31215a3537eSJed Brown PetscFunctionBeginUser; 3132b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &setup_context)); 3142b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &advection_ctx)); 3159529d636SJames Wright PetscCall(DMGetDimension(dm, &dim)); 316a515125bSLeila Ghaffari 317a515125bSLeila Ghaffari // ------------------------------------------------------ 318a515125bSLeila Ghaffari // SET UP ADVECTION 319a515125bSLeila Ghaffari // ------------------------------------------------------ 32028160fc2SJames Wright problem->print_info = PRINT_ADVECTION; 32128160fc2SJames Wright problem->jac_data_size_vol = 0; 32228160fc2SJames Wright problem->jac_data_size_sur = 0; 3239529d636SJames Wright switch (dim) { 3249529d636SJames Wright case 2: 325e07531f7SJames Wright problem->ics.qf_func_ptr = ICsAdvection2d; 326e07531f7SJames Wright problem->ics.qf_loc = ICsAdvection2d_loc; 327e07531f7SJames Wright problem->apply_vol_rhs.qf_func_ptr = RHS_Advection2d; 328e07531f7SJames Wright problem->apply_vol_rhs.qf_loc = RHS_Advection2d_loc; 329e07531f7SJames Wright problem->apply_vol_ifunction.qf_func_ptr = IFunction_Advection2d; 330e07531f7SJames Wright problem->apply_vol_ifunction.qf_loc = IFunction_Advection2d_loc; 331e07531f7SJames Wright problem->apply_inflow.qf_func_ptr = Advection2d_InOutFlow; 332e07531f7SJames Wright problem->apply_inflow.qf_loc = Advection2d_InOutFlow_loc; 33358ce1233SJames Wright problem->compute_exact_solution_error = PETSC_TRUE; 3349529d636SJames Wright break; 3359529d636SJames Wright case 3: 336e07531f7SJames Wright problem->ics.qf_func_ptr = ICsAdvection; 337e07531f7SJames Wright problem->ics.qf_loc = ICsAdvection_loc; 338e07531f7SJames Wright problem->apply_vol_rhs.qf_func_ptr = RHS_Advection; 339e07531f7SJames Wright problem->apply_vol_rhs.qf_loc = RHS_Advection_loc; 340e07531f7SJames Wright problem->apply_vol_ifunction.qf_func_ptr = IFunction_Advection; 341e07531f7SJames Wright problem->apply_vol_ifunction.qf_loc = IFunction_Advection_loc; 342e07531f7SJames Wright problem->apply_inflow.qf_func_ptr = Advection_InOutFlow; 343e07531f7SJames Wright problem->apply_inflow.qf_loc = Advection_InOutFlow_loc; 34458ce1233SJames Wright problem->compute_exact_solution_error = PETSC_FALSE; 3459529d636SJames Wright break; 3469529d636SJames Wright } 347a515125bSLeila Ghaffari 3480c373b74SJames Wright PetscCall(DivDiffFluxProjectionCreate(honee, 1, &honee->diff_flux_proj)); 3490c373b74SJames Wright if (honee->diff_flux_proj) { 3500c373b74SJames Wright DivDiffFluxProjectionData diff_flux_proj = honee->diff_flux_proj; 351c2d90829SJames Wright NodalProjectionData projection = diff_flux_proj->projection; 352c2d90829SJames Wright 353c2d90829SJames Wright diff_flux_proj->CreateRHSOperator_Direct = DivDiffFluxProjectionCreateRHS_Direct_AdvDif; 35440b78511SJames Wright diff_flux_proj->CreateRHSOperator_Indirect = DivDiffFluxProjectionCreateRHS_Indirect_AdvDif; 355c2d90829SJames Wright 3560c373b74SJames Wright switch (honee->diff_flux_proj->method) { 357c2d90829SJames Wright case DIV_DIFF_FLUX_PROJ_DIRECT: { 358c2d90829SJames Wright PetscSection section; 359c2d90829SJames Wright 360c2d90829SJames Wright PetscCall(DMGetLocalSection(projection->dm, §ion)); 361c2d90829SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 362c2d90829SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "DivDiffusiveFlux_Scalar")); 363c2d90829SJames Wright } break; 364c2d90829SJames Wright case DIV_DIFF_FLUX_PROJ_INDIRECT: { 365c2d90829SJames Wright PetscSection section; 366c2d90829SJames Wright 367c2d90829SJames Wright PetscCall(DMGetLocalSection(projection->dm, §ion)); 368c2d90829SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 369c2d90829SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "DiffusiveFlux_ScalarX")); 370c2d90829SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "DiffusiveFlux_ScalarY")); 37140b78511SJames Wright if (dim >= 3) PetscCall(PetscSectionSetComponentName(section, 0, 2, "DiffusiveFlux_ScalarZ")); 372c2d90829SJames Wright } break; 373c2d90829SJames Wright case DIV_DIFF_FLUX_PROJ_NONE: 3740c373b74SJames Wright SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s", 3750c373b74SJames Wright DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]); 376c2d90829SJames Wright break; 377c2d90829SJames Wright } 378c2d90829SJames Wright } 379c2d90829SJames Wright 380a515125bSLeila Ghaffari // ------------------------------------------------------ 381a515125bSLeila Ghaffari // Create the libCEED context 382a515125bSLeila Ghaffari // ------------------------------------------------------ 383a515125bSLeila Ghaffari CeedScalar rc = 1000.; // m (Radius of bubble) 384a515125bSLeila Ghaffari CeedScalar CtauS = 0.; // dimensionless 385059d1c40SJames Wright PetscBool strong_form = PETSC_FALSE; 386a515125bSLeila Ghaffari CeedScalar E_wind = 1.e6; // J 3870c373b74SJames Wright CeedScalar Ctau_a = PetscPowScalarInt(honee->app_ctx->degree, 2); 3880c373b74SJames Wright CeedScalar Ctau_d = PetscPowScalarInt(honee->app_ctx->degree, 4); 38957272ee0SJames Wright CeedScalar Ctau_t = 0.; 390a515125bSLeila Ghaffari PetscReal wind[3] = {1., 0, 0}; // m/s 391c8d249deSJames Wright CeedScalar diffusion_coeff = 0.; 392a62be6baSJames Wright CeedScalar wave_frequency = 2 * M_PI; 393a62be6baSJames Wright CeedScalar wave_phase = 0; 394108c6689SJames Wright AdvDifWaveType wave_type = -1; 395*b4fd18dfSJames Wright PetscScalar bl_height_factor = 1.; 39687d3884fSJeremy L Thompson PetscReal domain_min[3], domain_max[3], domain_size[3] = {0.}; 3972b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 39866d54740SJames Wright for (PetscInt i = 0; i < dim; i++) domain_size[i] = domain_max[i] - domain_min[i]; 39905a512bdSLeila Ghaffari 400a515125bSLeila Ghaffari // ------------------------------------------------------ 401a515125bSLeila Ghaffari // Create the PETSc context 402a515125bSLeila Ghaffari // ------------------------------------------------------ 403a515125bSLeila Ghaffari PetscScalar meter = 1e-2; // 1 meter in scaled length units 404a515125bSLeila Ghaffari PetscScalar kilogram = 1e-6; // 1 kilogram in scaled mass units 405a515125bSLeila Ghaffari PetscScalar second = 1e-2; // 1 second in scaled time units 406a515125bSLeila Ghaffari PetscScalar Joule; 407a515125bSLeila Ghaffari 408a515125bSLeila Ghaffari // ------------------------------------------------------ 409a515125bSLeila Ghaffari // Command line Options 410a515125bSLeila Ghaffari // ------------------------------------------------------ 4111485969bSJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for ADVECTION problem", NULL); 412a515125bSLeila Ghaffari // -- Physics 4132b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL)); 414a515125bSLeila Ghaffari PetscBool translation; 4155f636aeaSJames Wright PetscCall(PetscOptionsEnum("-wind_type", "Wind type in Advection", NULL, AdvDifWindTypes, (PetscEnum)(wind_type = ADVDIF_WIND_ROTATION), 4165f636aeaSJames Wright (PetscEnum *)&wind_type, &translation)); 41766d54740SJames Wright PetscInt n = dim; 418a515125bSLeila Ghaffari PetscBool user_wind; 4192b916ea7SJeremy L Thompson PetscCall(PetscOptionsRealArray("-wind_translation", "Constant wind vector", NULL, wind, &n, &user_wind)); 420c8d249deSJames Wright PetscCall(PetscOptionsScalar("-diffusion_coeff", "Diffusion coefficient", NULL, diffusion_coeff, &diffusion_coeff, NULL)); 4212b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-CtauS", "Scale coefficient for tau (nondimensional)", NULL, CtauS, &CtauS, NULL)); 422059d1c40SJames Wright PetscCall(PetscOptionsBool("-strong_form", "Strong (true) or weak/integrated by parts (false) advection residual", NULL, strong_form, &strong_form, 423059d1c40SJames Wright NULL)); 4242b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-E_wind", "Total energy of inflow wind", NULL, E_wind, &E_wind, NULL)); 4255f636aeaSJames Wright PetscCall(PetscOptionsEnum("-advection_ic_type", "Initial condition for Advection problem", NULL, AdvDifICTypes, 4265f636aeaSJames Wright (PetscEnum)(advectionic_type = ADVDIF_IC_BUBBLE_SPHERE), (PetscEnum *)&advectionic_type, NULL)); 4272b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 42857272ee0SJames Wright PetscCall(PetscOptionsEnum("-stab_tau", "Stabilization constant, tau", NULL, StabilizationTauTypes, (PetscEnum)(stab_tau = STAB_TAU_CTAU), 42957272ee0SJames Wright (PetscEnum *)&stab_tau, NULL)); 43057272ee0SJames Wright PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 431fbabb365SJames Wright PetscCall(PetscOptionsScalar("-Ctau_a", "Coefficient for the stabilization, advection component", NULL, Ctau_a, &Ctau_a, NULL)); 432fbabb365SJames Wright PetscCall(PetscOptionsScalar("-Ctau_d", "Coefficient for the stabilization, diffusion component", NULL, Ctau_d, &Ctau_d, NULL)); 4332b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 434a515125bSLeila Ghaffari 435108c6689SJames Wright if (advectionic_type == ADVDIF_IC_WAVE) { 436108c6689SJames Wright PetscCall(PetscOptionsEnum("-wave_type", "Type of wave", NULL, AdvDifWaveTypes, (PetscEnum)(wave_type = ADVDIF_WAVE_SINE), 437108c6689SJames Wright (PetscEnum *)&wave_type, NULL)); 438a62be6baSJames Wright PetscCall(PetscOptionsScalar("-wave_frequency", "Frequency of sine wave", NULL, wave_frequency, &wave_frequency, NULL)); 439a62be6baSJames Wright PetscCall(PetscOptionsScalar("-wave_phase", "Length correction", NULL, wave_phase, &wave_phase, NULL)); 440108c6689SJames Wright } 441108c6689SJames Wright 442*b4fd18dfSJames Wright if (advectionic_type == ADVDIF_IC_BOUNDARY_LAYER) 443*b4fd18dfSJames Wright PetscCall(PetscOptionsScalar("-bl_height_factor", "Height of boundary layer in IC", NULL, bl_height_factor, &bl_height_factor, NULL)); 444*b4fd18dfSJames Wright 445108c6689SJames Wright if (advectionic_type == ADVDIF_IC_BUBBLE_CYLINDER || advectionic_type == ADVDIF_IC_BUBBLE_SPHERE) { 446108c6689SJames Wright bubble_continuity_type = dim == 3 ? ADVDIF_BUBBLE_CONTINUITY_SMOOTH : ADVDIF_BUBBLE_CONTINUITY_COSINE; 447108c6689SJames Wright PetscCall(PetscOptionsEnum("-bubble_continuity", "Smooth, back_sharp, or thick", NULL, AdvDifBubbleContinuityTypes, 448108c6689SJames Wright (PetscEnum)bubble_continuity_type, (PetscEnum *)&bubble_continuity_type, NULL)); 449108c6689SJames Wright } 450a62be6baSJames Wright 451a515125bSLeila Ghaffari // -- Units 4522b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 453a515125bSLeila Ghaffari meter = fabs(meter); 4542b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL)); 455a515125bSLeila Ghaffari kilogram = fabs(kilogram); 4562b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 457a515125bSLeila Ghaffari second = fabs(second); 458a515125bSLeila Ghaffari 459a515125bSLeila Ghaffari // -- Warnings 4605f636aeaSJames Wright if (wind_type == ADVDIF_WIND_ROTATION && user_wind) { 4612b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Use -wind_translation only with -wind_type translation\n")); 462a515125bSLeila Ghaffari } 4635f636aeaSJames Wright if (wind_type == ADVDIF_WIND_TRANSLATION && advectionic_type == ADVDIF_IC_BUBBLE_CYLINDER && wind[2] != 0.) { 464a515125bSLeila Ghaffari wind[2] = 0; 465c51f031aSJames Wright PetscCall( 466c51f031aSJames Wright PetscPrintf(comm, "Warning! Background wind in the z direction should be zero (-wind_translation x,x,0) with -advection_ic_type cylinder\n")); 467a515125bSLeila Ghaffari } 468a515125bSLeila Ghaffari if (stab == STAB_NONE && CtauS != 0) { 4692b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Use -CtauS only with -stab su or -stab supg\n")); 470a515125bSLeila Ghaffari } 4711485969bSJeremy L Thompson PetscOptionsEnd(); 472a515125bSLeila Ghaffari 473a78efa86SJames Wright if (stab == STAB_SUPG) problem->create_mass_operator = CreateKSPMassOperator_AdvectionStabilized; 474a78efa86SJames Wright 475a515125bSLeila Ghaffari // ------------------------------------------------------ 476a515125bSLeila Ghaffari // Set up the PETSc context 477a515125bSLeila Ghaffari // ------------------------------------------------------ 478a515125bSLeila Ghaffari // -- Define derived units 479a515125bSLeila Ghaffari Joule = kilogram * PetscSqr(meter) / PetscSqr(second); 480a515125bSLeila Ghaffari 4810c373b74SJames Wright honee->units->meter = meter; 4820c373b74SJames Wright honee->units->kilogram = kilogram; 4830c373b74SJames Wright honee->units->second = second; 4840c373b74SJames Wright honee->units->Joule = Joule; 485a515125bSLeila Ghaffari 486a515125bSLeila Ghaffari // ------------------------------------------------------ 487a515125bSLeila Ghaffari // Set up the libCEED context 488a515125bSLeila Ghaffari // ------------------------------------------------------ 489a515125bSLeila Ghaffari // -- Scale variables to desired units 490a515125bSLeila Ghaffari E_wind *= Joule; 491a515125bSLeila Ghaffari rc = fabs(rc) * meter; 49266d54740SJames Wright for (PetscInt i = 0; i < dim; i++) { 49305a512bdSLeila Ghaffari wind[i] *= (meter / second); 49405a512bdSLeila Ghaffari domain_size[i] *= meter; 49505a512bdSLeila Ghaffari } 496a515125bSLeila Ghaffari 497a515125bSLeila Ghaffari // -- Setup Context 498a515125bSLeila Ghaffari setup_context->rc = rc; 49905a512bdSLeila Ghaffari setup_context->lx = domain_size[0]; 50005a512bdSLeila Ghaffari setup_context->ly = domain_size[1]; 50166d54740SJames Wright setup_context->lz = dim == 3 ? domain_size[2] : 0.; 502a515125bSLeila Ghaffari setup_context->wind[0] = wind[0]; 503a515125bSLeila Ghaffari setup_context->wind[1] = wind[1]; 50466d54740SJames Wright setup_context->wind[2] = dim == 3 ? wind[2] : 0.; 505a515125bSLeila Ghaffari setup_context->wind_type = wind_type; 506c51f031aSJames Wright setup_context->initial_condition_type = advectionic_type; 507a515125bSLeila Ghaffari setup_context->bubble_continuity_type = bubble_continuity_type; 508a515125bSLeila Ghaffari setup_context->time = 0; 509a62be6baSJames Wright setup_context->wave_frequency = wave_frequency; 510a62be6baSJames Wright setup_context->wave_phase = wave_phase; 511a62be6baSJames Wright setup_context->wave_type = wave_type; 512*b4fd18dfSJames Wright setup_context->bl_height_factor = bl_height_factor; 513a515125bSLeila Ghaffari 514a515125bSLeila Ghaffari // -- QFunction Context 5150c373b74SJames Wright honee->phys->implicit = implicit; 51615a3537eSJed Brown advection_ctx->CtauS = CtauS; 51715a3537eSJed Brown advection_ctx->E_wind = E_wind; 51815a3537eSJed Brown advection_ctx->implicit = implicit; 51915a3537eSJed Brown advection_ctx->strong_form = strong_form; 52015a3537eSJed Brown advection_ctx->stabilization = stab; 52157272ee0SJames Wright advection_ctx->stabilization_tau = stab_tau; 52257272ee0SJames Wright advection_ctx->Ctau_a = Ctau_a; 523fbabb365SJames Wright advection_ctx->Ctau_d = Ctau_d; 52457272ee0SJames Wright advection_ctx->Ctau_t = Ctau_t; 525c8d249deSJames Wright advection_ctx->diffusion_coeff = diffusion_coeff; 5260c373b74SJames Wright advection_ctx->divFdiff_method = honee->app_ctx->divFdiffproj_method; 527a515125bSLeila Ghaffari 5280c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &problem->ics.qfctx)); 529e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(problem->ics.qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 530e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfctx, CEED_MEM_HOST, FreeContextPetsc)); 531a515125bSLeila Ghaffari 5320c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &advection_qfctx)); 533e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(advection_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*advection_ctx), advection_ctx)); 534e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(advection_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 535e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(advection_qfctx, "timestep size", offsetof(struct AdvectionContext_, dt), 1, 53657272ee0SJames Wright "Size of timestep, delta t")); 537e07531f7SJames Wright problem->apply_vol_rhs.qfctx = advection_qfctx; 538e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &problem->apply_vol_ifunction.qfctx)); 539e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &problem->apply_inflow.qfctx)); 540d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 541a515125bSLeila Ghaffari } 542a515125bSLeila Ghaffari 5430c373b74SJames Wright PetscErrorCode PRINT_ADVECTION(Honee honee, ProblemData problem, AppCtx app_ctx) { 5440c373b74SJames Wright MPI_Comm comm = honee->comm; 5450c373b74SJames Wright Ceed ceed = honee->ceed; 5463636f6a4SJames Wright SetupContextAdv setup_ctx; 54715a3537eSJed Brown AdvectionContext advection_ctx; 54866d54740SJames Wright PetscInt dim; 549a515125bSLeila Ghaffari 55015a3537eSJed Brown PetscFunctionBeginUser; 5510c373b74SJames Wright PetscCall(DMGetDimension(honee->dm, &dim)); 552e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfctx, CEED_MEM_HOST, &setup_ctx)); 553e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &advection_ctx)); 5542b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, 555a515125bSLeila Ghaffari " Problem:\n" 556a515125bSLeila Ghaffari " Problem Name : %s\n" 557a515125bSLeila Ghaffari " Stabilization : %s\n" 5585bcc2007SJames Wright " Stabilization Tau : %s\n" 559a515125bSLeila Ghaffari " Wind Type : %s\n", 5605bcc2007SJames Wright app_ctx->problem_name, StabilizationTypes[advection_ctx->stabilization], 5615f636aeaSJames Wright StabilizationTauTypes[advection_ctx->stabilization_tau], AdvDifWindTypes[setup_ctx->wind_type])); 562a515125bSLeila Ghaffari 5635f636aeaSJames Wright if (setup_ctx->wind_type == ADVDIF_WIND_TRANSLATION) { 564c781aee8SJames Wright CeedScalar *wind = setup_ctx->wind; 56566d54740SJames Wright switch (dim) { 5669529d636SJames Wright case 2: 567c781aee8SJames Wright PetscCall(PetscPrintf(comm, " Background Wind : %f,%f\n", wind[0], wind[1])); 5689529d636SJames Wright break; 5699529d636SJames Wright case 3: 570c781aee8SJames Wright PetscCall(PetscPrintf(comm, " Background Wind : %f,%f,%f\n", wind[0], wind[1], wind[2])); 5719529d636SJames Wright break; 5729529d636SJames Wright } 573a515125bSLeila Ghaffari } 574c781aee8SJames Wright 5755f636aeaSJames Wright PetscCall(PetscPrintf(comm, " Initial Condition Type : %s\n", AdvDifICTypes[setup_ctx->initial_condition_type])); 576c781aee8SJames Wright switch (setup_ctx->initial_condition_type) { 5775f636aeaSJames Wright case ADVDIF_IC_SKEW: 5785f636aeaSJames Wright case ADVDIF_IC_COSINE_HILL: 5793d1afcc1SJames Wright case ADVDIF_IC_BOUNDARY_LAYER: 580c781aee8SJames Wright break; 5815f636aeaSJames Wright case ADVDIF_IC_BUBBLE_SPHERE: 5825f636aeaSJames Wright case ADVDIF_IC_BUBBLE_CYLINDER: 5835f636aeaSJames Wright PetscCall(PetscPrintf(comm, " Bubble Continuity : %s\n", AdvDifBubbleContinuityTypes[setup_ctx->bubble_continuity_type])); 584c781aee8SJames Wright break; 5855f636aeaSJames Wright case ADVDIF_IC_WAVE: 586c781aee8SJames Wright PetscCall(PetscPrintf(comm, " Wave Type : %s\n", AdvDifWaveTypes[setup_ctx->wave_type])); 587c781aee8SJames Wright break; 588c781aee8SJames Wright } 589c781aee8SJames Wright 590e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfctx, &setup_ctx)); 591e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &advection_ctx)); 592d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 593a515125bSLeila Ghaffari } 594