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 14*e5a8cae0SJames Wright const char *const AdvDifWindTypes[] = {"ROTATION", "TRANSLATION", "BOUNDARY_LAYER", "WindType", "ADVDIF_WIND_", NULL}; 15*e5a8cae0SJames Wright const char *const AdvDifICTypes[] = {"SPHERE", "CYLINDER", "COSINE_HILL", "SKEW", "WAVE", "BOUNDARY_LAYER", "AdvectionICType", "ADVDIF_IC_", NULL}; 16*e5a8cae0SJames Wright const char *const AdvDifWaveTypes[] = {"SINE", "SQUARE", "AdvDifWaveType", "ADVDIF_WAVE_", NULL}; 17*e5a8cae0SJames Wright const char *const AdvDifBubbleContinuityTypes[] = {"SMOOTH", "BACK_SHARP", "THICK", "COSINE", "BubbleContinuityType", "ADVDIF_BUBBLE_CONTINUITY_", 18*e5a8cae0SJames Wright NULL}; 19*e5a8cae0SJames Wright const char *const StabilizationTauTypes[] = {"CTAU", "ADVDIFF_SHAKIB", "StabilizationTauType", "STAB_TAU_", NULL}; 20*e5a8cae0SJames Wright 21a78efa86SJames Wright // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping 22a78efa86SJames Wright // 23a78efa86SJames Wright // Only used for SUPG stabilization 240c373b74SJames Wright PetscErrorCode CreateKSPMassOperator_AdvectionStabilized(Honee honee, CeedOperator *op_mass) { 250c373b74SJames Wright Ceed ceed = honee->ceed; 26a78efa86SJames Wright CeedInt num_comp_q, q_data_size; 2787d3884fSJeremy L Thompson CeedQFunction qf_mass = NULL; 2801e19bfaSJames Wright CeedElemRestriction elem_restr_q, elem_restr_qd; 29a78efa86SJames Wright CeedBasis basis_q; 30a78efa86SJames Wright CeedVector q_data; 31236a8ba6SJames Wright CeedQFunctionContext qfctx = NULL; 32a78efa86SJames Wright PetscInt dim; 33a78efa86SJames Wright 34a78efa86SJames Wright PetscFunctionBeginUser; 350c373b74SJames Wright PetscCall(DMGetDimension(honee->dm, &dim)); 36a78efa86SJames Wright { // Get restriction and basis from the RHS function 37a78efa86SJames Wright CeedOperator *sub_ops; 3801e19bfaSJames Wright CeedOperatorField op_field; 39a78efa86SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 40a78efa86SJames Wright 410c373b74SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_rhs_ctx->op, &sub_ops)); 4201e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &op_field)); 4301e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_q, &basis_q, NULL)); 4401e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &op_field)); 4501e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_qd, NULL, &q_data)); 46a78efa86SJames Wright 47236a8ba6SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qfctx)); 48a78efa86SJames Wright } 49a78efa86SJames Wright 50a78efa86SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q)); 5101e19bfaSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd, &q_data_size)); 52a78efa86SJames Wright 53a78efa86SJames Wright switch (dim) { 54a78efa86SJames Wright case 2: 55a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection2D, MassFunction_Advection2D_loc, &qf_mass)); 56a78efa86SJames Wright break; 57a78efa86SJames Wright case 3: 58a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection, MassFunction_Advection_loc, &qf_mass)); 59a78efa86SJames Wright break; 60a78efa86SJames Wright } 61a78efa86SJames Wright 62236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qfctx)); 63a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0)); 64a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP)); 65a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP)); 66a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE)); 67a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP)); 68a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD)); 69a78efa86SJames Wright 70a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass)); 7171f2ed29SJames Wright PetscCallCeed(ceed, CeedOperatorSetName(*op_mass, "RHS Mass Operator, Advection-Diffusion Stabilized")); 72a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 730c373b74SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, honee->q_ceed)); 7401e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 75a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 76a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 77a78efa86SJames Wright 78fff85bd3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 7901e19bfaSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 80fff85bd3SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 81fff85bd3SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 82236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&qfctx)); 83a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 84a78efa86SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 85a78efa86SJames Wright } 86a78efa86SJames Wright 87c2d90829SJames Wright /** 88c2d90829SJames Wright @brief Create RHS CeedOperator for direct projection of divergence of diffusive flux 89c2d90829SJames Wright 900c373b74SJames Wright @param[in] honee `Honee` context 91c2d90829SJames Wright @param[in] diff_flux_proj `DivDiffFluxProjectionData` object 92c2d90829SJames Wright @param[out] op_rhs Operator to calculate the RHS of the L^2 projection 93c2d90829SJames Wright **/ 94e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Direct_AdvDif(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) { 950c373b74SJames Wright Ceed ceed = honee->ceed; 96c2d90829SJames Wright NodalProjectionData projection = diff_flux_proj->projection; 97c2d90829SJames Wright CeedInt num_comp_q; 98c2d90829SJames Wright PetscInt dim, label_value = 0; 99c2d90829SJames Wright DMLabel domain_label = NULL; 100c2d90829SJames Wright CeedQFunctionContext advection_qfctx = NULL; 101c2d90829SJames Wright 102c2d90829SJames Wright PetscFunctionBeginUser; 103c2d90829SJames Wright // -- Get Pre-requisite things 104c2d90829SJames Wright PetscCall(DMGetDimension(projection->dm, &dim)); 105e3663b90SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q)); 106c2d90829SJames Wright 10740b78511SJames Wright { // Get advection-diffusion QF context 108c2d90829SJames Wright CeedOperator *sub_ops; 109c2d90829SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 110c2d90829SJames Wright 1110c373b74SJames Wright if (honee->op_ifunction) PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_ifunction, &sub_ops)); 1120c373b74SJames Wright else PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_rhs_ctx->op, &sub_ops)); 113c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &advection_qfctx)); 114c2d90829SJames Wright } 115c2d90829SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, op_rhs)); 116c2d90829SJames Wright { // Add the volume integral CeedOperator 117c2d90829SJames Wright CeedQFunction qf_rhs_volume = NULL; 118c2d90829SJames Wright CeedOperator op_rhs_volume; 119c2d90829SJames Wright CeedVector q_data; 120c2d90829SJames Wright CeedElemRestriction elem_restr_qd, elem_restr_diff_flux_volume = NULL; 121c2d90829SJames Wright CeedBasis basis_diff_flux = NULL; 122c2d90829SJames Wright CeedInt q_data_size; 123c2d90829SJames Wright 124c2d90829SJames Wright PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_diff_flux_volume, &basis_diff_flux, NULL, NULL)); 125e3663b90SJames 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, 126e3663b90SJames Wright &q_data_size)); 127c2d90829SJames Wright switch (dim) { 128c2d90829SJames Wright case 2: 129c2d90829SJames Wright PetscCallCeed( 130c2d90829SJames Wright ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_AdvDif_2D, DivDiffusiveFluxVolumeRHS_AdvDif_2D_loc, &qf_rhs_volume)); 131c2d90829SJames Wright break; 132c2d90829SJames Wright case 3: 133c2d90829SJames Wright PetscCallCeed( 134c2d90829SJames Wright ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_AdvDif_3D, DivDiffusiveFluxVolumeRHS_AdvDif_3D_loc, &qf_rhs_volume)); 135c2d90829SJames Wright break; 136c2d90829SJames Wright } 1370c373b74SJames Wright PetscCheck(qf_rhs_volume, honee->comm, PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim); 138c2d90829SJames Wright 139c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_volume, advection_qfctx)); 140c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 141c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "qdata", q_data_size, CEED_EVAL_NONE)); 142c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_volume, "diffusive flux RHS", projection->num_comp * dim, CEED_EVAL_GRAD)); 143c2d90829SJames Wright 144c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_volume, NULL, NULL, &op_rhs_volume)); 145e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 146c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 147c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "diffusive flux RHS", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE)); 148c2d90829SJames Wright 149c2d90829SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_volume)); 150c2d90829SJames Wright 151c2d90829SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 152c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 153c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_volume)); 154c2d90829SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 155c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_volume)); 156c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_volume)); 157c2d90829SJames Wright } 158c2d90829SJames Wright 159c2d90829SJames Wright { // Add the boundary integral CeedOperator 160c2d90829SJames Wright CeedQFunction qf_rhs_boundary; 161c2d90829SJames Wright DMLabel face_sets_label; 162c2d90829SJames Wright PetscInt num_face_set_values, *face_set_values; 163c2d90829SJames Wright CeedInt q_data_size; 164c2d90829SJames Wright 165c2d90829SJames Wright // -- Build RHS operator 166c2d90829SJames Wright switch (dim) { 167c2d90829SJames Wright case 2: 168c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_AdvDif_2D, DivDiffusiveFluxBoundaryRHS_AdvDif_2D_loc, 169c2d90829SJames Wright &qf_rhs_boundary)); 170c2d90829SJames Wright break; 171c2d90829SJames Wright case 3: 172c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_AdvDif_3D, DivDiffusiveFluxBoundaryRHS_AdvDif_3D_loc, 173c2d90829SJames Wright &qf_rhs_boundary)); 174c2d90829SJames Wright break; 175c2d90829SJames Wright } 176c2d90829SJames Wright 1770c373b74SJames Wright PetscCall(QDataBoundaryGradientGetNumComponents(honee->dm, &q_data_size)); 178c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_boundary, advection_qfctx)); 179c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 180c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "qdata", q_data_size, CEED_EVAL_NONE)); 181c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_boundary, "diffusive flux RHS", projection->num_comp, CEED_EVAL_INTERP)); 182c2d90829SJames Wright 183c2d90829SJames Wright PetscCall(DMGetLabel(projection->dm, "Face Sets", &face_sets_label)); 184c2d90829SJames Wright PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_sets_label, &num_face_set_values, &face_set_values)); 185c2d90829SJames Wright for (PetscInt f = 0; f < num_face_set_values; f++) { 186c2d90829SJames Wright DMLabel face_orientation_label; 187c2d90829SJames Wright PetscInt num_orientations_values, *orientation_values; 188c2d90829SJames Wright 189c2d90829SJames Wright { 190c2d90829SJames Wright char *face_orientation_label_name; 191c2d90829SJames Wright 192c2d90829SJames Wright PetscCall(DMPlexCreateFaceLabel(projection->dm, face_set_values[f], &face_orientation_label_name)); 193c2d90829SJames Wright PetscCall(DMGetLabel(projection->dm, face_orientation_label_name, &face_orientation_label)); 194c2d90829SJames Wright PetscCall(PetscFree(face_orientation_label_name)); 195c2d90829SJames Wright } 196c2d90829SJames Wright PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_orientation_label, &num_orientations_values, &orientation_values)); 197c2d90829SJames Wright for (PetscInt o = 0; o < num_orientations_values; o++) { 198c2d90829SJames Wright CeedOperator op_rhs_boundary; 199c2d90829SJames Wright CeedBasis basis_q, basis_diff_flux_boundary; 200c2d90829SJames Wright CeedElemRestriction elem_restr_qdata, elem_restr_q, elem_restr_diff_flux_boundary; 201c2d90829SJames Wright CeedVector q_data; 202c2d90829SJames Wright CeedInt q_data_size; 203c2d90829SJames Wright PetscInt orientation = orientation_values[o], dm_field_q = 0, height_cell = 0, height_face = 1; 204c2d90829SJames Wright 2050c373b74SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, face_orientation_label, orientation, height_cell, dm_field_q, &elem_restr_q)); 2060c373b74SJames Wright PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, honee->dm, face_orientation_label, orientation, orientation, dm_field_q, &basis_q)); 207c2d90829SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, 208c2d90829SJames Wright &elem_restr_diff_flux_boundary)); 209c2d90829SJames Wright PetscCall(CreateBasisFromPlex(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, &basis_diff_flux_boundary)); 210e3663b90SJames Wright PetscCall( 211e3663b90SJames Wright QDataBoundaryGradientGet(ceed, honee->dm, face_orientation_label, orientation, honee->x_coord, &elem_restr_qdata, &q_data, &q_data_size)); 212c2d90829SJames Wright 213c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_boundary, NULL, NULL, &op_rhs_boundary)); 214c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 215c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "qdata", elem_restr_qdata, CEED_BASIS_NONE, q_data)); 216c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "diffusive flux RHS", elem_restr_diff_flux_boundary, basis_diff_flux_boundary, 217c2d90829SJames Wright CEED_VECTOR_ACTIVE)); 218c2d90829SJames Wright 219c2d90829SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_boundary)); 220c2d90829SJames Wright 221c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_boundary)); 222c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qdata)); 223c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 224c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_boundary)); 225c2d90829SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 226c2d90829SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux_boundary)); 227c2d90829SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 228c2d90829SJames Wright } 229c2d90829SJames Wright PetscCall(PetscFree(orientation_values)); 230c2d90829SJames Wright } 231c2d90829SJames Wright PetscCall(PetscFree(face_set_values)); 232c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_boundary)); 233c2d90829SJames Wright } 234c2d90829SJames Wright 235c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&advection_qfctx)); 236c2d90829SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 237c2d90829SJames Wright } 238c2d90829SJames Wright 23940b78511SJames Wright /** 24040b78511SJames Wright @brief Create RHS CeedOperator for indirect projection of divergence of diffusive flux 24140b78511SJames Wright 2420c373b74SJames Wright @param[in] honee `Honee` context 24340b78511SJames Wright @param[in] diff_flux_proj `DivDiffFluxProjectionData` object 24440b78511SJames Wright @param[out] op_rhs Operator to calculate the RHS of the L^2 projection 24540b78511SJames Wright **/ 246e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Indirect_AdvDif(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) { 2470c373b74SJames Wright Ceed ceed = honee->ceed; 24840b78511SJames Wright NodalProjectionData projection = diff_flux_proj->projection; 24940b78511SJames Wright CeedBasis basis_diff_flux; 25040b78511SJames Wright CeedElemRestriction elem_restr_diff_flux, elem_restr_qd; 25140b78511SJames Wright CeedVector q_data; 25240b78511SJames Wright CeedInt num_comp_q, q_data_size; 25340b78511SJames Wright PetscInt dim; 25440b78511SJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 25540b78511SJames Wright DMLabel domain_label = NULL; 25640b78511SJames Wright CeedQFunction qf_rhs = NULL; 25740b78511SJames Wright CeedQFunctionContext advection_qfctx = NULL; 25840b78511SJames Wright 25940b78511SJames Wright PetscFunctionBeginUser; 26040b78511SJames Wright PetscCall(DMGetDimension(projection->dm, &dim)); 261e3663b90SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q)); 26240b78511SJames Wright 26340b78511SJames Wright { // Get advection-diffusion QF context 26440b78511SJames Wright CeedOperator *sub_ops; 26540b78511SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 26640b78511SJames Wright 2670c373b74SJames Wright if (honee->op_ifunction) PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_ifunction, &sub_ops)); 2680c373b74SJames Wright else PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_rhs_ctx->op, &sub_ops)); 26940b78511SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &advection_qfctx)); 27040b78511SJames Wright } 27140b78511SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, domain_label, label_value, height, dm_field, &elem_restr_diff_flux)); 27240b78511SJames Wright PetscCall(CreateBasisFromPlex(ceed, projection->dm, domain_label, label_value, height, dm_field, &basis_diff_flux)); 273e3663b90SJames 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, 274e3663b90SJames Wright &q_data_size)); 27540b78511SJames Wright 27640b78511SJames Wright switch (dim) { 27740b78511SJames Wright case 2: 27840b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_AdvDif_2D, DiffusiveFluxRHS_AdvDif_2D_loc, &qf_rhs)); 27940b78511SJames Wright break; 28040b78511SJames Wright case 3: 28140b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_AdvDif_3D, DiffusiveFluxRHS_AdvDif_3D_loc, &qf_rhs)); 28240b78511SJames Wright break; 28340b78511SJames Wright } 2840c373b74SJames Wright PetscCheck(qf_rhs, honee->comm, PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim); 28540b78511SJames Wright 28640b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs, advection_qfctx)); 28740b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 28840b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "qdata", q_data_size, CEED_EVAL_NONE)); 28940b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs, "F_diff RHS", projection->num_comp, CEED_EVAL_INTERP)); 29040b78511SJames Wright 29140b78511SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, op_rhs)); 292e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 29340b78511SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 29440b78511SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "F_diff RHS", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE)); 29540b78511SJames Wright 29640b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs)); 29740b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&advection_qfctx)); 29840b78511SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 29940b78511SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 30040b78511SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 30140b78511SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux)); 30240b78511SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 30340b78511SJames Wright } 30440b78511SJames Wright 305d6cac220SJames Wright static PetscErrorCode AdvectionInflowBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) { 306d6cac220SJames Wright HoneeBCStruct honee_bc; 307d6cac220SJames Wright DM dm; 308d6cac220SJames Wright PetscInt dim; 309d6cac220SJames Wright 310d6cac220SJames Wright PetscFunctionBeginUser; 311d6cac220SJames Wright PetscCall(BCDefinitionGetContext(bc_def, &honee_bc)); 312d6cac220SJames Wright PetscCall(BCDefinitionGetDM(bc_def, &dm)); 313d6cac220SJames Wright PetscCall(DMGetDimension(dm, &dim)); 314d6cac220SJames Wright switch (dim) { 315d6cac220SJames Wright case 2: 316d6cac220SJames Wright PetscCall(HoneeBCCreateIFunctionQF(bc_def, Advection2d_InOutFlow, Advection2d_InOutFlow_loc, honee_bc->qfctx, qf)); 317d6cac220SJames Wright break; 318d6cac220SJames Wright case 3: 319d6cac220SJames Wright PetscCall(HoneeBCCreateIFunctionQF(bc_def, Advection_InOutFlow, Advection_InOutFlow_loc, honee_bc->qfctx, qf)); 320d6cac220SJames Wright break; 321d6cac220SJames Wright } 322d6cac220SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 323d6cac220SJames Wright } 324d6cac220SJames Wright 325d3c60affSJames Wright PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx) { 3265f636aeaSJames Wright AdvDifWindType wind_type; 3275f636aeaSJames Wright AdvDifICType advectionic_type; 328108c6689SJames Wright AdvDifBubbleContinuityType bubble_continuity_type = -1; 329a515125bSLeila Ghaffari StabilizationType stab; 33057272ee0SJames Wright StabilizationTauType stab_tau; 3313636f6a4SJames Wright SetupContextAdv setup_context; 3320c373b74SJames Wright Honee honee = *(Honee *)ctx; 3330c373b74SJames Wright MPI_Comm comm = honee->comm; 3340c373b74SJames Wright Ceed ceed = honee->ceed; 335a515125bSLeila Ghaffari PetscBool implicit; 33615a3537eSJed Brown AdvectionContext advection_ctx; 337e07531f7SJames Wright CeedQFunctionContext advection_qfctx; 3389529d636SJames Wright PetscInt dim; 339a515125bSLeila Ghaffari 34015a3537eSJed Brown PetscFunctionBeginUser; 3412b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &setup_context)); 3422b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &advection_ctx)); 3439529d636SJames Wright PetscCall(DMGetDimension(dm, &dim)); 344a515125bSLeila Ghaffari 345a515125bSLeila Ghaffari // ------------------------------------------------------ 346a515125bSLeila Ghaffari // SET UP ADVECTION 347a515125bSLeila Ghaffari // ------------------------------------------------------ 34828160fc2SJames Wright problem->print_info = PRINT_ADVECTION; 3491abc2837SJames Wright problem->num_comps_jac_data = 0; 3509529d636SJames Wright switch (dim) { 3519529d636SJames Wright case 2: 352e07531f7SJames Wright problem->ics.qf_func_ptr = ICsAdvection2d; 353e07531f7SJames Wright problem->ics.qf_loc = ICsAdvection2d_loc; 354e07531f7SJames Wright problem->apply_vol_rhs.qf_func_ptr = RHS_Advection2d; 355e07531f7SJames Wright problem->apply_vol_rhs.qf_loc = RHS_Advection2d_loc; 356e07531f7SJames Wright problem->apply_vol_ifunction.qf_func_ptr = IFunction_Advection2d; 357e07531f7SJames Wright problem->apply_vol_ifunction.qf_loc = IFunction_Advection2d_loc; 35858ce1233SJames Wright problem->compute_exact_solution_error = PETSC_TRUE; 3599529d636SJames Wright break; 3609529d636SJames Wright case 3: 361e07531f7SJames Wright problem->ics.qf_func_ptr = ICsAdvection; 362e07531f7SJames Wright problem->ics.qf_loc = ICsAdvection_loc; 363e07531f7SJames Wright problem->apply_vol_rhs.qf_func_ptr = RHS_Advection; 364e07531f7SJames Wright problem->apply_vol_rhs.qf_loc = RHS_Advection_loc; 365e07531f7SJames Wright problem->apply_vol_ifunction.qf_func_ptr = IFunction_Advection; 366e07531f7SJames Wright problem->apply_vol_ifunction.qf_loc = IFunction_Advection_loc; 36758ce1233SJames Wright problem->compute_exact_solution_error = PETSC_FALSE; 3689529d636SJames Wright break; 3699529d636SJames Wright } 370a515125bSLeila Ghaffari 3710c373b74SJames Wright PetscCall(DivDiffFluxProjectionCreate(honee, 1, &honee->diff_flux_proj)); 3720c373b74SJames Wright if (honee->diff_flux_proj) { 3730c373b74SJames Wright DivDiffFluxProjectionData diff_flux_proj = honee->diff_flux_proj; 374c2d90829SJames Wright NodalProjectionData projection = diff_flux_proj->projection; 375c2d90829SJames Wright 376c2d90829SJames Wright diff_flux_proj->CreateRHSOperator_Direct = DivDiffFluxProjectionCreateRHS_Direct_AdvDif; 37740b78511SJames Wright diff_flux_proj->CreateRHSOperator_Indirect = DivDiffFluxProjectionCreateRHS_Indirect_AdvDif; 378c2d90829SJames Wright 3790c373b74SJames Wright switch (honee->diff_flux_proj->method) { 380c2d90829SJames Wright case DIV_DIFF_FLUX_PROJ_DIRECT: { 381c2d90829SJames Wright PetscSection section; 382c2d90829SJames Wright 383c2d90829SJames Wright PetscCall(DMGetLocalSection(projection->dm, §ion)); 384c2d90829SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 385c2d90829SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "DivDiffusiveFlux_Scalar")); 386c2d90829SJames Wright } break; 387c2d90829SJames Wright case DIV_DIFF_FLUX_PROJ_INDIRECT: { 388c2d90829SJames Wright PetscSection section; 389c2d90829SJames Wright 390c2d90829SJames Wright PetscCall(DMGetLocalSection(projection->dm, §ion)); 391c2d90829SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 392c2d90829SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "DiffusiveFlux_ScalarX")); 393c2d90829SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "DiffusiveFlux_ScalarY")); 39440b78511SJames Wright if (dim >= 3) PetscCall(PetscSectionSetComponentName(section, 0, 2, "DiffusiveFlux_ScalarZ")); 395c2d90829SJames Wright } break; 396c2d90829SJames Wright case DIV_DIFF_FLUX_PROJ_NONE: 3970c373b74SJames Wright SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s", 3980c373b74SJames Wright DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]); 399c2d90829SJames Wright break; 400c2d90829SJames Wright } 401c2d90829SJames Wright } 402c2d90829SJames Wright 403a515125bSLeila Ghaffari // ------------------------------------------------------ 404ea615d4cSJames Wright // Create the QFunction context 405a515125bSLeila Ghaffari // ------------------------------------------------------ 406a515125bSLeila Ghaffari CeedScalar rc = 1000.; // m (Radius of bubble) 407a515125bSLeila Ghaffari CeedScalar CtauS = 0.; // dimensionless 408059d1c40SJames Wright PetscBool strong_form = PETSC_FALSE; 409a515125bSLeila Ghaffari CeedScalar E_wind = 1.e6; // J 4100c373b74SJames Wright CeedScalar Ctau_a = PetscPowScalarInt(honee->app_ctx->degree, 2); 4110c373b74SJames Wright CeedScalar Ctau_d = PetscPowScalarInt(honee->app_ctx->degree, 4); 41257272ee0SJames Wright CeedScalar Ctau_t = 0.; 413a515125bSLeila Ghaffari PetscReal wind[3] = {1., 0, 0}; // m/s 414c8d249deSJames Wright CeedScalar diffusion_coeff = 0.; 415a62be6baSJames Wright CeedScalar wave_frequency = 2 * M_PI; 416a62be6baSJames Wright CeedScalar wave_phase = 0; 417108c6689SJames Wright AdvDifWaveType wave_type = -1; 418b4fd18dfSJames Wright PetscScalar bl_height_factor = 1.; 41987d3884fSJeremy L Thompson PetscReal domain_min[3], domain_max[3], domain_size[3] = {0.}; 4202b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 42166d54740SJames Wright for (PetscInt i = 0; i < dim; i++) domain_size[i] = domain_max[i] - domain_min[i]; 42205a512bdSLeila Ghaffari 423a515125bSLeila Ghaffari // ------------------------------------------------------ 424a515125bSLeila Ghaffari // Create the PETSc context 425a515125bSLeila Ghaffari // ------------------------------------------------------ 426a515125bSLeila Ghaffari PetscScalar meter = 1e-2; // 1 meter in scaled length units 427a515125bSLeila Ghaffari PetscScalar kilogram = 1e-6; // 1 kilogram in scaled mass units 428a515125bSLeila Ghaffari PetscScalar second = 1e-2; // 1 second in scaled time units 429a515125bSLeila Ghaffari PetscScalar Joule; 430a515125bSLeila Ghaffari 431a515125bSLeila Ghaffari // ------------------------------------------------------ 432a515125bSLeila Ghaffari // Command line Options 433a515125bSLeila Ghaffari // ------------------------------------------------------ 4341485969bSJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for ADVECTION problem", NULL); 435a515125bSLeila Ghaffari // -- Physics 436a515125bSLeila Ghaffari PetscBool translation; 4375f636aeaSJames Wright PetscCall(PetscOptionsEnum("-wind_type", "Wind type in Advection", NULL, AdvDifWindTypes, (PetscEnum)(wind_type = ADVDIF_WIND_ROTATION), 4385f636aeaSJames Wright (PetscEnum *)&wind_type, &translation)); 43966d54740SJames Wright PetscInt n = dim; 440a515125bSLeila Ghaffari PetscBool user_wind; 4412b916ea7SJeremy L Thompson PetscCall(PetscOptionsRealArray("-wind_translation", "Constant wind vector", NULL, wind, &n, &user_wind)); 442c8d249deSJames Wright PetscCall(PetscOptionsScalar("-diffusion_coeff", "Diffusion coefficient", NULL, diffusion_coeff, &diffusion_coeff, NULL)); 4432b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-CtauS", "Scale coefficient for tau (nondimensional)", NULL, CtauS, &CtauS, NULL)); 444059d1c40SJames Wright PetscCall(PetscOptionsBool("-strong_form", "Strong (true) or weak/integrated by parts (false) advection residual", NULL, strong_form, &strong_form, 445059d1c40SJames Wright NULL)); 4462b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-E_wind", "Total energy of inflow wind", NULL, E_wind, &E_wind, NULL)); 4472b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 44857272ee0SJames Wright PetscCall(PetscOptionsEnum("-stab_tau", "Stabilization constant, tau", NULL, StabilizationTauTypes, (PetscEnum)(stab_tau = STAB_TAU_CTAU), 44957272ee0SJames Wright (PetscEnum *)&stab_tau, NULL)); 45057272ee0SJames Wright PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 451fbabb365SJames Wright PetscCall(PetscOptionsScalar("-Ctau_a", "Coefficient for the stabilization, advection component", NULL, Ctau_a, &Ctau_a, NULL)); 452fbabb365SJames Wright PetscCall(PetscOptionsScalar("-Ctau_d", "Coefficient for the stabilization, diffusion component", NULL, Ctau_d, &Ctau_d, NULL)); 4532b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 45480e9ac5bSJames Wright PetscCall(PetscOptionsEnum("-advection_ic_type", "Initial condition for Advection problem", NULL, AdvDifICTypes, 45580e9ac5bSJames Wright (PetscEnum)(advectionic_type = ADVDIF_IC_BUBBLE_SPHERE), (PetscEnum *)&advectionic_type, NULL)); 45657bb6b22SJames Wright // IC-specific options 45757bb6b22SJames Wright switch (advectionic_type) { 45857bb6b22SJames Wright case ADVDIF_IC_WAVE: 45980e9ac5bSJames Wright PetscCall(PetscOptionsDeprecated("-wave_type", "-advection_ic_wave_type", "HONEE 0.0", NULL)); 46080e9ac5bSJames Wright PetscCall(PetscOptionsDeprecated("-wave_frequency", "-advection_ic_wave_frequency", "HONEE 0.0", NULL)); 46180e9ac5bSJames Wright PetscCall(PetscOptionsDeprecated("-wave_phase", "-advection_ic_wave_phase", "HONEE 0.0", NULL)); 46280e9ac5bSJames Wright PetscCall(PetscOptionsEnum("-advection_ic_wave_type", "Type of wave", NULL, AdvDifWaveTypes, (PetscEnum)(wave_type = ADVDIF_WAVE_SINE), 463108c6689SJames Wright (PetscEnum *)&wave_type, NULL)); 46480e9ac5bSJames Wright PetscCall(PetscOptionsScalar("-advection_ic_wave_frequency", "Frequency of sine wave", NULL, wave_frequency, &wave_frequency, NULL)); 46580e9ac5bSJames Wright PetscCall(PetscOptionsScalar("-advection_ic_wave_phase", "Length correction", NULL, wave_phase, &wave_phase, NULL)); 46657bb6b22SJames Wright break; 46757bb6b22SJames Wright case ADVDIF_IC_BOUNDARY_LAYER: 46880e9ac5bSJames Wright PetscCall( 46980e9ac5bSJames Wright PetscOptionsScalar("-advection_ic_bl_height_factor", "Height of boundary layer in IC", NULL, bl_height_factor, &bl_height_factor, NULL)); 47057bb6b22SJames Wright break; 47157bb6b22SJames Wright case ADVDIF_IC_BUBBLE_CYLINDER: 47257bb6b22SJames Wright case ADVDIF_IC_BUBBLE_SPHERE: 47380e9ac5bSJames Wright PetscCall(PetscOptionsDeprecated("-rc", "-advection_ic_bubble_rc", "HONEE 0.0", NULL)); 47480e9ac5bSJames Wright PetscCall(PetscOptionsDeprecated("-bubble_continuity", "-advection_ic_bubble_continuity", "HONEE 0.0", NULL)); 47580e9ac5bSJames Wright PetscCall(PetscOptionsScalar("-advection_ic_bubble_rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL)); 476108c6689SJames Wright bubble_continuity_type = dim == 3 ? ADVDIF_BUBBLE_CONTINUITY_SMOOTH : ADVDIF_BUBBLE_CONTINUITY_COSINE; 47780e9ac5bSJames Wright PetscCall(PetscOptionsEnum("-advection_ic_bubble_continuity", "Smooth, back_sharp, or thick", NULL, AdvDifBubbleContinuityTypes, 478108c6689SJames Wright (PetscEnum)bubble_continuity_type, (PetscEnum *)&bubble_continuity_type, NULL)); 47957bb6b22SJames Wright break; 48057bb6b22SJames Wright case ADVDIF_IC_SKEW: 48157bb6b22SJames Wright case ADVDIF_IC_COSINE_HILL: 48257bb6b22SJames Wright break; 483108c6689SJames Wright } 484a62be6baSJames Wright 485a515125bSLeila Ghaffari // -- Units 4862b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 487a515125bSLeila Ghaffari meter = fabs(meter); 4882b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL)); 489a515125bSLeila Ghaffari kilogram = fabs(kilogram); 4902b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 491a515125bSLeila Ghaffari second = fabs(second); 492a515125bSLeila Ghaffari 493a515125bSLeila Ghaffari // -- Warnings 4945f636aeaSJames Wright if (wind_type == ADVDIF_WIND_ROTATION && user_wind) { 4952b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Use -wind_translation only with -wind_type translation\n")); 496a515125bSLeila Ghaffari } 4975f636aeaSJames Wright if (wind_type == ADVDIF_WIND_TRANSLATION && advectionic_type == ADVDIF_IC_BUBBLE_CYLINDER && wind[2] != 0.) { 498a515125bSLeila Ghaffari wind[2] = 0; 499c51f031aSJames Wright PetscCall( 500c51f031aSJames Wright PetscPrintf(comm, "Warning! Background wind in the z direction should be zero (-wind_translation x,x,0) with -advection_ic_type cylinder\n")); 501a515125bSLeila Ghaffari } 502a515125bSLeila Ghaffari if (stab == STAB_NONE && CtauS != 0) { 5032b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Use -CtauS only with -stab su or -stab supg\n")); 504a515125bSLeila Ghaffari } 5051485969bSJeremy L Thompson PetscOptionsEnd(); 506a515125bSLeila Ghaffari 507a78efa86SJames Wright if (stab == STAB_SUPG) problem->create_mass_operator = CreateKSPMassOperator_AdvectionStabilized; 508a78efa86SJames Wright 509a515125bSLeila Ghaffari // ------------------------------------------------------ 510a515125bSLeila Ghaffari // Set up the PETSc context 511a515125bSLeila Ghaffari // ------------------------------------------------------ 512a515125bSLeila Ghaffari // -- Define derived units 513a515125bSLeila Ghaffari Joule = kilogram * PetscSqr(meter) / PetscSqr(second); 514a515125bSLeila Ghaffari 5150c373b74SJames Wright honee->units->meter = meter; 5160c373b74SJames Wright honee->units->kilogram = kilogram; 5170c373b74SJames Wright honee->units->second = second; 5180c373b74SJames Wright honee->units->Joule = Joule; 519a515125bSLeila Ghaffari 520a515125bSLeila Ghaffari // ------------------------------------------------------ 521ea615d4cSJames Wright // Set up the QFunction contexts 522a515125bSLeila Ghaffari // ------------------------------------------------------ 523a515125bSLeila Ghaffari // -- Scale variables to desired units 524a515125bSLeila Ghaffari E_wind *= Joule; 525a515125bSLeila Ghaffari rc = fabs(rc) * meter; 52666d54740SJames Wright for (PetscInt i = 0; i < dim; i++) { 52705a512bdSLeila Ghaffari wind[i] *= (meter / second); 52805a512bdSLeila Ghaffari domain_size[i] *= meter; 52905a512bdSLeila Ghaffari } 530a515125bSLeila Ghaffari 531a515125bSLeila Ghaffari // -- Setup Context 532a515125bSLeila Ghaffari setup_context->rc = rc; 53305a512bdSLeila Ghaffari setup_context->lx = domain_size[0]; 53405a512bdSLeila Ghaffari setup_context->ly = domain_size[1]; 53566d54740SJames Wright setup_context->lz = dim == 3 ? domain_size[2] : 0.; 536a515125bSLeila Ghaffari setup_context->wind[0] = wind[0]; 537a515125bSLeila Ghaffari setup_context->wind[1] = wind[1]; 53866d54740SJames Wright setup_context->wind[2] = dim == 3 ? wind[2] : 0.; 539a515125bSLeila Ghaffari setup_context->wind_type = wind_type; 540c51f031aSJames Wright setup_context->initial_condition_type = advectionic_type; 541a515125bSLeila Ghaffari setup_context->bubble_continuity_type = bubble_continuity_type; 542a515125bSLeila Ghaffari setup_context->time = 0; 543a62be6baSJames Wright setup_context->wave_frequency = wave_frequency; 544a62be6baSJames Wright setup_context->wave_phase = wave_phase; 545a62be6baSJames Wright setup_context->wave_type = wave_type; 546b4fd18dfSJames Wright setup_context->bl_height_factor = bl_height_factor; 547a515125bSLeila Ghaffari 548a515125bSLeila Ghaffari // -- QFunction Context 5490c373b74SJames Wright honee->phys->implicit = implicit; 55015a3537eSJed Brown advection_ctx->CtauS = CtauS; 55115a3537eSJed Brown advection_ctx->E_wind = E_wind; 55215a3537eSJed Brown advection_ctx->implicit = implicit; 55315a3537eSJed Brown advection_ctx->strong_form = strong_form; 55415a3537eSJed Brown advection_ctx->stabilization = stab; 55557272ee0SJames Wright advection_ctx->stabilization_tau = stab_tau; 55657272ee0SJames Wright advection_ctx->Ctau_a = Ctau_a; 557fbabb365SJames Wright advection_ctx->Ctau_d = Ctau_d; 55857272ee0SJames Wright advection_ctx->Ctau_t = Ctau_t; 559c8d249deSJames Wright advection_ctx->diffusion_coeff = diffusion_coeff; 5600c373b74SJames Wright advection_ctx->divFdiff_method = honee->app_ctx->divFdiffproj_method; 561a515125bSLeila Ghaffari 5620c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &problem->ics.qfctx)); 563e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(problem->ics.qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 564e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfctx, CEED_MEM_HOST, FreeContextPetsc)); 565a515125bSLeila Ghaffari 5660c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &advection_qfctx)); 567e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(advection_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*advection_ctx), advection_ctx)); 568e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(advection_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 569e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(advection_qfctx, "timestep size", offsetof(struct AdvectionContext_, dt), 1, 57057272ee0SJames Wright "Size of timestep, delta t")); 571e07531f7SJames Wright problem->apply_vol_rhs.qfctx = advection_qfctx; 572e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &problem->apply_vol_ifunction.qfctx)); 573d6cac220SJames Wright 574d6cac220SJames Wright for (PetscCount b = 0; b < problem->num_bc_defs; b++) { 575d6cac220SJames Wright BCDefinition bc_def = problem->bc_defs[b]; 576d6cac220SJames Wright const char *name; 577d6cac220SJames Wright 578d6cac220SJames Wright PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL)); 579d6cac220SJames Wright if (!strcmp(name, "inflow")) { 580d6cac220SJames Wright HoneeBCStruct honee_bc; 581d6cac220SJames Wright 582d6cac220SJames Wright PetscCall(PetscNew(&honee_bc)); 583d6cac220SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &honee_bc->qfctx)); 584d6cac220SJames Wright honee_bc->honee = honee; 5851abc2837SJames Wright honee_bc->num_comps_jac_data = 0; 586d6cac220SJames Wright PetscCall(BCDefinitionSetContext(bc_def, HoneeBCDestroy, honee_bc)); 587d6cac220SJames Wright 588d6cac220SJames Wright PetscCall(BCDefinitionSetIFunction(bc_def, AdvectionInflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp)); 589d6cac220SJames Wright PetscCall(BCDefinitionSetIJacobian(bc_def, NULL, NULL)); 590d6cac220SJames Wright } 591d6cac220SJames Wright } 592d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 593a515125bSLeila Ghaffari } 594a515125bSLeila Ghaffari 5950c373b74SJames Wright PetscErrorCode PRINT_ADVECTION(Honee honee, ProblemData problem, AppCtx app_ctx) { 5960c373b74SJames Wright MPI_Comm comm = honee->comm; 5970c373b74SJames Wright Ceed ceed = honee->ceed; 5983636f6a4SJames Wright SetupContextAdv setup_ctx; 59915a3537eSJed Brown AdvectionContext advection_ctx; 60066d54740SJames Wright PetscInt dim; 601a515125bSLeila Ghaffari 60215a3537eSJed Brown PetscFunctionBeginUser; 6030c373b74SJames Wright PetscCall(DMGetDimension(honee->dm, &dim)); 604e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfctx, CEED_MEM_HOST, &setup_ctx)); 605e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &advection_ctx)); 6062b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, 607a515125bSLeila Ghaffari " Problem:\n" 608a515125bSLeila Ghaffari " Problem Name : %s\n" 609a515125bSLeila Ghaffari " Stabilization : %s\n" 6105bcc2007SJames Wright " Stabilization Tau : %s\n" 611a515125bSLeila Ghaffari " Wind Type : %s\n", 6125bcc2007SJames Wright app_ctx->problem_name, StabilizationTypes[advection_ctx->stabilization], 6135f636aeaSJames Wright StabilizationTauTypes[advection_ctx->stabilization_tau], AdvDifWindTypes[setup_ctx->wind_type])); 614a515125bSLeila Ghaffari 6155f636aeaSJames Wright if (setup_ctx->wind_type == ADVDIF_WIND_TRANSLATION) { 616c781aee8SJames Wright CeedScalar *wind = setup_ctx->wind; 61766d54740SJames Wright switch (dim) { 6189529d636SJames Wright case 2: 619c781aee8SJames Wright PetscCall(PetscPrintf(comm, " Background Wind : %f,%f\n", wind[0], wind[1])); 6209529d636SJames Wright break; 6219529d636SJames Wright case 3: 622c781aee8SJames Wright PetscCall(PetscPrintf(comm, " Background Wind : %f,%f,%f\n", wind[0], wind[1], wind[2])); 6239529d636SJames Wright break; 6249529d636SJames Wright } 625a515125bSLeila Ghaffari } 626c781aee8SJames Wright 6275f636aeaSJames Wright PetscCall(PetscPrintf(comm, " Initial Condition Type : %s\n", AdvDifICTypes[setup_ctx->initial_condition_type])); 628c781aee8SJames Wright switch (setup_ctx->initial_condition_type) { 6295f636aeaSJames Wright case ADVDIF_IC_SKEW: 6305f636aeaSJames Wright case ADVDIF_IC_COSINE_HILL: 6313d1afcc1SJames Wright case ADVDIF_IC_BOUNDARY_LAYER: 632c781aee8SJames Wright break; 6335f636aeaSJames Wright case ADVDIF_IC_BUBBLE_SPHERE: 6345f636aeaSJames Wright case ADVDIF_IC_BUBBLE_CYLINDER: 6355f636aeaSJames Wright PetscCall(PetscPrintf(comm, " Bubble Continuity : %s\n", AdvDifBubbleContinuityTypes[setup_ctx->bubble_continuity_type])); 636c781aee8SJames Wright break; 6375f636aeaSJames Wright case ADVDIF_IC_WAVE: 638c781aee8SJames Wright PetscCall(PetscPrintf(comm, " Wave Type : %s\n", AdvDifWaveTypes[setup_ctx->wave_type])); 639c781aee8SJames Wright break; 640c781aee8SJames Wright } 641c781aee8SJames Wright 642e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfctx, &setup_ctx)); 643e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &advection_ctx)); 644d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 645a515125bSLeila Ghaffari } 646