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