1 // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2 // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3 // reserved. See files LICENSE and NOTICE for details. 4 // 5 // This file is part of CEED, a collection of benchmarks, miniapps, software 6 // libraries and APIs for efficient high-order finite element and spectral 7 // element discretizations for exascale applications. For more information and 8 // source code availability see http://github.com/ceed. 9 // 10 // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11 // a collaborative effort of two U.S. Department of Energy organizations (Office 12 // of Science and the National Nuclear Security Administration) responsible for 13 // the planning and preparation of a capable exascale ecosystem, including 14 // software, applications, hardware, advanced system engineering and early 15 // testbed platforms, in support of the nation's exascale computing imperative. 16 17 /// @file 18 /// Utility functions for setting up SHOCKTUBE 19 20 #include "../navierstokes.h" 21 #include "../qfunctions/setupgeo.h" 22 #include "../qfunctions/shocktube.h" 23 24 PetscErrorCode NS_SHOCKTUBE(ProblemData *problem, DM dm, void *setup_ctx, 25 void *ctx) { 26 SetupContext setup_context = *(SetupContext *)setup_ctx; 27 User user = *(User *)ctx; 28 MPI_Comm comm = PETSC_COMM_WORLD; 29 PetscBool implicit; 30 PetscBool yzb; 31 PetscInt stab; 32 PetscBool has_curr_time = PETSC_FALSE; 33 PetscInt ierr; 34 ShockTubeContext shocktube_ctx; 35 CeedQFunctionContext shocktube_context; 36 37 38 PetscFunctionBeginUser; 39 ierr = PetscCalloc1(1, &shocktube_ctx); CHKERRQ(ierr); 40 41 // ------------------------------------------------------ 42 // SET UP SHOCKTUBE 43 // ------------------------------------------------------ 44 problem->dim = 3; 45 problem->q_data_size_vol = 10; 46 problem->q_data_size_sur = 4; 47 problem->setup_vol.qfunction = Setup; 48 problem->setup_vol.qfunction_loc = Setup_loc; 49 problem->setup_sur.qfunction = SetupBoundary; 50 problem->setup_sur.qfunction_loc = SetupBoundary_loc; 51 problem->ics.qfunction = ICsShockTube; 52 problem->ics.qfunction_loc = ICsShockTube_loc; 53 problem->apply_vol_rhs.qfunction = EulerShockTube; 54 problem->apply_vol_rhs.qfunction_loc = EulerShockTube_loc; 55 problem->apply_vol_ifunction.qfunction = NULL; 56 problem->apply_vol_ifunction.qfunction_loc = NULL; 57 problem->bc = Exact_ShockTube; 58 problem->non_zero_time = PETSC_FALSE; 59 problem->print_info = PRINT_SHOCKTUBE; 60 61 // ------------------------------------------------------ 62 // Create the libCEED context 63 // ------------------------------------------------------ 64 // Driver section initial conditions 65 CeedScalar P_high = 1.0; // Pa 66 CeedScalar rho_high = 1.0; // kg/m^3 67 // Driven section initial conditions 68 CeedScalar P_low = 0.1; // Pa 69 CeedScalar rho_low = 0.125; // kg/m^3 70 // Stabilization parameter 71 CeedScalar c_tau = 0.5; // -, based on Hughes et al (2010) 72 // Tuning parameters for the YZB shock capturing 73 CeedScalar Cyzb = 0.1; // -, used in approximation of (Na),x 74 CeedScalar Byzb = 2.0; // -, 1 for smooth shocks 75 // 2 for sharp shocks 76 PetscReal domain_min[3], domain_max[3], domain_size[3]; 77 ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 78 for (int i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 79 80 // ------------------------------------------------------ 81 // Create the PETSc context 82 // ------------------------------------------------------ 83 PetscScalar meter = 1e-2; // 1 meter in scaled length units 84 PetscScalar second = 1e-2; // 1 second in scaled time units 85 86 // ------------------------------------------------------ 87 // Command line Options 88 // ------------------------------------------------------ 89 PetscOptionsBegin(comm, NULL, "Options for SHOCKTUBE problem", NULL); 90 91 // -- Numerical formulation options 92 ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", 93 NULL, implicit=PETSC_FALSE, &implicit, NULL); CHKERRQ(ierr); 94 ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL, 95 StabilizationTypes, (PetscEnum)(stab = STAB_NONE), 96 (PetscEnum *)&stab, NULL); CHKERRQ(ierr); 97 ierr = PetscOptionsScalar("-c_tau", "Stabilization constant", 98 NULL, c_tau, &c_tau, NULL); CHKERRQ(ierr); 99 ierr = PetscOptionsBool("-yzb", "Use YZB discontinuity capturing", 100 NULL, yzb=PETSC_FALSE, &yzb, NULL); CHKERRQ(ierr); 101 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) { 112 ierr = PetscPrintf(comm, 113 "Warning! -stab supg not implemented for the shocktube problem. \n"); 114 CHKERRQ(ierr); 115 } 116 if (yzb && implicit) { 117 ierr = PetscPrintf(comm, 118 "Warning! -yzb only implemented for explicit timestepping. \n"); 119 CHKERRQ(ierr); 120 } 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 domain_size[i] *= meter; 137 domain_min[i] *= meter; 138 } 139 problem->dm_scale = meter; 140 CeedScalar mid_point = 0.5*(domain_size[0]+domain_min[0]); 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->mid_point = mid_point; 147 setup_context->time = 0.0; 148 setup_context->P_high = P_high; 149 setup_context->rho_high = rho_high; 150 setup_context->P_low = P_low; 151 setup_context->rho_low = rho_low; 152 153 // -- QFunction Context 154 user->phys->implicit = implicit; 155 user->phys->has_curr_time = has_curr_time; 156 shocktube_ctx->implicit = implicit; 157 shocktube_ctx->stabilization = stab; 158 shocktube_ctx->yzb = yzb; 159 shocktube_ctx->Cyzb = Cyzb; 160 shocktube_ctx->Byzb = Byzb; 161 shocktube_ctx->c_tau = c_tau; 162 163 CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context); 164 CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, 165 CEED_USE_POINTER, sizeof(*setup_context), setup_context); 166 167 CeedQFunctionContextCreate(user->ceed, &shocktube_context); 168 CeedQFunctionContextSetData(shocktube_context, CEED_MEM_HOST, 169 CEED_USE_POINTER, 170 sizeof(*shocktube_ctx), shocktube_ctx); 171 CeedQFunctionContextSetDataDestroy(shocktube_context, CEED_MEM_HOST, 172 FreeContextPetsc); 173 problem->apply_vol_rhs.qfunction_context = shocktube_context; 174 PetscFunctionReturn(0); 175 } 176 177 PetscErrorCode PRINT_SHOCKTUBE(ProblemData *problem, SetupContext setup_ctx, 178 AppCtx app_ctx) { 179 MPI_Comm comm = PETSC_COMM_WORLD; 180 PetscErrorCode ierr; 181 PetscFunctionBeginUser; 182 183 ierr = PetscPrintf(comm, 184 " Problem:\n" 185 " Problem Name : %s\n", 186 app_ctx->problem_name); CHKERRQ(ierr); 187 188 PetscFunctionReturn(0); 189 } 190