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