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 14e5a8cae0SJames Wright const char *const AdvDifWindTypes[] = {"ROTATION", "TRANSLATION", "BOUNDARY_LAYER", "WindType", "ADVDIF_WIND_", NULL}; 15e5a8cae0SJames Wright const char *const AdvDifICTypes[] = {"SPHERE", "CYLINDER", "COSINE_HILL", "SKEW", "WAVE", "BOUNDARY_LAYER", "AdvectionICType", "ADVDIF_IC_", NULL}; 16e5a8cae0SJames Wright const char *const AdvDifWaveTypes[] = {"SINE", "SQUARE", "AdvDifWaveType", "ADVDIF_WAVE_", NULL}; 17e5a8cae0SJames Wright const char *const AdvDifBubbleContinuityTypes[] = {"SMOOTH", "BACK_SHARP", "THICK", "COSINE", "BubbleContinuityType", "ADVDIF_BUBBLE_CONTINUITY_", 18e5a8cae0SJames Wright NULL}; 19e5a8cae0SJames Wright const char *const StabilizationTauTypes[] = {"CTAU", "ADVDIFF_SHAKIB", "StabilizationTauType", "STAB_TAU_", NULL}; 20e5a8cae0SJames Wright 219b05e62eSJames Wright static PetscErrorCode PRINT_ADVECTION(Honee honee, ProblemData problem, AppCtx app_ctx) { 229b05e62eSJames Wright MPI_Comm comm = honee->comm; 239b05e62eSJames Wright Ceed ceed = honee->ceed; 249b05e62eSJames Wright SetupContextAdv setup_ctx; 259b05e62eSJames Wright AdvectionContext advection_ctx; 269b05e62eSJames Wright PetscInt dim; 279b05e62eSJames Wright 289b05e62eSJames Wright PetscFunctionBeginUser; 299b05e62eSJames Wright PetscCall(DMGetDimension(honee->dm, &dim)); 309b05e62eSJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfctx, CEED_MEM_HOST, &setup_ctx)); 319b05e62eSJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &advection_ctx)); 329b05e62eSJames Wright PetscCall(PetscPrintf(comm, 339b05e62eSJames Wright " Problem:\n" 349b05e62eSJames Wright " Problem Name : %s\n" 359b05e62eSJames Wright " Stabilization : %s\n" 369b05e62eSJames Wright " Stabilization Tau : %s\n" 379b05e62eSJames Wright " Wind Type : %s\n", 389b05e62eSJames Wright app_ctx->problem_name, StabilizationTypes[advection_ctx->stabilization], 399b05e62eSJames Wright StabilizationTauTypes[advection_ctx->stabilization_tau], AdvDifWindTypes[setup_ctx->wind_type])); 409b05e62eSJames Wright 419b05e62eSJames Wright if (setup_ctx->wind_type == ADVDIF_WIND_TRANSLATION) { 429b05e62eSJames Wright CeedScalar *wind = setup_ctx->wind; 439b05e62eSJames Wright switch (dim) { 449b05e62eSJames Wright case 2: 459b05e62eSJames Wright PetscCall(PetscPrintf(comm, " Background Wind : %f,%f\n", wind[0], wind[1])); 469b05e62eSJames Wright break; 479b05e62eSJames Wright case 3: 489b05e62eSJames Wright PetscCall(PetscPrintf(comm, " Background Wind : %f,%f,%f\n", wind[0], wind[1], wind[2])); 499b05e62eSJames Wright break; 509b05e62eSJames Wright } 519b05e62eSJames Wright } 529b05e62eSJames Wright 539b05e62eSJames Wright PetscCall(PetscPrintf(comm, " Initial Condition Type : %s\n", AdvDifICTypes[setup_ctx->initial_condition_type])); 549b05e62eSJames Wright switch (setup_ctx->initial_condition_type) { 559b05e62eSJames Wright case ADVDIF_IC_SKEW: 569b05e62eSJames Wright case ADVDIF_IC_COSINE_HILL: 579b05e62eSJames Wright case ADVDIF_IC_BOUNDARY_LAYER: 589b05e62eSJames Wright break; 599b05e62eSJames Wright case ADVDIF_IC_BUBBLE_SPHERE: 609b05e62eSJames Wright case ADVDIF_IC_BUBBLE_CYLINDER: 619b05e62eSJames Wright PetscCall(PetscPrintf(comm, " Bubble Continuity : %s\n", AdvDifBubbleContinuityTypes[setup_ctx->bubble_continuity_type])); 629b05e62eSJames Wright break; 639b05e62eSJames Wright case ADVDIF_IC_WAVE: 649b05e62eSJames Wright PetscCall(PetscPrintf(comm, " Wave Type : %s\n", AdvDifWaveTypes[setup_ctx->wave_type])); 659b05e62eSJames Wright break; 669b05e62eSJames Wright } 679b05e62eSJames Wright 689b05e62eSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfctx, &setup_ctx)); 699b05e62eSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &advection_ctx)); 709b05e62eSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 719b05e62eSJames Wright } 729b05e62eSJames Wright 73a78efa86SJames Wright // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping 74a78efa86SJames Wright // 75a78efa86SJames Wright // Only used for SUPG stabilization 760c373b74SJames Wright PetscErrorCode CreateKSPMassOperator_AdvectionStabilized(Honee honee, CeedOperator *op_mass) { 770c373b74SJames Wright Ceed ceed = honee->ceed; 78a78efa86SJames Wright CeedInt num_comp_q, q_data_size; 7987d3884fSJeremy L Thompson CeedQFunction qf_mass = NULL; 8001e19bfaSJames Wright CeedElemRestriction elem_restr_q, elem_restr_qd; 81a78efa86SJames Wright CeedBasis basis_q; 82a78efa86SJames Wright CeedVector q_data; 83236a8ba6SJames Wright CeedQFunctionContext qfctx = NULL; 84a78efa86SJames Wright PetscInt dim; 85a78efa86SJames Wright 86a78efa86SJames Wright PetscFunctionBeginUser; 870c373b74SJames Wright PetscCall(DMGetDimension(honee->dm, &dim)); 88a78efa86SJames Wright { // Get restriction and basis from the RHS function 89a78efa86SJames Wright CeedOperator *sub_ops; 9001e19bfaSJames Wright CeedOperatorField op_field; 91a78efa86SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 92a78efa86SJames Wright 930c373b74SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_rhs_ctx->op, &sub_ops)); 9401e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &op_field)); 9501e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_q, &basis_q, NULL)); 9601e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &op_field)); 9701e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_qd, NULL, &q_data)); 98a78efa86SJames Wright 99236a8ba6SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qfctx)); 100a78efa86SJames Wright } 101a78efa86SJames Wright 102a78efa86SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q)); 10301e19bfaSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd, &q_data_size)); 104a78efa86SJames Wright 105a78efa86SJames Wright switch (dim) { 106a78efa86SJames Wright case 2: 107a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection2D, MassFunction_Advection2D_loc, &qf_mass)); 108a78efa86SJames Wright break; 109a78efa86SJames Wright case 3: 110a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection, MassFunction_Advection_loc, &qf_mass)); 111a78efa86SJames Wright break; 112a78efa86SJames Wright } 113a78efa86SJames Wright 114236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qfctx)); 115a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0)); 116a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP)); 117a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP)); 118a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE)); 119a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP)); 120a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD)); 121a78efa86SJames Wright 122a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass)); 12371f2ed29SJames Wright PetscCallCeed(ceed, CeedOperatorSetName(*op_mass, "RHS Mass Operator, Advection-Diffusion Stabilized")); 124a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 1250c373b74SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, honee->q_ceed)); 12601e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 127a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 128a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 129a78efa86SJames Wright 130fff85bd3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 13101e19bfaSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 132fff85bd3SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 133fff85bd3SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 134236a8ba6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&qfctx)); 135a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 136a78efa86SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 137a78efa86SJames Wright } 138a78efa86SJames Wright 139c2d90829SJames Wright /** 140c2d90829SJames Wright @brief Create RHS CeedOperator for direct projection of divergence of diffusive flux 141c2d90829SJames Wright 1420c373b74SJames Wright @param[in] honee `Honee` context 143c2d90829SJames Wright @param[in] diff_flux_proj `DivDiffFluxProjectionData` object 144c2d90829SJames Wright @param[out] op_rhs Operator to calculate the RHS of the L^2 projection 145c2d90829SJames Wright **/ 146e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Direct_AdvDif(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) { 1470c373b74SJames Wright Ceed ceed = honee->ceed; 148c2d90829SJames Wright NodalProjectionData projection = diff_flux_proj->projection; 149c2d90829SJames Wright CeedInt num_comp_q; 150c2d90829SJames Wright PetscInt dim, label_value = 0; 151c2d90829SJames Wright DMLabel domain_label = NULL; 152c2d90829SJames Wright CeedQFunctionContext advection_qfctx = NULL; 153c2d90829SJames Wright 154c2d90829SJames Wright PetscFunctionBeginUser; 155c2d90829SJames Wright // -- Get Pre-requisite things 156c2d90829SJames Wright PetscCall(DMGetDimension(projection->dm, &dim)); 157e3663b90SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q)); 158c2d90829SJames Wright 15940b78511SJames Wright { // Get advection-diffusion QF context 160c2d90829SJames Wright CeedOperator *sub_ops; 161c2d90829SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 162c2d90829SJames Wright 1630c373b74SJames Wright if (honee->op_ifunction) PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_ifunction, &sub_ops)); 1640c373b74SJames Wright else PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_rhs_ctx->op, &sub_ops)); 165c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &advection_qfctx)); 166c2d90829SJames Wright } 167c2d90829SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, op_rhs)); 168c2d90829SJames Wright { // Add the volume integral CeedOperator 169c2d90829SJames Wright CeedQFunction qf_rhs_volume = NULL; 170c2d90829SJames Wright CeedOperator op_rhs_volume; 171c2d90829SJames Wright CeedVector q_data; 172c2d90829SJames Wright CeedElemRestriction elem_restr_qd, elem_restr_diff_flux_volume = NULL; 173c2d90829SJames Wright CeedBasis basis_diff_flux = NULL; 174c2d90829SJames Wright CeedInt q_data_size; 175c2d90829SJames Wright 176c2d90829SJames Wright PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_diff_flux_volume, &basis_diff_flux, NULL, NULL)); 177e3663b90SJames 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, 178e3663b90SJames Wright &q_data_size)); 179c2d90829SJames Wright switch (dim) { 180c2d90829SJames Wright case 2: 181c2d90829SJames Wright PetscCallCeed( 182c2d90829SJames Wright ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_AdvDif_2D, DivDiffusiveFluxVolumeRHS_AdvDif_2D_loc, &qf_rhs_volume)); 183c2d90829SJames Wright break; 184c2d90829SJames Wright case 3: 185c2d90829SJames Wright PetscCallCeed( 186c2d90829SJames Wright ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_AdvDif_3D, DivDiffusiveFluxVolumeRHS_AdvDif_3D_loc, &qf_rhs_volume)); 187c2d90829SJames Wright break; 188c2d90829SJames Wright } 1890c373b74SJames Wright PetscCheck(qf_rhs_volume, honee->comm, PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim); 190c2d90829SJames Wright 191c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_volume, advection_qfctx)); 192c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 193c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "qdata", q_data_size, CEED_EVAL_NONE)); 194c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_volume, "diffusive flux RHS", projection->num_comp * dim, CEED_EVAL_GRAD)); 195c2d90829SJames Wright 196c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_volume, NULL, NULL, &op_rhs_volume)); 197e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 198c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 199c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "diffusive flux RHS", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE)); 200c2d90829SJames Wright 201c2d90829SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_volume)); 202c2d90829SJames Wright 203c2d90829SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 204c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 205c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_volume)); 206c2d90829SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 207c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_volume)); 208c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_volume)); 209c2d90829SJames Wright } 210c2d90829SJames Wright 211c2d90829SJames Wright { // Add the boundary integral CeedOperator 212c2d90829SJames Wright CeedQFunction qf_rhs_boundary; 213c2d90829SJames Wright DMLabel face_sets_label; 214c2d90829SJames Wright PetscInt num_face_set_values, *face_set_values; 215c2d90829SJames Wright CeedInt q_data_size; 216c2d90829SJames Wright 217c2d90829SJames Wright // -- Build RHS operator 218c2d90829SJames Wright switch (dim) { 219c2d90829SJames Wright case 2: 220c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_AdvDif_2D, DivDiffusiveFluxBoundaryRHS_AdvDif_2D_loc, 221c2d90829SJames Wright &qf_rhs_boundary)); 222c2d90829SJames Wright break; 223c2d90829SJames Wright case 3: 224c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_AdvDif_3D, DivDiffusiveFluxBoundaryRHS_AdvDif_3D_loc, 225c2d90829SJames Wright &qf_rhs_boundary)); 226c2d90829SJames Wright break; 227c2d90829SJames Wright } 228c2d90829SJames Wright 2290c373b74SJames Wright PetscCall(QDataBoundaryGradientGetNumComponents(honee->dm, &q_data_size)); 230c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_boundary, advection_qfctx)); 231c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 232c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "qdata", q_data_size, CEED_EVAL_NONE)); 233c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_boundary, "diffusive flux RHS", projection->num_comp, CEED_EVAL_INTERP)); 234c2d90829SJames Wright 235c2d90829SJames Wright PetscCall(DMGetLabel(projection->dm, "Face Sets", &face_sets_label)); 236c2d90829SJames Wright PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_sets_label, &num_face_set_values, &face_set_values)); 237c2d90829SJames Wright for (PetscInt f = 0; f < num_face_set_values; f++) { 238c2d90829SJames Wright DMLabel face_orientation_label; 239c2d90829SJames Wright PetscInt num_orientations_values, *orientation_values; 240c2d90829SJames Wright 241c2d90829SJames Wright { 242c2d90829SJames Wright char *face_orientation_label_name; 243c2d90829SJames Wright 244c2d90829SJames Wright PetscCall(DMPlexCreateFaceLabel(projection->dm, face_set_values[f], &face_orientation_label_name)); 245c2d90829SJames Wright PetscCall(DMGetLabel(projection->dm, face_orientation_label_name, &face_orientation_label)); 246c2d90829SJames Wright PetscCall(PetscFree(face_orientation_label_name)); 247c2d90829SJames Wright } 248c2d90829SJames Wright PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_orientation_label, &num_orientations_values, &orientation_values)); 249c2d90829SJames Wright for (PetscInt o = 0; o < num_orientations_values; o++) { 250c2d90829SJames Wright CeedOperator op_rhs_boundary; 251c2d90829SJames Wright CeedBasis basis_q, basis_diff_flux_boundary; 252c2d90829SJames Wright CeedElemRestriction elem_restr_qdata, elem_restr_q, elem_restr_diff_flux_boundary; 253c2d90829SJames Wright CeedVector q_data; 254c2d90829SJames Wright CeedInt q_data_size; 255c2d90829SJames Wright PetscInt orientation = orientation_values[o], dm_field_q = 0, height_cell = 0, height_face = 1; 256c2d90829SJames Wright 2570c373b74SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, face_orientation_label, orientation, height_cell, dm_field_q, &elem_restr_q)); 2580c373b74SJames Wright PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, honee->dm, face_orientation_label, orientation, orientation, dm_field_q, &basis_q)); 259c2d90829SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, 260c2d90829SJames Wright &elem_restr_diff_flux_boundary)); 261c2d90829SJames Wright PetscCall(CreateBasisFromPlex(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, &basis_diff_flux_boundary)); 262e3663b90SJames Wright PetscCall( 263e3663b90SJames Wright QDataBoundaryGradientGet(ceed, honee->dm, face_orientation_label, orientation, honee->x_coord, &elem_restr_qdata, &q_data, &q_data_size)); 264c2d90829SJames Wright 265c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_boundary, NULL, NULL, &op_rhs_boundary)); 266c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 267c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "qdata", elem_restr_qdata, CEED_BASIS_NONE, q_data)); 268c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "diffusive flux RHS", elem_restr_diff_flux_boundary, basis_diff_flux_boundary, 269c2d90829SJames Wright CEED_VECTOR_ACTIVE)); 270c2d90829SJames Wright 271c2d90829SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_boundary)); 272c2d90829SJames Wright 273c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_boundary)); 274c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qdata)); 275c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 276c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_boundary)); 277c2d90829SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 278c2d90829SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux_boundary)); 279c2d90829SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 280c2d90829SJames Wright } 281c2d90829SJames Wright PetscCall(PetscFree(orientation_values)); 282c2d90829SJames Wright } 283c2d90829SJames Wright PetscCall(PetscFree(face_set_values)); 284c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_boundary)); 285c2d90829SJames Wright } 286c2d90829SJames Wright 287c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&advection_qfctx)); 288c2d90829SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 289c2d90829SJames Wright } 290c2d90829SJames Wright 29140b78511SJames Wright /** 29240b78511SJames Wright @brief Create RHS CeedOperator for indirect projection of divergence of diffusive flux 29340b78511SJames Wright 2940c373b74SJames Wright @param[in] honee `Honee` context 29540b78511SJames Wright @param[in] diff_flux_proj `DivDiffFluxProjectionData` object 29640b78511SJames Wright @param[out] op_rhs Operator to calculate the RHS of the L^2 projection 29740b78511SJames Wright **/ 298e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Indirect_AdvDif(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) { 2990c373b74SJames Wright Ceed ceed = honee->ceed; 30040b78511SJames Wright NodalProjectionData projection = diff_flux_proj->projection; 30140b78511SJames Wright CeedBasis basis_diff_flux; 30240b78511SJames Wright CeedElemRestriction elem_restr_diff_flux, elem_restr_qd; 30340b78511SJames Wright CeedVector q_data; 30440b78511SJames Wright CeedInt num_comp_q, q_data_size; 30540b78511SJames Wright PetscInt dim; 30640b78511SJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 30740b78511SJames Wright DMLabel domain_label = NULL; 30840b78511SJames Wright CeedQFunction qf_rhs = NULL; 30940b78511SJames Wright CeedQFunctionContext advection_qfctx = NULL; 31040b78511SJames Wright 31140b78511SJames Wright PetscFunctionBeginUser; 31240b78511SJames Wright PetscCall(DMGetDimension(projection->dm, &dim)); 313e3663b90SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q)); 31440b78511SJames Wright 31540b78511SJames Wright { // Get advection-diffusion QF context 31640b78511SJames Wright CeedOperator *sub_ops; 31740b78511SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 31840b78511SJames Wright 3190c373b74SJames Wright if (honee->op_ifunction) PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_ifunction, &sub_ops)); 3200c373b74SJames Wright else PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_rhs_ctx->op, &sub_ops)); 32140b78511SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &advection_qfctx)); 32240b78511SJames Wright } 32340b78511SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, domain_label, label_value, height, dm_field, &elem_restr_diff_flux)); 32440b78511SJames Wright PetscCall(CreateBasisFromPlex(ceed, projection->dm, domain_label, label_value, height, dm_field, &basis_diff_flux)); 325e3663b90SJames 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, 326e3663b90SJames Wright &q_data_size)); 32740b78511SJames Wright 32840b78511SJames Wright switch (dim) { 32940b78511SJames Wright case 2: 33040b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_AdvDif_2D, DiffusiveFluxRHS_AdvDif_2D_loc, &qf_rhs)); 33140b78511SJames Wright break; 33240b78511SJames Wright case 3: 33340b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_AdvDif_3D, DiffusiveFluxRHS_AdvDif_3D_loc, &qf_rhs)); 33440b78511SJames Wright break; 33540b78511SJames Wright } 3360c373b74SJames Wright PetscCheck(qf_rhs, honee->comm, PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim); 33740b78511SJames Wright 33840b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs, advection_qfctx)); 33940b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 34040b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "qdata", q_data_size, CEED_EVAL_NONE)); 34140b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs, "F_diff RHS", projection->num_comp, CEED_EVAL_INTERP)); 34240b78511SJames Wright 34340b78511SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, op_rhs)); 344e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 34540b78511SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 34640b78511SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "F_diff RHS", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE)); 34740b78511SJames Wright 34840b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs)); 34940b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&advection_qfctx)); 35040b78511SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 35140b78511SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 35240b78511SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 35340b78511SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux)); 35440b78511SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 35540b78511SJames Wright } 35640b78511SJames Wright 357d6cac220SJames Wright static PetscErrorCode AdvectionInflowBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) { 358d6cac220SJames Wright HoneeBCStruct honee_bc; 359d6cac220SJames Wright DM dm; 360d6cac220SJames Wright PetscInt dim; 361d6cac220SJames Wright 362d6cac220SJames Wright PetscFunctionBeginUser; 363d6cac220SJames Wright PetscCall(BCDefinitionGetContext(bc_def, &honee_bc)); 364d6cac220SJames Wright PetscCall(BCDefinitionGetDM(bc_def, &dm)); 365d6cac220SJames Wright PetscCall(DMGetDimension(dm, &dim)); 366d6cac220SJames Wright switch (dim) { 367d6cac220SJames Wright case 2: 368d6cac220SJames Wright PetscCall(HoneeBCCreateIFunctionQF(bc_def, Advection2d_InOutFlow, Advection2d_InOutFlow_loc, honee_bc->qfctx, qf)); 369d6cac220SJames Wright break; 370d6cac220SJames Wright case 3: 371d6cac220SJames Wright PetscCall(HoneeBCCreateIFunctionQF(bc_def, Advection_InOutFlow, Advection_InOutFlow_loc, honee_bc->qfctx, qf)); 372d6cac220SJames Wright break; 373d6cac220SJames Wright } 374d6cac220SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 375d6cac220SJames Wright } 376d6cac220SJames Wright 377f27c2204SJames Wright static const char *const component_names[] = {"Density", "MomentumX", "MomentumY", "MomentumZ", "TotalEnergy"}; 378f27c2204SJames Wright 379d3c60affSJames Wright PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx) { 3805f636aeaSJames Wright AdvDifWindType wind_type; 3815f636aeaSJames Wright AdvDifICType advectionic_type; 382108c6689SJames Wright AdvDifBubbleContinuityType bubble_continuity_type = -1; 383a515125bSLeila Ghaffari StabilizationType stab; 38457272ee0SJames Wright StabilizationTauType stab_tau; 3853636f6a4SJames Wright SetupContextAdv setup_context; 3860c373b74SJames Wright Honee honee = *(Honee *)ctx; 3870c373b74SJames Wright MPI_Comm comm = honee->comm; 3880c373b74SJames Wright Ceed ceed = honee->ceed; 389a515125bSLeila Ghaffari PetscBool implicit; 39015a3537eSJed Brown AdvectionContext advection_ctx; 391e07531f7SJames Wright CeedQFunctionContext advection_qfctx; 3929529d636SJames Wright PetscInt dim; 393a515125bSLeila Ghaffari 39415a3537eSJed Brown PetscFunctionBeginUser; 395*2d898fa6SJames Wright PetscCall(PetscNew(&setup_context)); 396*2d898fa6SJames Wright PetscCall(PetscNew(&advection_ctx)); 3979529d636SJames Wright PetscCall(DMGetDimension(dm, &dim)); 398a515125bSLeila Ghaffari 399a515125bSLeila Ghaffari // ------------------------------------------------------ 400a515125bSLeila Ghaffari // SET UP ADVECTION 401a515125bSLeila Ghaffari // ------------------------------------------------------ 40228160fc2SJames Wright problem->print_info = PRINT_ADVECTION; 4031abc2837SJames Wright problem->num_comps_jac_data = 0; 4049529d636SJames Wright switch (dim) { 4059529d636SJames Wright case 2: 406e07531f7SJames Wright problem->ics.qf_func_ptr = ICsAdvection2d; 407e07531f7SJames Wright problem->ics.qf_loc = ICsAdvection2d_loc; 408e07531f7SJames Wright problem->apply_vol_rhs.qf_func_ptr = RHS_Advection2d; 409e07531f7SJames Wright problem->apply_vol_rhs.qf_loc = RHS_Advection2d_loc; 410e07531f7SJames Wright problem->apply_vol_ifunction.qf_func_ptr = IFunction_Advection2d; 411e07531f7SJames Wright problem->apply_vol_ifunction.qf_loc = IFunction_Advection2d_loc; 41258ce1233SJames Wright problem->compute_exact_solution_error = PETSC_TRUE; 4139529d636SJames Wright break; 4149529d636SJames Wright case 3: 415e07531f7SJames Wright problem->ics.qf_func_ptr = ICsAdvection; 416e07531f7SJames Wright problem->ics.qf_loc = ICsAdvection_loc; 417e07531f7SJames Wright problem->apply_vol_rhs.qf_func_ptr = RHS_Advection; 418e07531f7SJames Wright problem->apply_vol_rhs.qf_loc = RHS_Advection_loc; 419e07531f7SJames Wright problem->apply_vol_ifunction.qf_func_ptr = IFunction_Advection; 420e07531f7SJames Wright problem->apply_vol_ifunction.qf_loc = IFunction_Advection_loc; 42158ce1233SJames Wright problem->compute_exact_solution_error = PETSC_FALSE; 4229529d636SJames Wright break; 4239529d636SJames Wright } 424a515125bSLeila Ghaffari 425f27c2204SJames Wright problem->num_components = 5; 426f27c2204SJames Wright PetscCall(PetscMalloc1(problem->num_components, &problem->component_names)); 427f27c2204SJames Wright for (PetscInt i = 0; i < 5; i++) PetscCall(PetscStrallocpy(component_names[i], &problem->component_names[i])); 428f27c2204SJames Wright 429d9e69621SJames Wright PetscCall(DivDiffFluxProjectionCreate(honee, honee->app_ctx->divFdiffproj_method, 1, &honee->diff_flux_proj)); 4300c373b74SJames Wright if (honee->diff_flux_proj) { 4310c373b74SJames Wright DivDiffFluxProjectionData diff_flux_proj = honee->diff_flux_proj; 432c2d90829SJames Wright NodalProjectionData projection = diff_flux_proj->projection; 433c2d90829SJames Wright 434c2d90829SJames Wright diff_flux_proj->CreateRHSOperator_Direct = DivDiffFluxProjectionCreateRHS_Direct_AdvDif; 43540b78511SJames Wright diff_flux_proj->CreateRHSOperator_Indirect = DivDiffFluxProjectionCreateRHS_Indirect_AdvDif; 436c2d90829SJames Wright 4370c373b74SJames Wright switch (honee->diff_flux_proj->method) { 438c2d90829SJames Wright case DIV_DIFF_FLUX_PROJ_DIRECT: { 439c2d90829SJames Wright PetscSection section; 440c2d90829SJames Wright 441c2d90829SJames Wright PetscCall(DMGetLocalSection(projection->dm, §ion)); 442c2d90829SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 443c2d90829SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "DivDiffusiveFlux_Scalar")); 444c2d90829SJames Wright } break; 445c2d90829SJames Wright case DIV_DIFF_FLUX_PROJ_INDIRECT: { 446c2d90829SJames Wright PetscSection section; 447c2d90829SJames Wright 448c2d90829SJames Wright PetscCall(DMGetLocalSection(projection->dm, §ion)); 449c2d90829SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 450c2d90829SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "DiffusiveFlux_ScalarX")); 451c2d90829SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "DiffusiveFlux_ScalarY")); 45240b78511SJames Wright if (dim >= 3) PetscCall(PetscSectionSetComponentName(section, 0, 2, "DiffusiveFlux_ScalarZ")); 453c2d90829SJames Wright } break; 454c2d90829SJames Wright case DIV_DIFF_FLUX_PROJ_NONE: 4550c373b74SJames Wright SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s", 4560c373b74SJames Wright DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]); 457c2d90829SJames Wright break; 458c2d90829SJames Wright } 459c2d90829SJames Wright } 460c2d90829SJames Wright 461a515125bSLeila Ghaffari // ------------------------------------------------------ 462ea615d4cSJames Wright // Create the QFunction context 463a515125bSLeila Ghaffari // ------------------------------------------------------ 464a515125bSLeila Ghaffari CeedScalar rc = 1000.; // m (Radius of bubble) 465a515125bSLeila Ghaffari CeedScalar CtauS = 0.; // dimensionless 466059d1c40SJames Wright PetscBool strong_form = PETSC_FALSE; 467a515125bSLeila Ghaffari CeedScalar E_wind = 1.e6; // J 4680c373b74SJames Wright CeedScalar Ctau_a = PetscPowScalarInt(honee->app_ctx->degree, 2); 4690c373b74SJames Wright CeedScalar Ctau_d = PetscPowScalarInt(honee->app_ctx->degree, 4); 47057272ee0SJames Wright CeedScalar Ctau_t = 0.; 471a515125bSLeila Ghaffari PetscReal wind[3] = {1., 0, 0}; // m/s 472c8d249deSJames Wright CeedScalar diffusion_coeff = 0.; 473a62be6baSJames Wright CeedScalar wave_frequency = 2 * M_PI; 474a62be6baSJames Wright CeedScalar wave_phase = 0; 475108c6689SJames Wright AdvDifWaveType wave_type = -1; 476b4fd18dfSJames Wright PetscScalar bl_height_factor = 1.; 47787d3884fSJeremy L Thompson PetscReal domain_min[3], domain_max[3], domain_size[3] = {0.}; 4782b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 47966d54740SJames Wright for (PetscInt i = 0; i < dim; i++) domain_size[i] = domain_max[i] - domain_min[i]; 48005a512bdSLeila Ghaffari 481a515125bSLeila Ghaffari // ------------------------------------------------------ 482a515125bSLeila Ghaffari // Command line Options 483a515125bSLeila Ghaffari // ------------------------------------------------------ 4841485969bSJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for ADVECTION problem", NULL); 485a515125bSLeila Ghaffari // -- Physics 486a515125bSLeila Ghaffari PetscBool translation; 4875f636aeaSJames Wright PetscCall(PetscOptionsEnum("-wind_type", "Wind type in Advection", NULL, AdvDifWindTypes, (PetscEnum)(wind_type = ADVDIF_WIND_ROTATION), 4885f636aeaSJames Wright (PetscEnum *)&wind_type, &translation)); 48966d54740SJames Wright PetscInt n = dim; 490a515125bSLeila Ghaffari PetscBool user_wind; 4912b916ea7SJeremy L Thompson PetscCall(PetscOptionsRealArray("-wind_translation", "Constant wind vector", NULL, wind, &n, &user_wind)); 492c8d249deSJames Wright PetscCall(PetscOptionsScalar("-diffusion_coeff", "Diffusion coefficient", NULL, diffusion_coeff, &diffusion_coeff, NULL)); 4932b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-CtauS", "Scale coefficient for tau (nondimensional)", NULL, CtauS, &CtauS, NULL)); 494059d1c40SJames Wright PetscCall(PetscOptionsBool("-strong_form", "Strong (true) or weak/integrated by parts (false) advection residual", NULL, strong_form, &strong_form, 495059d1c40SJames Wright NULL)); 4962b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-E_wind", "Total energy of inflow wind", NULL, E_wind, &E_wind, NULL)); 4972b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 49857272ee0SJames Wright PetscCall(PetscOptionsEnum("-stab_tau", "Stabilization constant, tau", NULL, StabilizationTauTypes, (PetscEnum)(stab_tau = STAB_TAU_CTAU), 49957272ee0SJames Wright (PetscEnum *)&stab_tau, NULL)); 50057272ee0SJames Wright PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 501fbabb365SJames Wright PetscCall(PetscOptionsScalar("-Ctau_a", "Coefficient for the stabilization, advection component", NULL, Ctau_a, &Ctau_a, NULL)); 502fbabb365SJames Wright PetscCall(PetscOptionsScalar("-Ctau_d", "Coefficient for the stabilization, diffusion component", NULL, Ctau_d, &Ctau_d, NULL)); 5032b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 50480e9ac5bSJames Wright PetscCall(PetscOptionsEnum("-advection_ic_type", "Initial condition for Advection problem", NULL, AdvDifICTypes, 50580e9ac5bSJames Wright (PetscEnum)(advectionic_type = ADVDIF_IC_BUBBLE_SPHERE), (PetscEnum *)&advectionic_type, NULL)); 50657bb6b22SJames Wright // IC-specific options 50757bb6b22SJames Wright switch (advectionic_type) { 50857bb6b22SJames Wright case ADVDIF_IC_WAVE: 50980e9ac5bSJames Wright PetscCall(PetscOptionsDeprecated("-wave_type", "-advection_ic_wave_type", "HONEE 0.0", NULL)); 51080e9ac5bSJames Wright PetscCall(PetscOptionsDeprecated("-wave_frequency", "-advection_ic_wave_frequency", "HONEE 0.0", NULL)); 51180e9ac5bSJames Wright PetscCall(PetscOptionsDeprecated("-wave_phase", "-advection_ic_wave_phase", "HONEE 0.0", NULL)); 51280e9ac5bSJames Wright PetscCall(PetscOptionsEnum("-advection_ic_wave_type", "Type of wave", NULL, AdvDifWaveTypes, (PetscEnum)(wave_type = ADVDIF_WAVE_SINE), 513108c6689SJames Wright (PetscEnum *)&wave_type, NULL)); 51480e9ac5bSJames Wright PetscCall(PetscOptionsScalar("-advection_ic_wave_frequency", "Frequency of sine wave", NULL, wave_frequency, &wave_frequency, NULL)); 51580e9ac5bSJames Wright PetscCall(PetscOptionsScalar("-advection_ic_wave_phase", "Length correction", NULL, wave_phase, &wave_phase, NULL)); 51657bb6b22SJames Wright break; 51757bb6b22SJames Wright case ADVDIF_IC_BOUNDARY_LAYER: 51880e9ac5bSJames Wright PetscCall( 51980e9ac5bSJames Wright PetscOptionsScalar("-advection_ic_bl_height_factor", "Height of boundary layer in IC", NULL, bl_height_factor, &bl_height_factor, NULL)); 52057bb6b22SJames Wright break; 52157bb6b22SJames Wright case ADVDIF_IC_BUBBLE_CYLINDER: 52257bb6b22SJames Wright case ADVDIF_IC_BUBBLE_SPHERE: 52380e9ac5bSJames Wright PetscCall(PetscOptionsDeprecated("-rc", "-advection_ic_bubble_rc", "HONEE 0.0", NULL)); 52480e9ac5bSJames Wright PetscCall(PetscOptionsDeprecated("-bubble_continuity", "-advection_ic_bubble_continuity", "HONEE 0.0", NULL)); 52580e9ac5bSJames Wright PetscCall(PetscOptionsScalar("-advection_ic_bubble_rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL)); 526108c6689SJames Wright bubble_continuity_type = dim == 3 ? ADVDIF_BUBBLE_CONTINUITY_SMOOTH : ADVDIF_BUBBLE_CONTINUITY_COSINE; 52780e9ac5bSJames Wright PetscCall(PetscOptionsEnum("-advection_ic_bubble_continuity", "Smooth, back_sharp, or thick", NULL, AdvDifBubbleContinuityTypes, 528108c6689SJames Wright (PetscEnum)bubble_continuity_type, (PetscEnum *)&bubble_continuity_type, NULL)); 52957bb6b22SJames Wright break; 53057bb6b22SJames Wright case ADVDIF_IC_SKEW: 53157bb6b22SJames Wright case ADVDIF_IC_COSINE_HILL: 53257bb6b22SJames Wright break; 533108c6689SJames Wright } 534a62be6baSJames Wright 535a515125bSLeila Ghaffari // -- Warnings 5365f636aeaSJames Wright if (wind_type == ADVDIF_WIND_ROTATION && user_wind) { 5372b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Use -wind_translation only with -wind_type translation\n")); 538a515125bSLeila Ghaffari } 5395f636aeaSJames Wright if (wind_type == ADVDIF_WIND_TRANSLATION && advectionic_type == ADVDIF_IC_BUBBLE_CYLINDER && wind[2] != 0.) { 540a515125bSLeila Ghaffari wind[2] = 0; 541c51f031aSJames Wright PetscCall( 542c51f031aSJames Wright PetscPrintf(comm, "Warning! Background wind in the z direction should be zero (-wind_translation x,x,0) with -advection_ic_type cylinder\n")); 543a515125bSLeila Ghaffari } 544a515125bSLeila Ghaffari if (stab == STAB_NONE && CtauS != 0) { 5452b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Use -CtauS only with -stab su or -stab supg\n")); 546a515125bSLeila Ghaffari } 5471485969bSJeremy L Thompson PetscOptionsEnd(); 548a515125bSLeila Ghaffari 549a78efa86SJames Wright if (stab == STAB_SUPG) problem->create_mass_operator = CreateKSPMassOperator_AdvectionStabilized; 550a78efa86SJames Wright 551a515125bSLeila Ghaffari // ------------------------------------------------------ 552ea615d4cSJames Wright // Set up the QFunction contexts 553a515125bSLeila Ghaffari // ------------------------------------------------------ 554a515125bSLeila Ghaffari // -- Scale variables to desired units 555c9f37605SMohammed Amin Units units = honee->units; 556c9f37605SMohammed Amin E_wind *= units->Joule; 557c9f37605SMohammed Amin rc = fabs(rc) * units->meter; 55866d54740SJames Wright for (PetscInt i = 0; i < dim; i++) { 559c9f37605SMohammed Amin wind[i] *= (units->meter / units->second); 560c9f37605SMohammed Amin domain_size[i] *= units->meter; 56105a512bdSLeila Ghaffari } 562a515125bSLeila Ghaffari 563a515125bSLeila Ghaffari // -- Setup Context 564a515125bSLeila Ghaffari setup_context->rc = rc; 56505a512bdSLeila Ghaffari setup_context->lx = domain_size[0]; 56605a512bdSLeila Ghaffari setup_context->ly = domain_size[1]; 56766d54740SJames Wright setup_context->lz = dim == 3 ? domain_size[2] : 0.; 568a515125bSLeila Ghaffari setup_context->wind[0] = wind[0]; 569a515125bSLeila Ghaffari setup_context->wind[1] = wind[1]; 57066d54740SJames Wright setup_context->wind[2] = dim == 3 ? wind[2] : 0.; 571a515125bSLeila Ghaffari setup_context->wind_type = wind_type; 572c51f031aSJames Wright setup_context->initial_condition_type = advectionic_type; 573a515125bSLeila Ghaffari setup_context->bubble_continuity_type = bubble_continuity_type; 574a515125bSLeila Ghaffari setup_context->time = 0; 575a62be6baSJames Wright setup_context->wave_frequency = wave_frequency; 576a62be6baSJames Wright setup_context->wave_phase = wave_phase; 577a62be6baSJames Wright setup_context->wave_type = wave_type; 578b4fd18dfSJames Wright setup_context->bl_height_factor = bl_height_factor; 579a515125bSLeila Ghaffari 580a515125bSLeila Ghaffari // -- QFunction Context 5810c373b74SJames Wright honee->phys->implicit = implicit; 58215a3537eSJed Brown advection_ctx->CtauS = CtauS; 58315a3537eSJed Brown advection_ctx->E_wind = E_wind; 58415a3537eSJed Brown advection_ctx->implicit = implicit; 58515a3537eSJed Brown advection_ctx->strong_form = strong_form; 58615a3537eSJed Brown advection_ctx->stabilization = stab; 58757272ee0SJames Wright advection_ctx->stabilization_tau = stab_tau; 58857272ee0SJames Wright advection_ctx->Ctau_a = Ctau_a; 589fbabb365SJames Wright advection_ctx->Ctau_d = Ctau_d; 59057272ee0SJames Wright advection_ctx->Ctau_t = Ctau_t; 591c8d249deSJames Wright advection_ctx->diffusion_coeff = diffusion_coeff; 5920c373b74SJames Wright advection_ctx->divFdiff_method = honee->app_ctx->divFdiffproj_method; 593a515125bSLeila Ghaffari 5940c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &problem->ics.qfctx)); 595e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(problem->ics.qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 596e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfctx, CEED_MEM_HOST, FreeContextPetsc)); 597a515125bSLeila Ghaffari 5980c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &advection_qfctx)); 599e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(advection_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*advection_ctx), advection_ctx)); 600e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(advection_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 601e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(advection_qfctx, "timestep size", offsetof(struct AdvectionContext_, dt), 1, 60257272ee0SJames Wright "Size of timestep, delta t")); 603e07531f7SJames Wright problem->apply_vol_rhs.qfctx = advection_qfctx; 604e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &problem->apply_vol_ifunction.qfctx)); 605d6cac220SJames Wright 606d6cac220SJames Wright for (PetscCount b = 0; b < problem->num_bc_defs; b++) { 607d6cac220SJames Wright BCDefinition bc_def = problem->bc_defs[b]; 608d6cac220SJames Wright const char *name; 609d6cac220SJames Wright 610d6cac220SJames Wright PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL)); 611d6cac220SJames Wright if (!strcmp(name, "inflow")) { 612d6cac220SJames Wright HoneeBCStruct honee_bc; 613d6cac220SJames Wright 614d6cac220SJames Wright PetscCall(PetscNew(&honee_bc)); 615d6cac220SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &honee_bc->qfctx)); 616d6cac220SJames Wright honee_bc->honee = honee; 6171abc2837SJames Wright honee_bc->num_comps_jac_data = 0; 618d6cac220SJames Wright PetscCall(BCDefinitionSetContext(bc_def, HoneeBCDestroy, honee_bc)); 619d6cac220SJames Wright 620d6cac220SJames Wright PetscCall(BCDefinitionSetIFunction(bc_def, AdvectionInflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp)); 621d6cac220SJames Wright PetscCall(BCDefinitionSetIJacobian(bc_def, NULL, NULL)); 622d6cac220SJames Wright } 623d6cac220SJames Wright } 624d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 625a515125bSLeila Ghaffari } 626