1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3 4 /// @file 5 /// Utility functions for setting up EULER_VORTEX 6 7 #include "../qfunctions/eulervortex.h" 8 9 #include <ceed.h> 10 #include <petscdm.h> 11 12 #include <navierstokes.h> 13 14 static PetscErrorCode PRINT_EULER_VORTEX(Honee honee, ProblemData problem, AppCtx app_ctx) { 15 MPI_Comm comm = honee->comm; 16 Ceed ceed = honee->ceed; 17 EulerContext euler_ctx; 18 19 PetscFunctionBeginUser; 20 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfctx, CEED_MEM_HOST, &euler_ctx)); 21 PetscCall(PetscPrintf(comm, 22 " Problem:\n" 23 " Problem Name : %s\n" 24 " Test Case : %s\n" 25 " Background Velocity : %f,%f,%f\n" 26 " Stabilization : %s\n", 27 app_ctx->problem_name, EulerTestTypes[euler_ctx->euler_test], euler_ctx->mean_velocity[0], euler_ctx->mean_velocity[1], 28 euler_ctx->mean_velocity[2], StabilizationTypes[euler_ctx->stabilization])); 29 30 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfctx, &euler_ctx)); 31 PetscFunctionReturn(PETSC_SUCCESS); 32 } 33 34 static PetscErrorCode EulerVortexOutflowBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) { 35 HoneeBCStruct honee_bc; 36 37 PetscFunctionBeginUser; 38 PetscCall(BCDefinitionGetContext(bc_def, &honee_bc)); 39 PetscCall(HoneeBCCreateIFunctionQF(bc_def, Euler_Outflow, Euler_Outflow_loc, honee_bc->qfctx, qf)); 40 PetscFunctionReturn(PETSC_SUCCESS); 41 } 42 43 static PetscErrorCode EulerVortexInflowBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) { 44 HoneeBCStruct honee_bc; 45 46 PetscFunctionBeginUser; 47 PetscCall(BCDefinitionGetContext(bc_def, &honee_bc)); 48 PetscCall(HoneeBCCreateIFunctionQF(bc_def, TravelingVortex_Inflow, TravelingVortex_Inflow_loc, honee_bc->qfctx, qf)); 49 PetscFunctionReturn(PETSC_SUCCESS); 50 } 51 52 static const char *const component_names[] = {"Density", "MomentumX", "MomentumY", "MomentumZ", "TotalEnergy"}; 53 54 PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx) { 55 EulerTestType euler_test; 56 Honee honee = *(Honee *)ctx; 57 StabilizationType stab; 58 MPI_Comm comm = honee->comm; 59 Ceed ceed = honee->ceed; 60 PetscBool implicit; 61 EulerContext euler_ctx; 62 CeedQFunctionContext euler_qfctx; 63 PetscInt dim; 64 65 PetscFunctionBeginUser; 66 PetscCall(PetscCalloc1(1, &euler_ctx)); 67 68 // ------------------------------------------------------ 69 // SET UP EULER VORTEX 70 // ------------------------------------------------------ 71 problem->ics.qf_func_ptr = ICsEuler; 72 problem->ics.qf_loc = ICsEuler_loc; 73 problem->apply_vol_rhs.qf_func_ptr = Euler; 74 problem->apply_vol_rhs.qf_loc = Euler_loc; 75 problem->apply_vol_ifunction.qf_func_ptr = IFunction_Euler; 76 problem->apply_vol_ifunction.qf_loc = IFunction_Euler_loc; 77 problem->num_comps_jac_data = 0; 78 problem->compute_exact_solution_error = PETSC_TRUE; 79 problem->print_info = PRINT_EULER_VORTEX; 80 81 problem->num_components = 5; 82 PetscCall(PetscMalloc1(problem->num_components, &problem->component_names)); 83 for (PetscInt i = 0; i < 5; i++) PetscCall(PetscStrallocpy(component_names[i], &problem->component_names[i])); 84 85 // ------------------------------------------------------ 86 // Create the libCEED context 87 // ------------------------------------------------------ 88 CeedScalar vortex_strength = 5.; // - 89 CeedScalar c_tau = 0.5; // - 90 // c_tau = 0.5 is reported as "optimal" in Hughes et al 2010 91 PetscReal center[3] = {0.}, // m 92 mean_velocity[3] = {1., 1., 0}; // m/s 93 PetscReal domain_min[3], domain_max[3], domain_size[3] = {0.}; 94 PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 95 PetscCall(DMGetDimension(dm, &dim)); 96 for (PetscInt i = 0; i < dim; i++) domain_size[i] = domain_max[i] - domain_min[i]; 97 98 // ------------------------------------------------------ 99 // Command line Options 100 // ------------------------------------------------------ 101 PetscOptionsBegin(comm, NULL, "Options for EULER_VORTEX problem", NULL); 102 // -- Physics 103 PetscCall(PetscOptionsScalar("-vortex_strength", "Strength of Vortex", NULL, vortex_strength, &vortex_strength, NULL)); 104 PetscInt n = dim; 105 PetscBool user_velocity; 106 PetscCall(PetscOptionsRealArray("-mean_velocity", "Background velocity vector", NULL, mean_velocity, &n, &user_velocity)); 107 for (PetscInt i = 0; i < dim; i++) center[i] = .5 * domain_size[i]; 108 n = dim; 109 PetscCall(PetscOptionsRealArray("-center", "Location of vortex center", NULL, center, &n, NULL)); 110 PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 111 PetscCall(PetscOptionsEnum("-euler_test", "Euler test option", NULL, EulerTestTypes, (PetscEnum)(euler_test = EULER_TEST_ISENTROPIC_VORTEX), 112 (PetscEnum *)&euler_test, NULL)); 113 PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 114 PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 115 116 // -- Warnings 117 if (stab == STAB_SUPG && !implicit) { 118 PetscCall(PetscPrintf(comm, "Warning! Use -stab supg only with -implicit\n")); 119 } 120 if (user_velocity && (euler_test == EULER_TEST_1 || euler_test == EULER_TEST_3)) { 121 PetscCall(PetscPrintf(comm, "Warning! Background velocity vector for -euler_test t1 and -euler_test t3 is (0,0,0)\n")); 122 } 123 124 PetscOptionsEnd(); 125 126 // ------------------------------------------------------ 127 // Set up the libCEED context 128 // ------------------------------------------------------ 129 // -- Scale variables to desired units 130 Units units = honee->units; 131 for (PetscInt i = 0; i < 3; i++) { 132 center[i] *= units->meter; 133 domain_size[i] *= units->meter; 134 mean_velocity[i] *= (units->meter / units->second); 135 } 136 137 // -- QFunction Context 138 honee->phys->implicit = implicit; 139 euler_ctx->curr_time = 0.; 140 euler_ctx->implicit = implicit; 141 euler_ctx->euler_test = euler_test; 142 euler_ctx->center[0] = center[0]; 143 euler_ctx->center[1] = center[1]; 144 euler_ctx->center[2] = center[2]; 145 euler_ctx->vortex_strength = vortex_strength; 146 euler_ctx->c_tau = c_tau; 147 euler_ctx->mean_velocity[0] = mean_velocity[0]; 148 euler_ctx->mean_velocity[1] = mean_velocity[1]; 149 euler_ctx->mean_velocity[2] = mean_velocity[2]; 150 euler_ctx->stabilization = stab; 151 152 PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &euler_qfctx)); 153 PetscCallCeed(ceed, CeedQFunctionContextSetData(euler_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*euler_ctx), euler_ctx)); 154 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(euler_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 155 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(euler_qfctx, "solution time", offsetof(struct EulerContext_, curr_time), 1, 156 "Physical time of the solution")); 157 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->ics.qfctx)); 158 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->apply_vol_rhs.qfctx)); 159 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->apply_vol_ifunction.qfctx)); 160 161 for (PetscCount b = 0; b < problem->num_bc_defs; b++) { 162 BCDefinition bc_def = problem->bc_defs[b]; 163 const char *name; 164 165 PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL)); 166 if (!strcmp(name, "outflow")) { 167 HoneeBCStruct honee_bc; 168 169 PetscCall(PetscNew(&honee_bc)); 170 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &honee_bc->qfctx)); 171 honee_bc->honee = honee; 172 honee_bc->num_comps_jac_data = 0; 173 PetscCall(BCDefinitionSetContext(bc_def, HoneeBCDestroy, honee_bc)); 174 175 PetscCall(BCDefinitionSetIFunction(bc_def, EulerVortexOutflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp)); 176 PetscCall(BCDefinitionSetIJacobian(bc_def, NULL, NULL)); 177 } else if (!strcmp(name, "inflow")) { 178 HoneeBCStruct honee_bc; 179 180 PetscCall(PetscNew(&honee_bc)); 181 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &honee_bc->qfctx)); 182 honee_bc->honee = honee; 183 honee_bc->num_comps_jac_data = 0; 184 PetscCall(BCDefinitionSetContext(bc_def, HoneeBCDestroy, honee_bc)); 185 186 PetscCall(BCDefinitionSetIFunction(bc_def, EulerVortexInflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp)); 187 PetscCall(BCDefinitionSetIJacobian(bc_def, NULL, NULL)); 188 } 189 } 190 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&euler_qfctx)); 191 PetscFunctionReturn(PETSC_SUCCESS); 192 } 193