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 PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 15 EulerTestType euler_test; 16 User user = *(User *)ctx; 17 StabilizationType stab; 18 MPI_Comm comm = user->comm; 19 Ceed ceed = user->ceed; 20 PetscBool implicit; 21 EulerContext euler_ctx; 22 CeedQFunctionContext euler_qfctx; 23 PetscInt dim; 24 25 PetscFunctionBeginUser; 26 PetscCall(PetscCalloc1(1, &euler_ctx)); 27 28 // ------------------------------------------------------ 29 // SET UP DENSITY_CURRENT 30 // ------------------------------------------------------ 31 problem->ics.qf_func_ptr = ICsEuler; 32 problem->ics.qf_loc = ICsEuler_loc; 33 problem->apply_vol_rhs.qf_func_ptr = Euler; 34 problem->apply_vol_rhs.qf_loc = Euler_loc; 35 problem->apply_vol_ifunction.qf_func_ptr = IFunction_Euler; 36 problem->apply_vol_ifunction.qf_loc = IFunction_Euler_loc; 37 problem->apply_inflow.qf_func_ptr = TravelingVortex_Inflow; 38 problem->apply_inflow.qf_loc = TravelingVortex_Inflow_loc; 39 problem->apply_outflow.qf_func_ptr = Euler_Outflow; 40 problem->apply_outflow.qf_loc = Euler_Outflow_loc; 41 problem->compute_exact_solution_error = PETSC_TRUE; 42 problem->print_info = PRINT_EULER_VORTEX; 43 44 // ------------------------------------------------------ 45 // Create the libCEED context 46 // ------------------------------------------------------ 47 CeedScalar vortex_strength = 5.; // - 48 CeedScalar c_tau = 0.5; // - 49 // c_tau = 0.5 is reported as "optimal" in Hughes et al 2010 50 PetscReal center[3] = {0.}, // m 51 mean_velocity[3] = {1., 1., 0}; // m/s 52 PetscReal domain_min[3], domain_max[3], domain_size[3] = {0.}; 53 PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 54 PetscCall(DMGetDimension(dm, &dim)); 55 for (PetscInt i = 0; i < dim; i++) domain_size[i] = domain_max[i] - domain_min[i]; 56 57 // ------------------------------------------------------ 58 // Create the PETSc context 59 // ------------------------------------------------------ 60 PetscScalar meter = 1e-2; // 1 meter in scaled length units 61 PetscScalar second = 1e-2; // 1 second in scaled time units 62 63 // ------------------------------------------------------ 64 // Command line Options 65 // ------------------------------------------------------ 66 PetscOptionsBegin(comm, NULL, "Options for EULER_VORTEX problem", NULL); 67 // -- Physics 68 PetscCall(PetscOptionsScalar("-vortex_strength", "Strength of Vortex", NULL, vortex_strength, &vortex_strength, NULL)); 69 PetscInt n = dim; 70 PetscBool user_velocity; 71 PetscCall(PetscOptionsRealArray("-mean_velocity", "Background velocity vector", NULL, mean_velocity, &n, &user_velocity)); 72 for (PetscInt i = 0; i < dim; i++) center[i] = .5 * domain_size[i]; 73 n = dim; 74 PetscCall(PetscOptionsRealArray("-center", "Location of vortex center", NULL, center, &n, NULL)); 75 PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 76 PetscCall(PetscOptionsEnum("-euler_test", "Euler test option", NULL, EulerTestTypes, (PetscEnum)(euler_test = EULER_TEST_ISENTROPIC_VORTEX), 77 (PetscEnum *)&euler_test, NULL)); 78 PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 79 PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 80 // -- Units 81 PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 82 meter = fabs(meter); 83 PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 84 second = fabs(second); 85 86 // -- Warnings 87 if (stab == STAB_SUPG && !implicit) { 88 PetscCall(PetscPrintf(comm, "Warning! Use -stab supg only with -implicit\n")); 89 } 90 if (user_velocity && (euler_test == EULER_TEST_1 || euler_test == EULER_TEST_3)) { 91 PetscCall(PetscPrintf(comm, "Warning! Background velocity vector for -euler_test t1 and -euler_test t3 is (0,0,0)\n")); 92 } 93 94 PetscOptionsEnd(); 95 96 // ------------------------------------------------------ 97 // Set up the PETSc context 98 // ------------------------------------------------------ 99 user->units->meter = meter; 100 user->units->second = second; 101 102 // ------------------------------------------------------ 103 // Set up the libCEED context 104 // ------------------------------------------------------ 105 // -- Scale variables to desired units 106 for (PetscInt i = 0; i < 3; i++) { 107 center[i] *= meter; 108 domain_size[i] *= meter; 109 mean_velocity[i] *= (meter / second); 110 } 111 problem->dm_scale = meter; 112 113 // -- QFunction Context 114 user->phys->implicit = implicit; 115 euler_ctx->curr_time = 0.; 116 euler_ctx->implicit = implicit; 117 euler_ctx->euler_test = euler_test; 118 euler_ctx->center[0] = center[0]; 119 euler_ctx->center[1] = center[1]; 120 euler_ctx->center[2] = center[2]; 121 euler_ctx->vortex_strength = vortex_strength; 122 euler_ctx->c_tau = c_tau; 123 euler_ctx->mean_velocity[0] = mean_velocity[0]; 124 euler_ctx->mean_velocity[1] = mean_velocity[1]; 125 euler_ctx->mean_velocity[2] = mean_velocity[2]; 126 euler_ctx->stabilization = stab; 127 128 PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &euler_qfctx)); 129 PetscCallCeed(ceed, CeedQFunctionContextSetData(euler_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*euler_ctx), euler_ctx)); 130 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(euler_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 131 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(euler_qfctx, "solution time", offsetof(struct EulerContext_, curr_time), 1, 132 "Physical time of the solution")); 133 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->ics.qfctx)); 134 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->apply_vol_rhs.qfctx)); 135 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->apply_vol_ifunction.qfctx)); 136 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->apply_inflow.qfctx)); 137 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->apply_outflow.qfctx)); 138 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&euler_qfctx)); 139 PetscFunctionReturn(PETSC_SUCCESS); 140 } 141 142 PetscErrorCode PRINT_EULER_VORTEX(User user, ProblemData problem, AppCtx app_ctx) { 143 MPI_Comm comm = user->comm; 144 Ceed ceed = user->ceed; 145 EulerContext euler_ctx; 146 147 PetscFunctionBeginUser; 148 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfctx, CEED_MEM_HOST, &euler_ctx)); 149 PetscCall(PetscPrintf(comm, 150 " Problem:\n" 151 " Problem Name : %s\n" 152 " Test Case : %s\n" 153 " Background Velocity : %f,%f,%f\n" 154 " Stabilization : %s\n", 155 app_ctx->problem_name, EulerTestTypes[euler_ctx->euler_test], euler_ctx->mean_velocity[0], euler_ctx->mean_velocity[1], 156 euler_ctx->mean_velocity[2], StabilizationTypes[euler_ctx->stabilization])); 157 158 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfctx, &euler_ctx)); 159 PetscFunctionReturn(PETSC_SUCCESS); 160 } 161