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; 21a78efa86SJames Wright CeedElemRestriction elem_restr_q, elem_restr_qd_i; 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; 31a78efa86SJames Wright CeedOperatorField 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)); 35a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field)); 36a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q)); 37a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetBasis(field, &basis_q)); 38a78efa86SJames Wright 39a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field)); 40a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_qd_i)); 41a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetVector(field, &q_data)); 42a78efa86SJames Wright 43236a8ba6SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qfctx)); 44a78efa86SJames Wright } 45a78efa86SJames Wright 46a78efa86SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q)); 47a78efa86SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size)); 48a78efa86SJames Wright 49a78efa86SJames Wright switch (dim) { 50a78efa86SJames Wright case 2: 51a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection2D, MassFunction_Advection2D_loc, &qf_mass)); 52a78efa86SJames Wright break; 53a78efa86SJames Wright case 3: 54a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection, MassFunction_Advection_loc, &qf_mass)); 55a78efa86SJames Wright break; 56a78efa86SJames Wright } 57a78efa86SJames Wright 58236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qfctx)); 59a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0)); 60a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP)); 61a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP)); 62a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE)); 63a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP)); 64a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD)); 65a78efa86SJames Wright 66a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass)); 67a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 680c373b74SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, honee->q_ceed)); 69a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data)); 70a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 71a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 72a78efa86SJames Wright 73*fff85bd3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 74*fff85bd3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_i)); 75*fff85bd3SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 76*fff85bd3SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 77236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&qfctx)); 78a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 79a78efa86SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 80a78efa86SJames Wright } 81a78efa86SJames Wright 82c2d90829SJames Wright /** 83c2d90829SJames Wright @brief Create RHS CeedOperator for direct projection of divergence of diffusive flux 84c2d90829SJames Wright 850c373b74SJames Wright @param[in] honee `Honee` context 86c2d90829SJames Wright @param[in] diff_flux_proj `DivDiffFluxProjectionData` object 87c2d90829SJames Wright @param[out] op_rhs Operator to calculate the RHS of the L^2 projection 88c2d90829SJames Wright **/ 89e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Direct_AdvDif(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) { 900c373b74SJames Wright Ceed ceed = honee->ceed; 91c2d90829SJames Wright NodalProjectionData projection = diff_flux_proj->projection; 92c2d90829SJames Wright CeedInt num_comp_q; 93c2d90829SJames Wright PetscInt dim, label_value = 0; 94c2d90829SJames Wright DMLabel domain_label = NULL; 95c2d90829SJames Wright CeedQFunctionContext advection_qfctx = NULL; 96c2d90829SJames Wright 97c2d90829SJames Wright PetscFunctionBeginUser; 98c2d90829SJames Wright // -- Get Pre-requisite things 99c2d90829SJames Wright PetscCall(DMGetDimension(projection->dm, &dim)); 100e3663b90SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q)); 101c2d90829SJames Wright 10240b78511SJames Wright { // Get advection-diffusion QF context 103c2d90829SJames Wright CeedOperator *sub_ops; 104c2d90829SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 105c2d90829SJames Wright 1060c373b74SJames Wright if (honee->op_ifunction) PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_ifunction, &sub_ops)); 1070c373b74SJames Wright else PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_rhs_ctx->op, &sub_ops)); 108c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &advection_qfctx)); 109c2d90829SJames Wright } 110c2d90829SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, op_rhs)); 111c2d90829SJames Wright { // Add the volume integral CeedOperator 112c2d90829SJames Wright CeedQFunction qf_rhs_volume = NULL; 113c2d90829SJames Wright CeedOperator op_rhs_volume; 114c2d90829SJames Wright CeedVector q_data; 115c2d90829SJames Wright CeedElemRestriction elem_restr_qd, elem_restr_diff_flux_volume = NULL; 116c2d90829SJames Wright CeedBasis basis_diff_flux = NULL; 117c2d90829SJames Wright CeedInt q_data_size; 118c2d90829SJames Wright 119c2d90829SJames Wright PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_diff_flux_volume, &basis_diff_flux, NULL, NULL)); 120e3663b90SJames 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, 121e3663b90SJames Wright &q_data_size)); 122c2d90829SJames Wright switch (dim) { 123c2d90829SJames Wright case 2: 124c2d90829SJames Wright PetscCallCeed( 125c2d90829SJames Wright ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_AdvDif_2D, DivDiffusiveFluxVolumeRHS_AdvDif_2D_loc, &qf_rhs_volume)); 126c2d90829SJames Wright break; 127c2d90829SJames Wright case 3: 128c2d90829SJames Wright PetscCallCeed( 129c2d90829SJames Wright ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_AdvDif_3D, DivDiffusiveFluxVolumeRHS_AdvDif_3D_loc, &qf_rhs_volume)); 130c2d90829SJames Wright break; 131c2d90829SJames Wright } 1320c373b74SJames Wright PetscCheck(qf_rhs_volume, honee->comm, PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim); 133c2d90829SJames Wright 134c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_volume, advection_qfctx)); 135c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 136c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "qdata", q_data_size, CEED_EVAL_NONE)); 137c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_volume, "diffusive flux RHS", projection->num_comp * dim, CEED_EVAL_GRAD)); 138c2d90829SJames Wright 139c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_volume, NULL, NULL, &op_rhs_volume)); 140e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 141c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 142c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "diffusive flux RHS", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE)); 143c2d90829SJames Wright 144c2d90829SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_volume)); 145c2d90829SJames Wright 146c2d90829SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 147c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 148c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_volume)); 149c2d90829SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 150c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_volume)); 151c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_volume)); 152c2d90829SJames Wright } 153c2d90829SJames Wright 154c2d90829SJames Wright { // Add the boundary integral CeedOperator 155c2d90829SJames Wright CeedQFunction qf_rhs_boundary; 156c2d90829SJames Wright DMLabel face_sets_label; 157c2d90829SJames Wright PetscInt num_face_set_values, *face_set_values; 158c2d90829SJames Wright CeedInt q_data_size; 159c2d90829SJames Wright 160c2d90829SJames Wright // -- Build RHS operator 161c2d90829SJames Wright switch (dim) { 162c2d90829SJames Wright case 2: 163c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_AdvDif_2D, DivDiffusiveFluxBoundaryRHS_AdvDif_2D_loc, 164c2d90829SJames Wright &qf_rhs_boundary)); 165c2d90829SJames Wright break; 166c2d90829SJames Wright case 3: 167c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_AdvDif_3D, DivDiffusiveFluxBoundaryRHS_AdvDif_3D_loc, 168c2d90829SJames Wright &qf_rhs_boundary)); 169c2d90829SJames Wright break; 170c2d90829SJames Wright } 171c2d90829SJames Wright 1720c373b74SJames Wright PetscCall(QDataBoundaryGradientGetNumComponents(honee->dm, &q_data_size)); 173c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_boundary, advection_qfctx)); 174c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 175c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "qdata", q_data_size, CEED_EVAL_NONE)); 176c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_boundary, "diffusive flux RHS", projection->num_comp, CEED_EVAL_INTERP)); 177c2d90829SJames Wright 178c2d90829SJames Wright PetscCall(DMGetLabel(projection->dm, "Face Sets", &face_sets_label)); 179c2d90829SJames Wright PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_sets_label, &num_face_set_values, &face_set_values)); 180c2d90829SJames Wright for (PetscInt f = 0; f < num_face_set_values; f++) { 181c2d90829SJames Wright DMLabel face_orientation_label; 182c2d90829SJames Wright PetscInt num_orientations_values, *orientation_values; 183c2d90829SJames Wright 184c2d90829SJames Wright { 185c2d90829SJames Wright char *face_orientation_label_name; 186c2d90829SJames Wright 187c2d90829SJames Wright PetscCall(DMPlexCreateFaceLabel(projection->dm, face_set_values[f], &face_orientation_label_name)); 188c2d90829SJames Wright PetscCall(DMGetLabel(projection->dm, face_orientation_label_name, &face_orientation_label)); 189c2d90829SJames Wright PetscCall(PetscFree(face_orientation_label_name)); 190c2d90829SJames Wright } 191c2d90829SJames Wright PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_orientation_label, &num_orientations_values, &orientation_values)); 192c2d90829SJames Wright for (PetscInt o = 0; o < num_orientations_values; o++) { 193c2d90829SJames Wright CeedOperator op_rhs_boundary; 194c2d90829SJames Wright CeedBasis basis_q, basis_diff_flux_boundary; 195c2d90829SJames Wright CeedElemRestriction elem_restr_qdata, elem_restr_q, elem_restr_diff_flux_boundary; 196c2d90829SJames Wright CeedVector q_data; 197c2d90829SJames Wright CeedInt q_data_size; 198c2d90829SJames Wright PetscInt orientation = orientation_values[o], dm_field_q = 0, height_cell = 0, height_face = 1; 199c2d90829SJames Wright 2000c373b74SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, face_orientation_label, orientation, height_cell, dm_field_q, &elem_restr_q)); 2010c373b74SJames Wright PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, honee->dm, face_orientation_label, orientation, orientation, dm_field_q, &basis_q)); 202c2d90829SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, 203c2d90829SJames Wright &elem_restr_diff_flux_boundary)); 204c2d90829SJames Wright PetscCall(CreateBasisFromPlex(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, &basis_diff_flux_boundary)); 205e3663b90SJames Wright PetscCall( 206e3663b90SJames Wright QDataBoundaryGradientGet(ceed, honee->dm, face_orientation_label, orientation, honee->x_coord, &elem_restr_qdata, &q_data, &q_data_size)); 207c2d90829SJames Wright 208c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_boundary, NULL, NULL, &op_rhs_boundary)); 209c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 210c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "qdata", elem_restr_qdata, CEED_BASIS_NONE, q_data)); 211c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "diffusive flux RHS", elem_restr_diff_flux_boundary, basis_diff_flux_boundary, 212c2d90829SJames Wright CEED_VECTOR_ACTIVE)); 213c2d90829SJames Wright 214c2d90829SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_boundary)); 215c2d90829SJames Wright 216c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_boundary)); 217c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qdata)); 218c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 219c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_boundary)); 220c2d90829SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 221c2d90829SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux_boundary)); 222c2d90829SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 223c2d90829SJames Wright } 224c2d90829SJames Wright PetscCall(PetscFree(orientation_values)); 225c2d90829SJames Wright } 226c2d90829SJames Wright PetscCall(PetscFree(face_set_values)); 227c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_boundary)); 228c2d90829SJames Wright } 229c2d90829SJames Wright 230c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&advection_qfctx)); 231c2d90829SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 232c2d90829SJames Wright } 233c2d90829SJames Wright 23440b78511SJames Wright /** 23540b78511SJames Wright @brief Create RHS CeedOperator for indirect projection of divergence of diffusive flux 23640b78511SJames Wright 2370c373b74SJames Wright @param[in] honee `Honee` context 23840b78511SJames Wright @param[in] diff_flux_proj `DivDiffFluxProjectionData` object 23940b78511SJames Wright @param[out] op_rhs Operator to calculate the RHS of the L^2 projection 24040b78511SJames Wright **/ 241e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Indirect_AdvDif(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) { 2420c373b74SJames Wright Ceed ceed = honee->ceed; 24340b78511SJames Wright NodalProjectionData projection = diff_flux_proj->projection; 24440b78511SJames Wright CeedBasis basis_diff_flux; 24540b78511SJames Wright CeedElemRestriction elem_restr_diff_flux, elem_restr_qd; 24640b78511SJames Wright CeedVector q_data; 24740b78511SJames Wright CeedInt num_comp_q, q_data_size; 24840b78511SJames Wright PetscInt dim; 24940b78511SJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 25040b78511SJames Wright DMLabel domain_label = NULL; 25140b78511SJames Wright CeedQFunction qf_rhs = NULL; 25240b78511SJames Wright CeedQFunctionContext advection_qfctx = NULL; 25340b78511SJames Wright 25440b78511SJames Wright PetscFunctionBeginUser; 25540b78511SJames Wright PetscCall(DMGetDimension(projection->dm, &dim)); 256e3663b90SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q)); 25740b78511SJames Wright 25840b78511SJames Wright { // Get advection-diffusion QF context 25940b78511SJames Wright CeedOperator *sub_ops; 26040b78511SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 26140b78511SJames Wright 2620c373b74SJames Wright if (honee->op_ifunction) PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_ifunction, &sub_ops)); 2630c373b74SJames Wright else PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_rhs_ctx->op, &sub_ops)); 26440b78511SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &advection_qfctx)); 26540b78511SJames Wright } 26640b78511SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, domain_label, label_value, height, dm_field, &elem_restr_diff_flux)); 26740b78511SJames Wright PetscCall(CreateBasisFromPlex(ceed, projection->dm, domain_label, label_value, height, dm_field, &basis_diff_flux)); 268e3663b90SJames 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, 269e3663b90SJames Wright &q_data_size)); 27040b78511SJames Wright 27140b78511SJames Wright switch (dim) { 27240b78511SJames Wright case 2: 27340b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_AdvDif_2D, DiffusiveFluxRHS_AdvDif_2D_loc, &qf_rhs)); 27440b78511SJames Wright break; 27540b78511SJames Wright case 3: 27640b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_AdvDif_3D, DiffusiveFluxRHS_AdvDif_3D_loc, &qf_rhs)); 27740b78511SJames Wright break; 27840b78511SJames Wright } 2790c373b74SJames Wright PetscCheck(qf_rhs, honee->comm, PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim); 28040b78511SJames Wright 28140b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs, advection_qfctx)); 28240b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 28340b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "qdata", q_data_size, CEED_EVAL_NONE)); 28440b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs, "F_diff RHS", projection->num_comp, CEED_EVAL_INTERP)); 28540b78511SJames Wright 28640b78511SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, op_rhs)); 287e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 28840b78511SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 28940b78511SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "F_diff RHS", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE)); 29040b78511SJames Wright 29140b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs)); 29240b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&advection_qfctx)); 29340b78511SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 29440b78511SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 29540b78511SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 29640b78511SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux)); 29740b78511SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 29840b78511SJames Wright } 29940b78511SJames Wright 300991aef52SJames Wright PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 3015f636aeaSJames Wright AdvDifWindType wind_type; 3025f636aeaSJames Wright AdvDifICType advectionic_type; 303108c6689SJames Wright AdvDifBubbleContinuityType bubble_continuity_type = -1; 304a515125bSLeila Ghaffari StabilizationType stab; 30557272ee0SJames Wright StabilizationTauType stab_tau; 3063636f6a4SJames Wright SetupContextAdv setup_context; 3070c373b74SJames Wright Honee honee = *(Honee *)ctx; 3080c373b74SJames Wright MPI_Comm comm = honee->comm; 3090c373b74SJames Wright Ceed ceed = honee->ceed; 310a515125bSLeila Ghaffari PetscBool implicit; 31115a3537eSJed Brown AdvectionContext advection_ctx; 312e07531f7SJames Wright CeedQFunctionContext advection_qfctx; 3139529d636SJames Wright PetscInt dim; 314a515125bSLeila Ghaffari 31515a3537eSJed Brown PetscFunctionBeginUser; 3162b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &setup_context)); 3172b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &advection_ctx)); 3189529d636SJames Wright PetscCall(DMGetDimension(dm, &dim)); 319a515125bSLeila Ghaffari 320a515125bSLeila Ghaffari // ------------------------------------------------------ 321a515125bSLeila Ghaffari // SET UP ADVECTION 322a515125bSLeila Ghaffari // ------------------------------------------------------ 32328160fc2SJames Wright problem->print_info = PRINT_ADVECTION; 32428160fc2SJames Wright problem->jac_data_size_vol = 0; 32528160fc2SJames Wright problem->jac_data_size_sur = 0; 3269529d636SJames Wright switch (dim) { 3279529d636SJames Wright case 2: 328e07531f7SJames Wright problem->ics.qf_func_ptr = ICsAdvection2d; 329e07531f7SJames Wright problem->ics.qf_loc = ICsAdvection2d_loc; 330e07531f7SJames Wright problem->apply_vol_rhs.qf_func_ptr = RHS_Advection2d; 331e07531f7SJames Wright problem->apply_vol_rhs.qf_loc = RHS_Advection2d_loc; 332e07531f7SJames Wright problem->apply_vol_ifunction.qf_func_ptr = IFunction_Advection2d; 333e07531f7SJames Wright problem->apply_vol_ifunction.qf_loc = IFunction_Advection2d_loc; 334e07531f7SJames Wright problem->apply_inflow.qf_func_ptr = Advection2d_InOutFlow; 335e07531f7SJames Wright problem->apply_inflow.qf_loc = Advection2d_InOutFlow_loc; 33658ce1233SJames Wright problem->compute_exact_solution_error = PETSC_TRUE; 3379529d636SJames Wright break; 3389529d636SJames Wright case 3: 339e07531f7SJames Wright problem->ics.qf_func_ptr = ICsAdvection; 340e07531f7SJames Wright problem->ics.qf_loc = ICsAdvection_loc; 341e07531f7SJames Wright problem->apply_vol_rhs.qf_func_ptr = RHS_Advection; 342e07531f7SJames Wright problem->apply_vol_rhs.qf_loc = RHS_Advection_loc; 343e07531f7SJames Wright problem->apply_vol_ifunction.qf_func_ptr = IFunction_Advection; 344e07531f7SJames Wright problem->apply_vol_ifunction.qf_loc = IFunction_Advection_loc; 345e07531f7SJames Wright problem->apply_inflow.qf_func_ptr = Advection_InOutFlow; 346e07531f7SJames Wright problem->apply_inflow.qf_loc = Advection_InOutFlow_loc; 34758ce1233SJames Wright problem->compute_exact_solution_error = PETSC_FALSE; 3489529d636SJames Wright break; 3499529d636SJames Wright } 350a515125bSLeila Ghaffari 3510c373b74SJames Wright PetscCall(DivDiffFluxProjectionCreate(honee, 1, &honee->diff_flux_proj)); 3520c373b74SJames Wright if (honee->diff_flux_proj) { 3530c373b74SJames Wright DivDiffFluxProjectionData diff_flux_proj = honee->diff_flux_proj; 354c2d90829SJames Wright NodalProjectionData projection = diff_flux_proj->projection; 355c2d90829SJames Wright 356c2d90829SJames Wright diff_flux_proj->CreateRHSOperator_Direct = DivDiffFluxProjectionCreateRHS_Direct_AdvDif; 35740b78511SJames Wright diff_flux_proj->CreateRHSOperator_Indirect = DivDiffFluxProjectionCreateRHS_Indirect_AdvDif; 358c2d90829SJames Wright 3590c373b74SJames Wright switch (honee->diff_flux_proj->method) { 360c2d90829SJames Wright case DIV_DIFF_FLUX_PROJ_DIRECT: { 361c2d90829SJames Wright PetscSection section; 362c2d90829SJames Wright 363c2d90829SJames Wright PetscCall(DMGetLocalSection(projection->dm, §ion)); 364c2d90829SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 365c2d90829SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "DivDiffusiveFlux_Scalar")); 366c2d90829SJames Wright } break; 367c2d90829SJames Wright case DIV_DIFF_FLUX_PROJ_INDIRECT: { 368c2d90829SJames Wright PetscSection section; 369c2d90829SJames Wright 370c2d90829SJames Wright PetscCall(DMGetLocalSection(projection->dm, §ion)); 371c2d90829SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 372c2d90829SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "DiffusiveFlux_ScalarX")); 373c2d90829SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "DiffusiveFlux_ScalarY")); 37440b78511SJames Wright if (dim >= 3) PetscCall(PetscSectionSetComponentName(section, 0, 2, "DiffusiveFlux_ScalarZ")); 375c2d90829SJames Wright } break; 376c2d90829SJames Wright case DIV_DIFF_FLUX_PROJ_NONE: 3770c373b74SJames Wright SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s", 3780c373b74SJames Wright DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]); 379c2d90829SJames Wright break; 380c2d90829SJames Wright } 381c2d90829SJames Wright } 382c2d90829SJames Wright 383a515125bSLeila Ghaffari // ------------------------------------------------------ 384a515125bSLeila Ghaffari // Create the libCEED context 385a515125bSLeila Ghaffari // ------------------------------------------------------ 386a515125bSLeila Ghaffari CeedScalar rc = 1000.; // m (Radius of bubble) 387a515125bSLeila Ghaffari CeedScalar CtauS = 0.; // dimensionless 388059d1c40SJames Wright PetscBool strong_form = PETSC_FALSE; 389a515125bSLeila Ghaffari CeedScalar E_wind = 1.e6; // J 3900c373b74SJames Wright CeedScalar Ctau_a = PetscPowScalarInt(honee->app_ctx->degree, 2); 3910c373b74SJames Wright CeedScalar Ctau_d = PetscPowScalarInt(honee->app_ctx->degree, 4); 39257272ee0SJames Wright CeedScalar Ctau_t = 0.; 393a515125bSLeila Ghaffari PetscReal wind[3] = {1., 0, 0}; // m/s 394c8d249deSJames Wright CeedScalar diffusion_coeff = 0.; 395a62be6baSJames Wright CeedScalar wave_frequency = 2 * M_PI; 396a62be6baSJames Wright CeedScalar wave_phase = 0; 397108c6689SJames Wright AdvDifWaveType wave_type = -1; 39887d3884fSJeremy L Thompson PetscReal domain_min[3], domain_max[3], domain_size[3] = {0.}; 3992b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 40066d54740SJames Wright for (PetscInt i = 0; i < dim; i++) domain_size[i] = domain_max[i] - domain_min[i]; 40105a512bdSLeila Ghaffari 402a515125bSLeila Ghaffari // ------------------------------------------------------ 403a515125bSLeila Ghaffari // Create the PETSc context 404a515125bSLeila Ghaffari // ------------------------------------------------------ 405a515125bSLeila Ghaffari PetscScalar meter = 1e-2; // 1 meter in scaled length units 406a515125bSLeila Ghaffari PetscScalar kilogram = 1e-6; // 1 kilogram in scaled mass units 407a515125bSLeila Ghaffari PetscScalar second = 1e-2; // 1 second in scaled time units 408a515125bSLeila Ghaffari PetscScalar Joule; 409a515125bSLeila Ghaffari 410a515125bSLeila Ghaffari // ------------------------------------------------------ 411a515125bSLeila Ghaffari // Command line Options 412a515125bSLeila Ghaffari // ------------------------------------------------------ 4131485969bSJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for ADVECTION problem", NULL); 414a515125bSLeila Ghaffari // -- Physics 4152b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL)); 416a515125bSLeila Ghaffari PetscBool translation; 4175f636aeaSJames Wright PetscCall(PetscOptionsEnum("-wind_type", "Wind type in Advection", NULL, AdvDifWindTypes, (PetscEnum)(wind_type = ADVDIF_WIND_ROTATION), 4185f636aeaSJames Wright (PetscEnum *)&wind_type, &translation)); 41966d54740SJames Wright PetscInt n = dim; 420a515125bSLeila Ghaffari PetscBool user_wind; 4212b916ea7SJeremy L Thompson PetscCall(PetscOptionsRealArray("-wind_translation", "Constant wind vector", NULL, wind, &n, &user_wind)); 422c8d249deSJames Wright PetscCall(PetscOptionsScalar("-diffusion_coeff", "Diffusion coefficient", NULL, diffusion_coeff, &diffusion_coeff, NULL)); 4232b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-CtauS", "Scale coefficient for tau (nondimensional)", NULL, CtauS, &CtauS, NULL)); 424059d1c40SJames Wright PetscCall(PetscOptionsBool("-strong_form", "Strong (true) or weak/integrated by parts (false) advection residual", NULL, strong_form, &strong_form, 425059d1c40SJames Wright NULL)); 4262b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-E_wind", "Total energy of inflow wind", NULL, E_wind, &E_wind, NULL)); 4275f636aeaSJames Wright PetscCall(PetscOptionsEnum("-advection_ic_type", "Initial condition for Advection problem", NULL, AdvDifICTypes, 4285f636aeaSJames Wright (PetscEnum)(advectionic_type = ADVDIF_IC_BUBBLE_SPHERE), (PetscEnum *)&advectionic_type, NULL)); 4292b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 43057272ee0SJames Wright PetscCall(PetscOptionsEnum("-stab_tau", "Stabilization constant, tau", NULL, StabilizationTauTypes, (PetscEnum)(stab_tau = STAB_TAU_CTAU), 43157272ee0SJames Wright (PetscEnum *)&stab_tau, NULL)); 43257272ee0SJames Wright PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 433fbabb365SJames Wright PetscCall(PetscOptionsScalar("-Ctau_a", "Coefficient for the stabilization, advection component", NULL, Ctau_a, &Ctau_a, NULL)); 434fbabb365SJames Wright PetscCall(PetscOptionsScalar("-Ctau_d", "Coefficient for the stabilization, diffusion component", NULL, Ctau_d, &Ctau_d, NULL)); 4352b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 436a515125bSLeila Ghaffari 437108c6689SJames Wright if (advectionic_type == ADVDIF_IC_WAVE) { 438108c6689SJames Wright PetscCall(PetscOptionsEnum("-wave_type", "Type of wave", NULL, AdvDifWaveTypes, (PetscEnum)(wave_type = ADVDIF_WAVE_SINE), 439108c6689SJames Wright (PetscEnum *)&wave_type, NULL)); 440a62be6baSJames Wright PetscCall(PetscOptionsScalar("-wave_frequency", "Frequency of sine wave", NULL, wave_frequency, &wave_frequency, NULL)); 441a62be6baSJames Wright PetscCall(PetscOptionsScalar("-wave_phase", "Length correction", NULL, wave_phase, &wave_phase, NULL)); 442108c6689SJames Wright } 443108c6689SJames Wright 444108c6689SJames Wright if (advectionic_type == ADVDIF_IC_BUBBLE_CYLINDER || advectionic_type == ADVDIF_IC_BUBBLE_SPHERE) { 445108c6689SJames Wright bubble_continuity_type = dim == 3 ? ADVDIF_BUBBLE_CONTINUITY_SMOOTH : ADVDIF_BUBBLE_CONTINUITY_COSINE; 446108c6689SJames Wright PetscCall(PetscOptionsEnum("-bubble_continuity", "Smooth, back_sharp, or thick", NULL, AdvDifBubbleContinuityTypes, 447108c6689SJames Wright (PetscEnum)bubble_continuity_type, (PetscEnum *)&bubble_continuity_type, NULL)); 448108c6689SJames Wright } 449a62be6baSJames Wright 450a515125bSLeila Ghaffari // -- Units 4512b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 452a515125bSLeila Ghaffari meter = fabs(meter); 4532b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL)); 454a515125bSLeila Ghaffari kilogram = fabs(kilogram); 4552b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 456a515125bSLeila Ghaffari second = fabs(second); 457a515125bSLeila Ghaffari 458a515125bSLeila Ghaffari // -- Warnings 4595f636aeaSJames Wright if (wind_type == ADVDIF_WIND_ROTATION && user_wind) { 4602b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Use -wind_translation only with -wind_type translation\n")); 461a515125bSLeila Ghaffari } 4625f636aeaSJames Wright if (wind_type == ADVDIF_WIND_TRANSLATION && advectionic_type == ADVDIF_IC_BUBBLE_CYLINDER && wind[2] != 0.) { 463a515125bSLeila Ghaffari wind[2] = 0; 464c51f031aSJames Wright PetscCall( 465c51f031aSJames Wright PetscPrintf(comm, "Warning! Background wind in the z direction should be zero (-wind_translation x,x,0) with -advection_ic_type cylinder\n")); 466a515125bSLeila Ghaffari } 467a515125bSLeila Ghaffari if (stab == STAB_NONE && CtauS != 0) { 4682b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Use -CtauS only with -stab su or -stab supg\n")); 469a515125bSLeila Ghaffari } 4701485969bSJeremy L Thompson PetscOptionsEnd(); 471a515125bSLeila Ghaffari 472a78efa86SJames Wright if (stab == STAB_SUPG) problem->create_mass_operator = CreateKSPMassOperator_AdvectionStabilized; 473a78efa86SJames Wright 474a515125bSLeila Ghaffari // ------------------------------------------------------ 475a515125bSLeila Ghaffari // Set up the PETSc context 476a515125bSLeila Ghaffari // ------------------------------------------------------ 477a515125bSLeila Ghaffari // -- Define derived units 478a515125bSLeila Ghaffari Joule = kilogram * PetscSqr(meter) / PetscSqr(second); 479a515125bSLeila Ghaffari 4800c373b74SJames Wright honee->units->meter = meter; 4810c373b74SJames Wright honee->units->kilogram = kilogram; 4820c373b74SJames Wright honee->units->second = second; 4830c373b74SJames Wright honee->units->Joule = Joule; 484a515125bSLeila Ghaffari 485a515125bSLeila Ghaffari // ------------------------------------------------------ 486a515125bSLeila Ghaffari // Set up the libCEED context 487a515125bSLeila Ghaffari // ------------------------------------------------------ 488a515125bSLeila Ghaffari // -- Scale variables to desired units 489a515125bSLeila Ghaffari E_wind *= Joule; 490a515125bSLeila Ghaffari rc = fabs(rc) * meter; 49166d54740SJames Wright for (PetscInt i = 0; i < dim; i++) { 49205a512bdSLeila Ghaffari wind[i] *= (meter / second); 49305a512bdSLeila Ghaffari domain_size[i] *= meter; 49405a512bdSLeila Ghaffari } 495a515125bSLeila Ghaffari 496a515125bSLeila Ghaffari // -- Setup Context 497a515125bSLeila Ghaffari setup_context->rc = rc; 49805a512bdSLeila Ghaffari setup_context->lx = domain_size[0]; 49905a512bdSLeila Ghaffari setup_context->ly = domain_size[1]; 50066d54740SJames Wright setup_context->lz = dim == 3 ? domain_size[2] : 0.; 501a515125bSLeila Ghaffari setup_context->wind[0] = wind[0]; 502a515125bSLeila Ghaffari setup_context->wind[1] = wind[1]; 50366d54740SJames Wright setup_context->wind[2] = dim == 3 ? wind[2] : 0.; 504a515125bSLeila Ghaffari setup_context->wind_type = wind_type; 505c51f031aSJames Wright setup_context->initial_condition_type = advectionic_type; 506a515125bSLeila Ghaffari setup_context->bubble_continuity_type = bubble_continuity_type; 507a515125bSLeila Ghaffari setup_context->time = 0; 508a62be6baSJames Wright setup_context->wave_frequency = wave_frequency; 509a62be6baSJames Wright setup_context->wave_phase = wave_phase; 510a62be6baSJames Wright setup_context->wave_type = wave_type; 511a515125bSLeila Ghaffari 512a515125bSLeila Ghaffari // -- QFunction Context 5130c373b74SJames Wright honee->phys->implicit = implicit; 51415a3537eSJed Brown advection_ctx->CtauS = CtauS; 51515a3537eSJed Brown advection_ctx->E_wind = E_wind; 51615a3537eSJed Brown advection_ctx->implicit = implicit; 51715a3537eSJed Brown advection_ctx->strong_form = strong_form; 51815a3537eSJed Brown advection_ctx->stabilization = stab; 51957272ee0SJames Wright advection_ctx->stabilization_tau = stab_tau; 52057272ee0SJames Wright advection_ctx->Ctau_a = Ctau_a; 521fbabb365SJames Wright advection_ctx->Ctau_d = Ctau_d; 52257272ee0SJames Wright advection_ctx->Ctau_t = Ctau_t; 523c8d249deSJames Wright advection_ctx->diffusion_coeff = diffusion_coeff; 5240c373b74SJames Wright advection_ctx->divFdiff_method = honee->app_ctx->divFdiffproj_method; 525a515125bSLeila Ghaffari 5260c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &problem->ics.qfctx)); 527e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(problem->ics.qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 528e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfctx, CEED_MEM_HOST, FreeContextPetsc)); 529a515125bSLeila Ghaffari 5300c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &advection_qfctx)); 531e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(advection_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*advection_ctx), advection_ctx)); 532e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(advection_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 533e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(advection_qfctx, "timestep size", offsetof(struct AdvectionContext_, dt), 1, 53457272ee0SJames Wright "Size of timestep, delta t")); 535e07531f7SJames Wright problem->apply_vol_rhs.qfctx = advection_qfctx; 536e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &problem->apply_vol_ifunction.qfctx)); 537e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &problem->apply_inflow.qfctx)); 538d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 539a515125bSLeila Ghaffari } 540a515125bSLeila Ghaffari 5410c373b74SJames Wright PetscErrorCode PRINT_ADVECTION(Honee honee, ProblemData problem, AppCtx app_ctx) { 5420c373b74SJames Wright MPI_Comm comm = honee->comm; 5430c373b74SJames Wright Ceed ceed = honee->ceed; 5443636f6a4SJames Wright SetupContextAdv setup_ctx; 54515a3537eSJed Brown AdvectionContext advection_ctx; 54666d54740SJames Wright PetscInt dim; 547a515125bSLeila Ghaffari 54815a3537eSJed Brown PetscFunctionBeginUser; 5490c373b74SJames Wright PetscCall(DMGetDimension(honee->dm, &dim)); 550e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfctx, CEED_MEM_HOST, &setup_ctx)); 551e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &advection_ctx)); 5522b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, 553a515125bSLeila Ghaffari " Problem:\n" 554a515125bSLeila Ghaffari " Problem Name : %s\n" 555a515125bSLeila Ghaffari " Stabilization : %s\n" 5565bcc2007SJames Wright " Stabilization Tau : %s\n" 557a515125bSLeila Ghaffari " Wind Type : %s\n", 5585bcc2007SJames Wright app_ctx->problem_name, StabilizationTypes[advection_ctx->stabilization], 5595f636aeaSJames Wright StabilizationTauTypes[advection_ctx->stabilization_tau], AdvDifWindTypes[setup_ctx->wind_type])); 560a515125bSLeila Ghaffari 5615f636aeaSJames Wright if (setup_ctx->wind_type == ADVDIF_WIND_TRANSLATION) { 562c781aee8SJames Wright CeedScalar *wind = setup_ctx->wind; 56366d54740SJames Wright switch (dim) { 5649529d636SJames Wright case 2: 565c781aee8SJames Wright PetscCall(PetscPrintf(comm, " Background Wind : %f,%f\n", wind[0], wind[1])); 5669529d636SJames Wright break; 5679529d636SJames Wright case 3: 568c781aee8SJames Wright PetscCall(PetscPrintf(comm, " Background Wind : %f,%f,%f\n", wind[0], wind[1], wind[2])); 5699529d636SJames Wright break; 5709529d636SJames Wright } 571a515125bSLeila Ghaffari } 572c781aee8SJames Wright 5735f636aeaSJames Wright PetscCall(PetscPrintf(comm, " Initial Condition Type : %s\n", AdvDifICTypes[setup_ctx->initial_condition_type])); 574c781aee8SJames Wright switch (setup_ctx->initial_condition_type) { 5755f636aeaSJames Wright case ADVDIF_IC_SKEW: 5765f636aeaSJames Wright case ADVDIF_IC_COSINE_HILL: 5773d1afcc1SJames Wright case ADVDIF_IC_BOUNDARY_LAYER: 578c781aee8SJames Wright break; 5795f636aeaSJames Wright case ADVDIF_IC_BUBBLE_SPHERE: 5805f636aeaSJames Wright case ADVDIF_IC_BUBBLE_CYLINDER: 5815f636aeaSJames Wright PetscCall(PetscPrintf(comm, " Bubble Continuity : %s\n", AdvDifBubbleContinuityTypes[setup_ctx->bubble_continuity_type])); 582c781aee8SJames Wright break; 5835f636aeaSJames Wright case ADVDIF_IC_WAVE: 584c781aee8SJames Wright PetscCall(PetscPrintf(comm, " Wave Type : %s\n", AdvDifWaveTypes[setup_ctx->wave_type])); 585c781aee8SJames Wright break; 586c781aee8SJames Wright } 587c781aee8SJames Wright 588e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfctx, &setup_ctx)); 589e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &advection_ctx)); 590d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 591a515125bSLeila Ghaffari } 592