xref: /libCEED/examples/fluids/problems/shocktube.c (revision 8ec64e9ae9d5df169dba8c8ee61d8ec8907b8f80)
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 SHOCKTUBE
10 
11 #include "../qfunctions/shocktube.h"
12 
13 #include "../navierstokes.h"
14 #include "../qfunctions/setupgeo.h"
15 
16 PetscErrorCode NS_SHOCKTUBE(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) {
17   SetupContextShock    setup_context;
18   User                 user = *(User *)ctx;
19   MPI_Comm             comm = PETSC_COMM_WORLD;
20   PetscBool            implicit;
21   PetscBool            yzb;
22   PetscInt             stab;
23   PetscBool            has_curr_time = PETSC_FALSE;
24   ShockTubeContext     shocktube_ctx;
25   CeedQFunctionContext shocktube_context;
26 
27   PetscFunctionBeginUser;
28   PetscCall(PetscCalloc1(1, &setup_context));
29   PetscCall(PetscCalloc1(1, &shocktube_ctx));
30 
31   // ------------------------------------------------------
32   //               SET UP SHOCKTUBE
33   // ------------------------------------------------------
34   problem->dim                               = 3;
35   problem->q_data_size_vol                   = 10;
36   problem->q_data_size_sur                   = 4;
37   problem->setup_vol.qfunction               = Setup;
38   problem->setup_vol.qfunction_loc           = Setup_loc;
39   problem->setup_sur.qfunction               = SetupBoundary;
40   problem->setup_sur.qfunction_loc           = SetupBoundary_loc;
41   problem->ics.qfunction                     = ICsShockTube;
42   problem->ics.qfunction_loc                 = ICsShockTube_loc;
43   problem->apply_vol_rhs.qfunction           = EulerShockTube;
44   problem->apply_vol_rhs.qfunction_loc       = EulerShockTube_loc;
45   problem->apply_vol_ifunction.qfunction     = NULL;
46   problem->apply_vol_ifunction.qfunction_loc = NULL;
47   problem->bc                                = Exact_ShockTube;
48   problem->bc_ctx                            = setup_context;
49   problem->non_zero_time                     = PETSC_FALSE;
50   problem->print_info                        = PRINT_SHOCKTUBE;
51 
52   // ------------------------------------------------------
53   //             Create the libCEED context
54   // ------------------------------------------------------
55   // Driver section initial conditions
56   CeedScalar P_high   = 1.0;  // Pa
57   CeedScalar rho_high = 1.0;  // kg/m^3
58   // Driven section initial conditions
59   CeedScalar P_low   = 0.1;    // Pa
60   CeedScalar rho_low = 0.125;  // kg/m^3
61   // Stabilization parameter
62   CeedScalar c_tau = 0.5;  // -, based on Hughes et al (2010)
63   // Tuning parameters for the YZB shock capturing
64   CeedScalar Cyzb = 0.1;  // -, used in approximation of (Na),x
65   CeedScalar Byzb = 2.0;  // -, 1 for smooth shocks
66   //                                          2 for sharp shocks
67   PetscReal domain_min[3], domain_max[3], domain_size[3];
68   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
69   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
70 
71   // ------------------------------------------------------
72   //             Create the PETSc context
73   // ------------------------------------------------------
74   PetscScalar meter  = 1e-2;  // 1 meter in scaled length units
75   PetscScalar second = 1e-2;  // 1 second in scaled time units
76 
77   // ------------------------------------------------------
78   //              Command line Options
79   // ------------------------------------------------------
80   PetscOptionsBegin(comm, NULL, "Options for SHOCKTUBE problem", NULL);
81 
82   // -- Numerical formulation options
83   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
84   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
85   PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL));
86   PetscCall(PetscOptionsBool("-yzb", "Use YZB discontinuity capturing", NULL, yzb = PETSC_FALSE, &yzb, NULL));
87 
88   // -- Units
89   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
90   meter = fabs(meter);
91   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
92   second = fabs(second);
93 
94   // -- Warnings
95   if (stab == STAB_SUPG) {
96     PetscCall(PetscPrintf(comm, "Warning! -stab supg not implemented for the shocktube problem. \n"));
97   }
98   if (yzb && implicit) {
99     PetscCall(PetscPrintf(comm, "Warning! -yzb only implemented for explicit timestepping. \n"));
100   }
101 
102   PetscOptionsEnd();
103 
104   // ------------------------------------------------------
105   //           Set up the PETSc context
106   // ------------------------------------------------------
107   user->units->meter  = meter;
108   user->units->second = second;
109 
110   // ------------------------------------------------------
111   //           Set up the libCEED context
112   // ------------------------------------------------------
113   // -- Scale variables to desired units
114   for (PetscInt i = 0; i < 3; i++) {
115     domain_size[i] *= meter;
116     domain_min[i] *= meter;
117   }
118   problem->dm_scale    = meter;
119   CeedScalar mid_point = 0.5 * (domain_size[0] + domain_min[0]);
120 
121   // -- Setup Context
122   setup_context->mid_point = mid_point;
123   setup_context->time      = 0.0;
124   setup_context->P_high    = P_high;
125   setup_context->rho_high  = rho_high;
126   setup_context->P_low     = P_low;
127   setup_context->rho_low   = rho_low;
128 
129   // -- QFunction Context
130   user->phys->implicit         = implicit;
131   user->phys->has_curr_time    = has_curr_time;
132   shocktube_ctx->implicit      = implicit;
133   shocktube_ctx->stabilization = stab;
134   shocktube_ctx->yzb           = yzb;
135   shocktube_ctx->Cyzb          = Cyzb;
136   shocktube_ctx->Byzb          = Byzb;
137   shocktube_ctx->c_tau         = c_tau;
138 
139   CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context);
140   CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context);
141   CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc);
142 
143   CeedQFunctionContextCreate(user->ceed, &shocktube_context);
144   CeedQFunctionContextSetData(shocktube_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*shocktube_ctx), shocktube_ctx);
145   CeedQFunctionContextSetDataDestroy(shocktube_context, CEED_MEM_HOST, FreeContextPetsc);
146   problem->apply_vol_rhs.qfunction_context = shocktube_context;
147   PetscFunctionReturn(0);
148 }
149 
150 PetscErrorCode PRINT_SHOCKTUBE(ProblemData *problem, AppCtx app_ctx) {
151   MPI_Comm comm = PETSC_COMM_WORLD;
152   PetscFunctionBeginUser;
153 
154   PetscCall(PetscPrintf(comm,
155                         "  Problem:\n"
156                         "    Problem Name                       : %s\n",
157                         app_ctx->problem_name));
158 
159   PetscFunctionReturn(0);
160 }
161