1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 /// @file 9 /// Utility functions for setting up EULER_VORTEX 10 11 #include "../qfunctions/eulervortex.h" 12 13 #include <ceed.h> 14 #include <petscdm.h> 15 16 #include "../navierstokes.h" 17 #include "../qfunctions/setupgeo.h" 18 19 PetscErrorCode NS_EULER_VORTEX(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) { 20 EulerTestType euler_test; 21 User user = *(User *)ctx; 22 StabilizationType stab; 23 MPI_Comm comm = PETSC_COMM_WORLD; 24 PetscBool implicit; 25 PetscBool has_curr_time = PETSC_TRUE; 26 PetscBool has_neumann = PETSC_TRUE; 27 EulerContext euler_ctx; 28 CeedQFunctionContext euler_context; 29 30 PetscFunctionBeginUser; 31 PetscCall(PetscCalloc1(1, &euler_ctx)); 32 33 // ------------------------------------------------------ 34 // SET UP DENSITY_CURRENT 35 // ------------------------------------------------------ 36 problem->dim = 3; 37 problem->q_data_size_vol = 10; 38 problem->q_data_size_sur = 10; 39 problem->setup_vol.qfunction = Setup; 40 problem->setup_vol.qfunction_loc = Setup_loc; 41 problem->setup_sur.qfunction = SetupBoundary; 42 problem->setup_sur.qfunction_loc = SetupBoundary_loc; 43 problem->ics.qfunction = ICsEuler; 44 problem->ics.qfunction_loc = ICsEuler_loc; 45 problem->apply_vol_rhs.qfunction = Euler; 46 problem->apply_vol_rhs.qfunction_loc = Euler_loc; 47 problem->apply_vol_ifunction.qfunction = IFunction_Euler; 48 problem->apply_vol_ifunction.qfunction_loc = IFunction_Euler_loc; 49 problem->apply_inflow.qfunction = TravelingVortex_Inflow; 50 problem->apply_inflow.qfunction_loc = TravelingVortex_Inflow_loc; 51 problem->apply_outflow.qfunction = Euler_Outflow; 52 problem->apply_outflow.qfunction_loc = Euler_Outflow_loc; 53 problem->bc = Exact_Euler; 54 problem->bc_ctx = euler_ctx; 55 problem->non_zero_time = PETSC_TRUE; 56 problem->print_info = PRINT_EULER_VORTEX; 57 58 // ------------------------------------------------------ 59 // Create the libCEED context 60 // ------------------------------------------------------ 61 CeedScalar vortex_strength = 5.; // - 62 CeedScalar c_tau = 0.5; // - 63 // c_tau = 0.5 is reported as "optimal" in Hughes et al 2010 64 PetscReal center[3], // m 65 mean_velocity[3] = {1., 1., 0}; // m/s 66 PetscReal domain_min[3], domain_max[3], domain_size[3]; 67 PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 68 for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 69 70 // ------------------------------------------------------ 71 // Create the PETSc context 72 // ------------------------------------------------------ 73 PetscScalar meter = 1e-2; // 1 meter in scaled length units 74 PetscScalar second = 1e-2; // 1 second in scaled time units 75 76 // ------------------------------------------------------ 77 // Command line Options 78 // ------------------------------------------------------ 79 PetscOptionsBegin(comm, NULL, "Options for EULER_VORTEX problem", NULL); 80 // -- Physics 81 PetscCall(PetscOptionsScalar("-vortex_strength", "Strength of Vortex", NULL, vortex_strength, &vortex_strength, NULL)); 82 PetscInt n = problem->dim; 83 PetscBool user_velocity; 84 PetscCall(PetscOptionsRealArray("-mean_velocity", "Background velocity vector", NULL, mean_velocity, &n, &user_velocity)); 85 for (PetscInt i = 0; i < 3; i++) center[i] = .5 * domain_size[i]; 86 n = problem->dim; 87 PetscCall(PetscOptionsRealArray("-center", "Location of vortex center", NULL, center, &n, NULL)); 88 PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 89 PetscCall(PetscOptionsEnum("-euler_test", "Euler test option", NULL, EulerTestTypes, (PetscEnum)(euler_test = EULER_TEST_ISENTROPIC_VORTEX), 90 (PetscEnum *)&euler_test, NULL)); 91 PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 92 PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 93 // -- Units 94 PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 95 meter = fabs(meter); 96 PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 97 second = fabs(second); 98 99 // -- Warnings 100 if (stab == STAB_SUPG && !implicit) { 101 PetscCall(PetscPrintf(comm, "Warning! Use -stab supg only with -implicit\n")); 102 } 103 if (user_velocity && (euler_test == EULER_TEST_1 || euler_test == EULER_TEST_3)) { 104 PetscCall(PetscPrintf(comm, "Warning! Background velocity vector for -euler_test t1 and -euler_test t3 is (0,0,0)\n")); 105 } 106 107 PetscOptionsEnd(); 108 109 // ------------------------------------------------------ 110 // Set up the PETSc context 111 // ------------------------------------------------------ 112 user->units->meter = meter; 113 user->units->second = second; 114 115 // ------------------------------------------------------ 116 // Set up the libCEED context 117 // ------------------------------------------------------ 118 // -- Scale variables to desired units 119 for (PetscInt i = 0; i < 3; i++) { 120 center[i] *= meter; 121 domain_size[i] *= meter; 122 mean_velocity[i] *= (meter / second); 123 } 124 problem->dm_scale = meter; 125 126 // -- QFunction Context 127 user->phys->stab = stab; 128 user->phys->euler_test = euler_test; 129 user->phys->implicit = implicit; 130 user->phys->has_curr_time = has_curr_time; 131 user->phys->has_neumann = has_neumann; 132 euler_ctx->curr_time = 0.; 133 euler_ctx->implicit = implicit; 134 euler_ctx->euler_test = euler_test; 135 euler_ctx->center[0] = center[0]; 136 euler_ctx->center[1] = center[1]; 137 euler_ctx->center[2] = center[2]; 138 euler_ctx->vortex_strength = vortex_strength; 139 euler_ctx->c_tau = c_tau; 140 euler_ctx->mean_velocity[0] = mean_velocity[0]; 141 euler_ctx->mean_velocity[1] = mean_velocity[1]; 142 euler_ctx->mean_velocity[2] = mean_velocity[2]; 143 euler_ctx->stabilization = stab; 144 145 CeedQFunctionContextCreate(user->ceed, &euler_context); 146 CeedQFunctionContextSetData(euler_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*euler_ctx), euler_ctx); 147 CeedQFunctionContextSetDataDestroy(euler_context, CEED_MEM_HOST, FreeContextPetsc); 148 CeedQFunctionContextRegisterDouble(euler_context, "solution time", offsetof(struct EulerContext_, curr_time), 1, "Physical time of the solution"); 149 CeedQFunctionContextReferenceCopy(euler_context, &problem->ics.qfunction_context); 150 CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_vol_rhs.qfunction_context); 151 CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_vol_ifunction.qfunction_context); 152 CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_inflow.qfunction_context); 153 CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_outflow.qfunction_context); 154 PetscFunctionReturn(0); 155 } 156 157 PetscErrorCode PRINT_EULER_VORTEX(ProblemData *problem, AppCtx app_ctx) { 158 MPI_Comm comm = PETSC_COMM_WORLD; 159 EulerContext euler_ctx; 160 161 PetscFunctionBeginUser; 162 CeedQFunctionContextGetData(problem->ics.qfunction_context, CEED_MEM_HOST, &euler_ctx); 163 PetscCall(PetscPrintf(comm, 164 " Problem:\n" 165 " Problem Name : %s\n" 166 " Test Case : %s\n" 167 " Background Velocity : %f,%f,%f\n" 168 " Stabilization : %s\n", 169 app_ctx->problem_name, EulerTestTypes[euler_ctx->euler_test], euler_ctx->mean_velocity[0], euler_ctx->mean_velocity[1], 170 euler_ctx->mean_velocity[2], StabilizationTypes[euler_ctx->stabilization])); 171 172 CeedQFunctionContextRestoreData(problem->ics.qfunction_context, &euler_ctx); 173 PetscFunctionReturn(0); 174 } 175