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