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 93da7f3a07SJames Wright PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(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; 150*e3db12f8SJames Wright PetscInt dim; 151c2d90829SJames Wright CeedQFunctionContext advection_qfctx = NULL; 152c2d90829SJames Wright 153c2d90829SJames Wright PetscFunctionBeginUser; 154c2d90829SJames Wright // -- Get Pre-requisite things 155c2d90829SJames Wright PetscCall(DMGetDimension(projection->dm, &dim)); 156e3663b90SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q)); 157c2d90829SJames Wright 15840b78511SJames Wright { // Get advection-diffusion QF context 159c2d90829SJames Wright CeedOperator *sub_ops; 160c2d90829SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 161c2d90829SJames Wright 162da7f3a07SJames Wright if (honee->op_ifunction) PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_ifunction, &sub_ops)); 163da7f3a07SJames Wright else PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_rhs_ctx->op, &sub_ops)); 164c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &advection_qfctx)); 165c2d90829SJames Wright } 166da7f3a07SJames Wright PetscCallCeed(ceed, CeedOperatorCreateComposite(ceed, op_rhs)); 167c2d90829SJames Wright { // Add the volume integral CeedOperator 168c2d90829SJames Wright CeedQFunction qf_rhs_volume = NULL; 169c2d90829SJames Wright CeedOperator op_rhs_volume; 170c2d90829SJames Wright CeedVector q_data; 171c2d90829SJames Wright CeedElemRestriction elem_restr_qd, elem_restr_diff_flux_volume = NULL; 172c2d90829SJames Wright CeedBasis basis_diff_flux = NULL; 173c2d90829SJames Wright CeedInt q_data_size; 174c2d90829SJames Wright 175c2d90829SJames Wright PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_diff_flux_volume, &basis_diff_flux, NULL, NULL)); 176*e3db12f8SJames Wright PetscCall(QDataGet(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, honee->elem_restr_x, honee->basis_x, honee->x_coord, 177*e3db12f8SJames Wright &elem_restr_qd, &q_data, &q_data_size)); 178c2d90829SJames Wright switch (dim) { 179c2d90829SJames Wright case 2: 1809eadbee4SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_AdvDif_2D, DivDiffusiveFluxVolumeRHS_AdvDif_2D_loc, 1819eadbee4SJames Wright &qf_rhs_volume)); 182c2d90829SJames Wright break; 183c2d90829SJames Wright case 3: 1849eadbee4SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_AdvDif_3D, DivDiffusiveFluxVolumeRHS_AdvDif_3D_loc, 1859eadbee4SJames Wright &qf_rhs_volume)); 186c2d90829SJames Wright break; 187c2d90829SJames Wright } 1880c373b74SJames Wright PetscCheck(qf_rhs_volume, honee->comm, PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim); 189c2d90829SJames Wright 190c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_volume, advection_qfctx)); 191c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 192c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "qdata", q_data_size, CEED_EVAL_NONE)); 193c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_volume, "diffusive flux RHS", projection->num_comp * dim, CEED_EVAL_GRAD)); 194c2d90829SJames Wright 195c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_volume, NULL, NULL, &op_rhs_volume)); 196e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 197c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 198c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "diffusive flux RHS", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE)); 199c2d90829SJames Wright 200da7f3a07SJames Wright PetscCallCeed(ceed, CeedOperatorCompositeAddSub(*op_rhs, op_rhs_volume)); 201c2d90829SJames Wright 202c2d90829SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 203c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 204c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_volume)); 205c2d90829SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 206c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_volume)); 207c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_volume)); 208c2d90829SJames Wright } 209c2d90829SJames Wright 210c2d90829SJames Wright { // Add the boundary integral CeedOperator 211c2d90829SJames Wright CeedQFunction qf_rhs_boundary; 212c2d90829SJames Wright DMLabel face_sets_label; 213c2d90829SJames Wright PetscInt num_face_set_values, *face_set_values; 214c2d90829SJames Wright CeedInt q_data_size; 215c2d90829SJames Wright 216c2d90829SJames Wright // -- Build RHS operator 217c2d90829SJames Wright switch (dim) { 218c2d90829SJames Wright case 2: 219c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_AdvDif_2D, DivDiffusiveFluxBoundaryRHS_AdvDif_2D_loc, 220c2d90829SJames Wright &qf_rhs_boundary)); 221c2d90829SJames Wright break; 222c2d90829SJames Wright case 3: 223c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_AdvDif_3D, DivDiffusiveFluxBoundaryRHS_AdvDif_3D_loc, 224c2d90829SJames Wright &qf_rhs_boundary)); 225c2d90829SJames Wright break; 226c2d90829SJames Wright } 227c2d90829SJames Wright 2280c373b74SJames Wright PetscCall(QDataBoundaryGradientGetNumComponents(honee->dm, &q_data_size)); 229c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_boundary, advection_qfctx)); 230c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 231c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "qdata", q_data_size, CEED_EVAL_NONE)); 232c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_boundary, "diffusive flux RHS", projection->num_comp, CEED_EVAL_INTERP)); 233c2d90829SJames Wright 234c2d90829SJames Wright PetscCall(DMGetLabel(projection->dm, "Face Sets", &face_sets_label)); 235c2d90829SJames Wright PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_sets_label, &num_face_set_values, &face_set_values)); 236c2d90829SJames Wright for (PetscInt f = 0; f < num_face_set_values; f++) { 237c2d90829SJames Wright DMLabel face_orientation_label; 238c2d90829SJames Wright PetscInt num_orientations_values, *orientation_values; 239c2d90829SJames Wright 240c2d90829SJames Wright { 241c2d90829SJames Wright char *face_orientation_label_name; 242c2d90829SJames Wright 243c2d90829SJames Wright PetscCall(DMPlexCreateFaceLabel(projection->dm, face_set_values[f], &face_orientation_label_name)); 244c2d90829SJames Wright PetscCall(DMGetLabel(projection->dm, face_orientation_label_name, &face_orientation_label)); 245c2d90829SJames Wright PetscCall(PetscFree(face_orientation_label_name)); 246c2d90829SJames Wright } 247c2d90829SJames Wright PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_orientation_label, &num_orientations_values, &orientation_values)); 248c2d90829SJames Wright for (PetscInt o = 0; o < num_orientations_values; o++) { 249c2d90829SJames Wright CeedOperator op_rhs_boundary; 250c2d90829SJames Wright CeedBasis basis_q, basis_diff_flux_boundary; 251c2d90829SJames Wright CeedElemRestriction elem_restr_qdata, elem_restr_q, elem_restr_diff_flux_boundary; 252c2d90829SJames Wright CeedVector q_data; 253c2d90829SJames Wright CeedInt q_data_size; 254c2d90829SJames Wright PetscInt orientation = orientation_values[o], dm_field_q = 0, height_cell = 0, height_face = 1; 255c2d90829SJames Wright 2560c373b74SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, face_orientation_label, orientation, height_cell, dm_field_q, &elem_restr_q)); 2570c373b74SJames Wright PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, honee->dm, face_orientation_label, orientation, orientation, dm_field_q, &basis_q)); 258c2d90829SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, 259c2d90829SJames Wright &elem_restr_diff_flux_boundary)); 2606b9fd993SJames Wright PetscCall(DMPlexCeedBasisCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, &basis_diff_flux_boundary)); 2619eadbee4SJames Wright PetscCall(QDataBoundaryGradientGet(ceed, honee->dm, face_orientation_label, orientation, honee->x_coord, &elem_restr_qdata, &q_data, 2629eadbee4SJames Wright &q_data_size)); 263c2d90829SJames Wright 264c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_boundary, NULL, NULL, &op_rhs_boundary)); 265c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 266c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "qdata", elem_restr_qdata, CEED_BASIS_NONE, q_data)); 267c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "diffusive flux RHS", elem_restr_diff_flux_boundary, basis_diff_flux_boundary, 268c2d90829SJames Wright CEED_VECTOR_ACTIVE)); 269c2d90829SJames Wright 270da7f3a07SJames Wright PetscCallCeed(ceed, CeedOperatorCompositeAddSub(*op_rhs, op_rhs_boundary)); 271c2d90829SJames Wright 272c2d90829SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_boundary)); 273c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qdata)); 274c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 275c2d90829SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_boundary)); 276c2d90829SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 277c2d90829SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux_boundary)); 278c2d90829SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 279c2d90829SJames Wright } 280c2d90829SJames Wright PetscCall(PetscFree(orientation_values)); 281c2d90829SJames Wright } 282c2d90829SJames Wright PetscCall(PetscFree(face_set_values)); 283c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_boundary)); 284c2d90829SJames Wright } 285c2d90829SJames Wright 286c2d90829SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&advection_qfctx)); 287c2d90829SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 288c2d90829SJames Wright } 289c2d90829SJames Wright 29040b78511SJames Wright /** 29140b78511SJames Wright @brief Create RHS CeedOperator for indirect projection of divergence of diffusive flux 29240b78511SJames Wright 2930c373b74SJames Wright @param[in] honee `Honee` context 29440b78511SJames Wright @param[in] diff_flux_proj `DivDiffFluxProjectionData` object 29540b78511SJames Wright @param[out] op_rhs Operator to calculate the RHS of the L^2 projection 29640b78511SJames Wright **/ 297e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Indirect_AdvDif(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) { 2980c373b74SJames Wright Ceed ceed = honee->ceed; 29940b78511SJames Wright NodalProjectionData projection = diff_flux_proj->projection; 30040b78511SJames Wright CeedBasis basis_diff_flux; 30140b78511SJames Wright CeedElemRestriction elem_restr_diff_flux, elem_restr_qd; 30240b78511SJames Wright CeedVector q_data; 30340b78511SJames Wright CeedInt num_comp_q, q_data_size; 30440b78511SJames Wright PetscInt dim; 305*e3db12f8SJames Wright PetscInt height = 0, dm_field = 0; 30640b78511SJames Wright CeedQFunction qf_rhs = NULL; 30740b78511SJames Wright CeedQFunctionContext advection_qfctx = NULL; 30840b78511SJames Wright 30940b78511SJames Wright PetscFunctionBeginUser; 31040b78511SJames Wright PetscCall(DMGetDimension(projection->dm, &dim)); 311e3663b90SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q)); 31240b78511SJames Wright 31340b78511SJames Wright { // Get advection-diffusion QF context 31440b78511SJames Wright CeedOperator *sub_ops; 31540b78511SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 31640b78511SJames Wright 317da7f3a07SJames Wright if (honee->op_ifunction) PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_ifunction, &sub_ops)); 318da7f3a07SJames Wright else PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_rhs_ctx->op, &sub_ops)); 31940b78511SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &advection_qfctx)); 32040b78511SJames Wright } 321*e3db12f8SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &elem_restr_diff_flux)); 322*e3db12f8SJames Wright PetscCall(DMPlexCeedBasisCreate(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &basis_diff_flux)); 323*e3db12f8SJames Wright PetscCall(QDataGet(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, honee->elem_restr_x, honee->basis_x, honee->x_coord, 324*e3db12f8SJames Wright &elem_restr_qd, &q_data, &q_data_size)); 32540b78511SJames Wright 32640b78511SJames Wright switch (dim) { 32740b78511SJames Wright case 2: 32840b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_AdvDif_2D, DiffusiveFluxRHS_AdvDif_2D_loc, &qf_rhs)); 32940b78511SJames Wright break; 33040b78511SJames Wright case 3: 33140b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_AdvDif_3D, DiffusiveFluxRHS_AdvDif_3D_loc, &qf_rhs)); 33240b78511SJames Wright break; 33340b78511SJames Wright } 3340c373b74SJames Wright PetscCheck(qf_rhs, honee->comm, PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim); 33540b78511SJames Wright 33640b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs, advection_qfctx)); 33740b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD)); 33840b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "qdata", q_data_size, CEED_EVAL_NONE)); 33940b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs, "F_diff RHS", projection->num_comp, CEED_EVAL_INTERP)); 34040b78511SJames Wright 34140b78511SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, op_rhs)); 342e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 34340b78511SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 34440b78511SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "F_diff RHS", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE)); 34540b78511SJames Wright 34640b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs)); 34740b78511SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&advection_qfctx)); 34840b78511SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux)); 34940b78511SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 35040b78511SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 35140b78511SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux)); 35240b78511SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 35340b78511SJames Wright } 35440b78511SJames Wright 355d6cac220SJames Wright static PetscErrorCode AdvectionInflowBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) { 356d6cac220SJames Wright HoneeBCStruct honee_bc; 357d6cac220SJames Wright DM dm; 358d6cac220SJames Wright PetscInt dim; 359d6cac220SJames Wright 360d6cac220SJames Wright PetscFunctionBeginUser; 361d6cac220SJames Wright PetscCall(BCDefinitionGetContext(bc_def, &honee_bc)); 362d6cac220SJames Wright PetscCall(BCDefinitionGetDM(bc_def, &dm)); 363d6cac220SJames Wright PetscCall(DMGetDimension(dm, &dim)); 364d6cac220SJames Wright switch (dim) { 365d6cac220SJames Wright case 2: 366d6cac220SJames Wright PetscCall(HoneeBCCreateIFunctionQF(bc_def, Advection2d_InOutFlow, Advection2d_InOutFlow_loc, honee_bc->qfctx, qf)); 367d6cac220SJames Wright break; 368d6cac220SJames Wright case 3: 369d6cac220SJames Wright PetscCall(HoneeBCCreateIFunctionQF(bc_def, Advection_InOutFlow, Advection_InOutFlow_loc, honee_bc->qfctx, qf)); 370d6cac220SJames Wright break; 371d6cac220SJames Wright } 372d6cac220SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 373d6cac220SJames Wright } 374d6cac220SJames Wright 375f27c2204SJames Wright static const char *const component_names[] = {"Density", "MomentumX", "MomentumY", "MomentumZ", "TotalEnergy"}; 376f27c2204SJames Wright 377d3c60affSJames Wright PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx) { 3785f636aeaSJames Wright AdvDifWindType wind_type; 3795f636aeaSJames Wright AdvDifICType advectionic_type; 380108c6689SJames Wright AdvDifBubbleContinuityType bubble_continuity_type = -1; 381a515125bSLeila Ghaffari StabilizationType stab; 38257272ee0SJames Wright StabilizationTauType stab_tau; 3833636f6a4SJames Wright SetupContextAdv setup_context; 3840c373b74SJames Wright Honee honee = *(Honee *)ctx; 3850c373b74SJames Wright MPI_Comm comm = honee->comm; 3860c373b74SJames Wright Ceed ceed = honee->ceed; 387a515125bSLeila Ghaffari PetscBool implicit; 38815a3537eSJed Brown AdvectionContext advection_ctx; 389f5dc303cSJames Wright CeedQFunctionContext advection_qfctx, ics_qfctx; 3909529d636SJames Wright PetscInt dim; 391a515125bSLeila Ghaffari 39215a3537eSJed Brown PetscFunctionBeginUser; 3932d898fa6SJames Wright PetscCall(PetscNew(&setup_context)); 3942d898fa6SJames Wright PetscCall(PetscNew(&advection_ctx)); 3959529d636SJames Wright PetscCall(DMGetDimension(dm, &dim)); 396a515125bSLeila Ghaffari 397a515125bSLeila Ghaffari // ------------------------------------------------------ 398ea615d4cSJames Wright // Create the QFunction context 399a515125bSLeila Ghaffari // ------------------------------------------------------ 400a515125bSLeila Ghaffari CeedScalar rc = 1000.; // m (Radius of bubble) 401a515125bSLeila Ghaffari CeedScalar CtauS = 0.; // dimensionless 402059d1c40SJames Wright PetscBool strong_form = PETSC_FALSE; 403a515125bSLeila Ghaffari CeedScalar E_wind = 1.e6; // J 4040c373b74SJames Wright CeedScalar Ctau_a = PetscPowScalarInt(honee->app_ctx->degree, 2); 4050c373b74SJames Wright CeedScalar Ctau_d = PetscPowScalarInt(honee->app_ctx->degree, 4); 40657272ee0SJames Wright CeedScalar Ctau_t = 0.; 407a515125bSLeila Ghaffari PetscReal wind[3] = {1., 0, 0}; // m/s 408c8d249deSJames Wright CeedScalar diffusion_coeff = 0.; 409a62be6baSJames Wright CeedScalar wave_frequency = 2 * M_PI; 410a62be6baSJames Wright CeedScalar wave_phase = 0; 411108c6689SJames Wright AdvDifWaveType wave_type = -1; 412b4fd18dfSJames Wright PetscScalar bl_height_factor = 1.; 41387d3884fSJeremy L Thompson PetscReal domain_min[3], domain_max[3], domain_size[3] = {0.}; 4142b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 41566d54740SJames Wright for (PetscInt i = 0; i < dim; i++) domain_size[i] = domain_max[i] - domain_min[i]; 41605a512bdSLeila Ghaffari 417a515125bSLeila Ghaffari // ------------------------------------------------------ 418a515125bSLeila Ghaffari // Command line Options 419a515125bSLeila Ghaffari // ------------------------------------------------------ 4201485969bSJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for ADVECTION problem", NULL); 421a515125bSLeila Ghaffari // -- Physics 422a515125bSLeila Ghaffari PetscBool translation; 4235f636aeaSJames Wright PetscCall(PetscOptionsEnum("-wind_type", "Wind type in Advection", NULL, AdvDifWindTypes, (PetscEnum)(wind_type = ADVDIF_WIND_ROTATION), 4245f636aeaSJames Wright (PetscEnum *)&wind_type, &translation)); 42566d54740SJames Wright PetscInt n = dim; 426a515125bSLeila Ghaffari PetscBool user_wind; 4272b916ea7SJeremy L Thompson PetscCall(PetscOptionsRealArray("-wind_translation", "Constant wind vector", NULL, wind, &n, &user_wind)); 428c8d249deSJames Wright PetscCall(PetscOptionsScalar("-diffusion_coeff", "Diffusion coefficient", NULL, diffusion_coeff, &diffusion_coeff, NULL)); 4292b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-CtauS", "Scale coefficient for tau (nondimensional)", NULL, CtauS, &CtauS, NULL)); 430059d1c40SJames Wright PetscCall(PetscOptionsBool("-strong_form", "Strong (true) or weak/integrated by parts (false) advection residual", NULL, strong_form, &strong_form, 431059d1c40SJames Wright NULL)); 4322b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-E_wind", "Total energy of inflow wind", NULL, E_wind, &E_wind, NULL)); 4332b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 43457272ee0SJames Wright PetscCall(PetscOptionsEnum("-stab_tau", "Stabilization constant, tau", NULL, StabilizationTauTypes, (PetscEnum)(stab_tau = STAB_TAU_CTAU), 43557272ee0SJames Wright (PetscEnum *)&stab_tau, NULL)); 43657272ee0SJames Wright PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 437fbabb365SJames Wright PetscCall(PetscOptionsScalar("-Ctau_a", "Coefficient for the stabilization, advection component", NULL, Ctau_a, &Ctau_a, NULL)); 438fbabb365SJames Wright PetscCall(PetscOptionsScalar("-Ctau_d", "Coefficient for the stabilization, diffusion component", NULL, Ctau_d, &Ctau_d, NULL)); 4392b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 44080e9ac5bSJames Wright PetscCall(PetscOptionsEnum("-advection_ic_type", "Initial condition for Advection problem", NULL, AdvDifICTypes, 44180e9ac5bSJames Wright (PetscEnum)(advectionic_type = ADVDIF_IC_BUBBLE_SPHERE), (PetscEnum *)&advectionic_type, NULL)); 44257bb6b22SJames Wright // IC-specific options 44357bb6b22SJames Wright switch (advectionic_type) { 44457bb6b22SJames Wright case ADVDIF_IC_WAVE: 44580e9ac5bSJames Wright PetscCall(PetscOptionsDeprecated("-wave_type", "-advection_ic_wave_type", "HONEE 0.0", NULL)); 44680e9ac5bSJames Wright PetscCall(PetscOptionsDeprecated("-wave_frequency", "-advection_ic_wave_frequency", "HONEE 0.0", NULL)); 44780e9ac5bSJames Wright PetscCall(PetscOptionsDeprecated("-wave_phase", "-advection_ic_wave_phase", "HONEE 0.0", NULL)); 44880e9ac5bSJames Wright PetscCall(PetscOptionsEnum("-advection_ic_wave_type", "Type of wave", NULL, AdvDifWaveTypes, (PetscEnum)(wave_type = ADVDIF_WAVE_SINE), 449108c6689SJames Wright (PetscEnum *)&wave_type, NULL)); 45080e9ac5bSJames Wright PetscCall(PetscOptionsScalar("-advection_ic_wave_frequency", "Frequency of sine wave", NULL, wave_frequency, &wave_frequency, NULL)); 45180e9ac5bSJames Wright PetscCall(PetscOptionsScalar("-advection_ic_wave_phase", "Length correction", NULL, wave_phase, &wave_phase, NULL)); 45257bb6b22SJames Wright break; 45357bb6b22SJames Wright case ADVDIF_IC_BOUNDARY_LAYER: 4549eadbee4SJames Wright PetscCall(PetscOptionsScalar("-advection_ic_bl_height_factor", "Height of boundary layer in IC", NULL, bl_height_factor, &bl_height_factor, 4559eadbee4SJames Wright NULL)); 45657bb6b22SJames Wright break; 45757bb6b22SJames Wright case ADVDIF_IC_BUBBLE_CYLINDER: 45857bb6b22SJames Wright case ADVDIF_IC_BUBBLE_SPHERE: 45980e9ac5bSJames Wright PetscCall(PetscOptionsDeprecated("-rc", "-advection_ic_bubble_rc", "HONEE 0.0", NULL)); 46080e9ac5bSJames Wright PetscCall(PetscOptionsDeprecated("-bubble_continuity", "-advection_ic_bubble_continuity", "HONEE 0.0", NULL)); 46180e9ac5bSJames Wright PetscCall(PetscOptionsScalar("-advection_ic_bubble_rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL)); 462108c6689SJames Wright bubble_continuity_type = dim == 3 ? ADVDIF_BUBBLE_CONTINUITY_SMOOTH : ADVDIF_BUBBLE_CONTINUITY_COSINE; 46380e9ac5bSJames Wright PetscCall(PetscOptionsEnum("-advection_ic_bubble_continuity", "Smooth, back_sharp, or thick", NULL, AdvDifBubbleContinuityTypes, 464108c6689SJames Wright (PetscEnum)bubble_continuity_type, (PetscEnum *)&bubble_continuity_type, NULL)); 46557bb6b22SJames Wright break; 46657bb6b22SJames Wright case ADVDIF_IC_SKEW: 46757bb6b22SJames Wright case ADVDIF_IC_COSINE_HILL: 46857bb6b22SJames Wright break; 469108c6689SJames Wright } 470a62be6baSJames Wright 471a515125bSLeila Ghaffari // -- Warnings 4725f636aeaSJames Wright if (wind_type == ADVDIF_WIND_ROTATION && user_wind) { 4732b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Use -wind_translation only with -wind_type translation\n")); 474a515125bSLeila Ghaffari } 4755f636aeaSJames Wright if (wind_type == ADVDIF_WIND_TRANSLATION && advectionic_type == ADVDIF_IC_BUBBLE_CYLINDER && wind[2] != 0.) { 476a515125bSLeila Ghaffari wind[2] = 0; 4779eadbee4SJames Wright PetscCall(PetscPrintf(comm, 4789eadbee4SJames Wright "Warning! Background wind in the z direction should be zero (-wind_translation x,x,0) with -advection_ic_type cylinder\n")); 479a515125bSLeila Ghaffari } 480a515125bSLeila Ghaffari if (stab == STAB_NONE && CtauS != 0) { 4812b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Use -CtauS only with -stab su or -stab supg\n")); 482a515125bSLeila Ghaffari } 4831485969bSJeremy L Thompson PetscOptionsEnd(); 484a515125bSLeila Ghaffari 485a78efa86SJames Wright if (stab == STAB_SUPG) problem->create_mass_operator = CreateKSPMassOperator_AdvectionStabilized; 486a78efa86SJames Wright 487a515125bSLeila Ghaffari // ------------------------------------------------------ 488ea615d4cSJames Wright // Set up the QFunction contexts 489a515125bSLeila Ghaffari // ------------------------------------------------------ 490a515125bSLeila Ghaffari // -- Scale variables to desired units 491c9f37605SMohammed Amin Units units = honee->units; 492c9f37605SMohammed Amin E_wind *= units->Joule; 493c9f37605SMohammed Amin rc = fabs(rc) * units->meter; 49466d54740SJames Wright for (PetscInt i = 0; i < dim; i++) { 495c9f37605SMohammed Amin wind[i] *= (units->meter / units->second); 496c9f37605SMohammed Amin domain_size[i] *= units->meter; 49705a512bdSLeila Ghaffari } 498a515125bSLeila Ghaffari 499a515125bSLeila Ghaffari // -- Setup Context 500a515125bSLeila Ghaffari setup_context->rc = rc; 50105a512bdSLeila Ghaffari setup_context->lx = domain_size[0]; 50205a512bdSLeila Ghaffari setup_context->ly = domain_size[1]; 50366d54740SJames Wright setup_context->lz = dim == 3 ? domain_size[2] : 0.; 504a515125bSLeila Ghaffari setup_context->wind[0] = wind[0]; 505a515125bSLeila Ghaffari setup_context->wind[1] = wind[1]; 50666d54740SJames Wright setup_context->wind[2] = dim == 3 ? wind[2] : 0.; 507a515125bSLeila Ghaffari setup_context->wind_type = wind_type; 508c51f031aSJames Wright setup_context->initial_condition_type = advectionic_type; 509a515125bSLeila Ghaffari setup_context->bubble_continuity_type = bubble_continuity_type; 510a515125bSLeila Ghaffari setup_context->time = 0; 511a62be6baSJames Wright setup_context->wave_frequency = wave_frequency; 512a62be6baSJames Wright setup_context->wave_phase = wave_phase; 513a62be6baSJames Wright setup_context->wave_type = wave_type; 514b4fd18dfSJames Wright setup_context->bl_height_factor = bl_height_factor; 515a515125bSLeila Ghaffari 516a515125bSLeila Ghaffari // -- QFunction Context 5170c373b74SJames Wright honee->phys->implicit = implicit; 51815a3537eSJed Brown advection_ctx->CtauS = CtauS; 51915a3537eSJed Brown advection_ctx->E_wind = E_wind; 52015a3537eSJed Brown advection_ctx->implicit = implicit; 52115a3537eSJed Brown advection_ctx->strong_form = strong_form; 52215a3537eSJed Brown advection_ctx->stabilization = stab; 52357272ee0SJames Wright advection_ctx->stabilization_tau = stab_tau; 52457272ee0SJames Wright advection_ctx->Ctau_a = Ctau_a; 525fbabb365SJames Wright advection_ctx->Ctau_d = Ctau_d; 52657272ee0SJames Wright advection_ctx->Ctau_t = Ctau_t; 527c8d249deSJames Wright advection_ctx->diffusion_coeff = diffusion_coeff; 5280c373b74SJames Wright advection_ctx->divFdiff_method = honee->app_ctx->divFdiffproj_method; 529a515125bSLeila Ghaffari 530f5dc303cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &ics_qfctx)); 531f5dc303cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(ics_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 532f5dc303cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(ics_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 533a515125bSLeila Ghaffari 5340c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &advection_qfctx)); 535e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(advection_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*advection_ctx), advection_ctx)); 536e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(advection_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 537e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(advection_qfctx, "timestep size", offsetof(struct AdvectionContext_, dt), 1, 53857272ee0SJames Wright "Size of timestep, delta t")); 539f5dc303cSJames Wright 540f5dc303cSJames Wright // ------------------------------------------------------ 541f5dc303cSJames Wright // SET UP ADVECTION 542f5dc303cSJames Wright // ------------------------------------------------------ 543f5dc303cSJames Wright problem->print_info = PRINT_ADVECTION; 544f5dc303cSJames Wright problem->num_comps_jac_data = 0; 545f5dc303cSJames Wright switch (dim) { 546f5dc303cSJames Wright case 2: 547f5dc303cSJames Wright problem->ics = (HoneeQFSpec){.qf_func_ptr = ICsAdvection2d, .qf_loc = ICsAdvection2d_loc}; 548f5dc303cSJames Wright problem->apply_vol_rhs = (HoneeQFSpec){.qf_func_ptr = RHS_Advection2d, .qf_loc = RHS_Advection2d_loc}; 549f5dc303cSJames Wright problem->apply_vol_ifunction = (HoneeQFSpec){.qf_func_ptr = IFunction_Advection2d, .qf_loc = IFunction_Advection2d_loc}; 550f5dc303cSJames Wright problem->compute_exact_solution_error = PETSC_TRUE; 551f5dc303cSJames Wright break; 552f5dc303cSJames Wright case 3: 553f5dc303cSJames Wright problem->ics = (HoneeQFSpec){.qf_func_ptr = ICsAdvection, .qf_loc = ICsAdvection_loc}; 554f5dc303cSJames Wright problem->apply_vol_rhs = (HoneeQFSpec){.qf_func_ptr = RHS_Advection, .qf_loc = RHS_Advection_loc}; 555f5dc303cSJames Wright problem->apply_vol_ifunction = (HoneeQFSpec){.qf_func_ptr = IFunction_Advection, .qf_loc = IFunction_Advection_loc}; 556f5dc303cSJames Wright problem->compute_exact_solution_error = PETSC_FALSE; 557f5dc303cSJames Wright break; 558f5dc303cSJames Wright } 559f5dc303cSJames Wright problem->ics.qfctx = ics_qfctx; 560e07531f7SJames Wright problem->apply_vol_rhs.qfctx = advection_qfctx; 561e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &problem->apply_vol_ifunction.qfctx)); 562d6cac220SJames Wright 563f5dc303cSJames Wright problem->num_components = 5; 564f5dc303cSJames Wright PetscCall(PetscMalloc1(problem->num_components, &problem->component_names)); 565f5dc303cSJames Wright for (PetscInt i = 0; i < 5; i++) PetscCall(PetscStrallocpy(component_names[i], &problem->component_names[i])); 566f5dc303cSJames Wright 567f5dc303cSJames Wright PetscCall(DivDiffFluxProjectionCreate(honee, honee->app_ctx->divFdiffproj_method, 1, &honee->diff_flux_proj)); 568f5dc303cSJames Wright if (honee->diff_flux_proj) { 569f5dc303cSJames Wright DivDiffFluxProjectionData diff_flux_proj = honee->diff_flux_proj; 570f5dc303cSJames Wright NodalProjectionData projection = diff_flux_proj->projection; 571ae2e5212SJames Wright PetscSection section; 572f5dc303cSJames Wright 573f5dc303cSJames Wright diff_flux_proj->CreateRHSOperator_Direct = DivDiffFluxProjectionCreateRHS_Direct_AdvDif; 574f5dc303cSJames Wright diff_flux_proj->CreateRHSOperator_Indirect = DivDiffFluxProjectionCreateRHS_Indirect_AdvDif; 575ae2e5212SJames Wright PetscCall(DMGetLocalSection(projection->dm, §ion)); 576f5dc303cSJames Wright switch (honee->diff_flux_proj->method) { 577f5dc303cSJames Wright case DIV_DIFF_FLUX_PROJ_DIRECT: { 578f5dc303cSJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 579f5dc303cSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "DivDiffusiveFlux_Scalar")); 580f5dc303cSJames Wright } break; 581f5dc303cSJames Wright case DIV_DIFF_FLUX_PROJ_INDIRECT: { 582f5dc303cSJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 583f5dc303cSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "DiffusiveFlux_ScalarX")); 584f5dc303cSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "DiffusiveFlux_ScalarY")); 585f5dc303cSJames Wright if (dim >= 3) PetscCall(PetscSectionSetComponentName(section, 0, 2, "DiffusiveFlux_ScalarZ")); 586f5dc303cSJames Wright } break; 587f5dc303cSJames Wright case DIV_DIFF_FLUX_PROJ_NONE: 588f5dc303cSJames Wright SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s", 589f5dc303cSJames Wright DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]); 590f5dc303cSJames Wright break; 591f5dc303cSJames Wright } 592f5dc303cSJames Wright } 593f5dc303cSJames Wright 594d6cac220SJames Wright for (PetscCount b = 0; b < problem->num_bc_defs; b++) { 595d6cac220SJames Wright BCDefinition bc_def = problem->bc_defs[b]; 596d6cac220SJames Wright const char *name; 597d6cac220SJames Wright 598d6cac220SJames Wright PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL)); 599d6cac220SJames Wright if (!strcmp(name, "inflow")) { 600d6cac220SJames Wright HoneeBCStruct honee_bc; 601d6cac220SJames Wright 602d6cac220SJames Wright PetscCall(PetscNew(&honee_bc)); 603d6cac220SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &honee_bc->qfctx)); 604d6cac220SJames Wright honee_bc->honee = honee; 6051abc2837SJames Wright honee_bc->num_comps_jac_data = 0; 606d6cac220SJames Wright PetscCall(BCDefinitionSetContext(bc_def, HoneeBCDestroy, honee_bc)); 607d6cac220SJames Wright 608d6cac220SJames Wright PetscCall(BCDefinitionSetIFunction(bc_def, AdvectionInflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp)); 609d6cac220SJames Wright PetscCall(BCDefinitionSetIJacobian(bc_def, NULL, NULL)); 610d6cac220SJames Wright } 611d6cac220SJames Wright } 612d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 613a515125bSLeila Ghaffari } 614