xref: /honee/problems/shocktube.c (revision 9b05e62e4a81b29e40cfcd7da7a7fab1171c8eb3)
1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3af8870a9STimothy Aiken 
4af8870a9STimothy Aiken /// @file
5af8870a9STimothy Aiken /// Utility functions for setting up SHOCKTUBE
6af8870a9STimothy Aiken 
7af8870a9STimothy Aiken #include "../qfunctions/shocktube.h"
8af8870a9STimothy Aiken 
9e419654dSJeremy L Thompson #include <ceed.h>
10e419654dSJeremy L Thompson #include <petscdm.h>
11e419654dSJeremy L Thompson 
12149fb536SJames Wright #include <navierstokes.h>
136dd99bedSLeila Ghaffari 
14*9b05e62eSJames Wright static PetscErrorCode PRINT_SHOCKTUBE(Honee honee, ProblemData problem, AppCtx app_ctx) {
15*9b05e62eSJames Wright   MPI_Comm comm = honee->comm;
16*9b05e62eSJames Wright 
17*9b05e62eSJames Wright   PetscFunctionBeginUser;
18*9b05e62eSJames Wright   PetscCall(PetscPrintf(comm,
19*9b05e62eSJames Wright                         "  Problem:\n"
20*9b05e62eSJames Wright                         "    Problem Name                       : %s\n",
21*9b05e62eSJames Wright                         app_ctx->problem_name));
22*9b05e62eSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
23*9b05e62eSJames Wright }
24*9b05e62eSJames Wright 
25d3c60affSJames Wright PetscErrorCode NS_SHOCKTUBE(ProblemData problem, DM dm, void *ctx) {
263636f6a4SJames Wright   SetupContextShock    setup_context;
270c373b74SJames Wright   Honee                honee = *(Honee *)ctx;
280c373b74SJames Wright   MPI_Comm             comm  = honee->comm;
290c373b74SJames Wright   Ceed                 ceed  = honee->ceed;
30af8870a9STimothy Aiken   PetscBool            implicit;
31af8870a9STimothy Aiken   PetscBool            yzb;
32af8870a9STimothy Aiken   PetscInt             stab;
3315a3537eSJed Brown   ShockTubeContext     shocktube_ctx;
34e07531f7SJames Wright   CeedQFunctionContext shocktube_qfctx;
35af8870a9STimothy Aiken 
3615a3537eSJed Brown   PetscFunctionBeginUser;
372b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &setup_context));
382b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &shocktube_ctx));
39af8870a9STimothy Aiken 
40af8870a9STimothy Aiken   // ------------------------------------------------------
41af8870a9STimothy Aiken   //               SET UP SHOCKTUBE
42af8870a9STimothy Aiken   // ------------------------------------------------------
43e07531f7SJames Wright   problem->ics.qf_func_ptr                 = ICsShockTube;
44e07531f7SJames Wright   problem->ics.qf_loc                      = ICsShockTube_loc;
45e07531f7SJames Wright   problem->apply_vol_rhs.qf_func_ptr       = EulerShockTube;
46e07531f7SJames Wright   problem->apply_vol_rhs.qf_loc            = EulerShockTube_loc;
47e07531f7SJames Wright   problem->apply_vol_ifunction.qf_func_ptr = NULL;
48e07531f7SJames Wright   problem->apply_vol_ifunction.qf_loc      = NULL;
491abc2837SJames Wright   problem->num_comps_jac_data              = 0;
5058ce1233SJames Wright   problem->compute_exact_solution_error    = PETSC_FALSE;
51af8870a9STimothy Aiken   problem->print_info                      = PRINT_SHOCKTUBE;
52af8870a9STimothy Aiken 
53af8870a9STimothy Aiken   // ------------------------------------------------------
54ea615d4cSJames Wright   //             Create the QFunction context
55af8870a9STimothy Aiken   // ------------------------------------------------------
56af8870a9STimothy Aiken   // Driver section initial conditions
57af8870a9STimothy Aiken   CeedScalar P_high   = 1.0;  // Pa
58af8870a9STimothy Aiken   CeedScalar rho_high = 1.0;  // kg/m^3
59af8870a9STimothy Aiken   // Driven section initial conditions
60af8870a9STimothy Aiken   CeedScalar P_low   = 0.1;    // Pa
61af8870a9STimothy Aiken   CeedScalar rho_low = 0.125;  // kg/m^3
62af8870a9STimothy Aiken   // Stabilization parameter
63af8870a9STimothy Aiken   CeedScalar c_tau = 0.5;  // -, based on Hughes et al (2010)
64af8870a9STimothy Aiken   // Tuning parameters for the YZB shock capturing
65af8870a9STimothy Aiken   CeedScalar Cyzb = 0.1;  // -, used in approximation of (Na),x
66af8870a9STimothy Aiken   CeedScalar Byzb = 2.0;  // -, 1 for smooth shocks
67af8870a9STimothy Aiken   //                                          2 for sharp shocks
68af8870a9STimothy Aiken   PetscReal domain_min[3], domain_max[3], domain_size[3];
692b916ea7SJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
70493642f1SJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
71af8870a9STimothy Aiken 
72af8870a9STimothy Aiken   // ------------------------------------------------------
73af8870a9STimothy Aiken   //             Create the PETSc context
74af8870a9STimothy Aiken   // ------------------------------------------------------
75af8870a9STimothy Aiken   PetscScalar meter  = 1e-2;  // 1 meter in scaled length units
76af8870a9STimothy Aiken   PetscScalar second = 1e-2;  // 1 second in scaled time units
77af8870a9STimothy Aiken 
78af8870a9STimothy Aiken   // ------------------------------------------------------
79af8870a9STimothy Aiken   //              Command line Options
80af8870a9STimothy Aiken   // ------------------------------------------------------
81af8870a9STimothy Aiken   PetscOptionsBegin(comm, NULL, "Options for SHOCKTUBE problem", NULL);
82af8870a9STimothy Aiken 
83af8870a9STimothy Aiken   // -- Numerical formulation options
842b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
852b916ea7SJeremy L Thompson   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
862b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL));
872b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-yzb", "Use YZB discontinuity capturing", NULL, yzb = PETSC_FALSE, &yzb, NULL));
88af8870a9STimothy Aiken 
89af8870a9STimothy Aiken   // -- Units
902b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
91af8870a9STimothy Aiken   meter = fabs(meter);
922b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
93af8870a9STimothy Aiken   second = fabs(second);
94af8870a9STimothy Aiken 
95af8870a9STimothy Aiken   // -- Warnings
96af8870a9STimothy Aiken   if (stab == STAB_SUPG) {
972b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! -stab supg not implemented for the shocktube problem. \n"));
98af8870a9STimothy Aiken   }
99af8870a9STimothy Aiken   if (yzb && implicit) {
1002b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! -yzb only implemented for explicit timestepping. \n"));
101af8870a9STimothy Aiken   }
102af8870a9STimothy Aiken 
103af8870a9STimothy Aiken   PetscOptionsEnd();
104af8870a9STimothy Aiken 
105af8870a9STimothy Aiken   // ------------------------------------------------------
106af8870a9STimothy Aiken   //           Set up the PETSc context
107af8870a9STimothy Aiken   // ------------------------------------------------------
1080c373b74SJames Wright   honee->units->meter  = meter;
1090c373b74SJames Wright   honee->units->second = second;
110af8870a9STimothy Aiken 
111af8870a9STimothy Aiken   // ------------------------------------------------------
112ea615d4cSJames Wright   //           Set up the QFunction context
113af8870a9STimothy Aiken   // ------------------------------------------------------
114af8870a9STimothy Aiken   // -- Scale variables to desired units
115493642f1SJames Wright   for (PetscInt i = 0; i < 3; i++) {
116af8870a9STimothy Aiken     domain_size[i] *= meter;
117af8870a9STimothy Aiken     domain_min[i] *= meter;
118af8870a9STimothy Aiken   }
119af8870a9STimothy Aiken   CeedScalar mid_point = 0.5 * (domain_size[0] + domain_min[0]);
120af8870a9STimothy Aiken 
121af8870a9STimothy Aiken   // -- Setup Context
122af8870a9STimothy Aiken   setup_context->mid_point = mid_point;
123af8870a9STimothy Aiken   setup_context->time      = 0.0;
124af8870a9STimothy Aiken   setup_context->P_high    = P_high;
125af8870a9STimothy Aiken   setup_context->rho_high  = rho_high;
126af8870a9STimothy Aiken   setup_context->P_low     = P_low;
127af8870a9STimothy Aiken   setup_context->rho_low   = rho_low;
128af8870a9STimothy Aiken 
129af8870a9STimothy Aiken   // -- QFunction Context
1300c373b74SJames Wright   honee->phys->implicit        = implicit;
13115a3537eSJed Brown   shocktube_ctx->implicit      = implicit;
13215a3537eSJed Brown   shocktube_ctx->stabilization = stab;
13315a3537eSJed Brown   shocktube_ctx->yzb           = yzb;
13415a3537eSJed Brown   shocktube_ctx->Cyzb          = Cyzb;
13515a3537eSJed Brown   shocktube_ctx->Byzb          = Byzb;
13615a3537eSJed Brown   shocktube_ctx->c_tau         = c_tau;
137af8870a9STimothy Aiken 
1380c373b74SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &problem->ics.qfctx));
139e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(problem->ics.qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context));
140e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfctx, CEED_MEM_HOST, FreeContextPetsc));
141af8870a9STimothy Aiken 
1420c373b74SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &shocktube_qfctx));
143e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(shocktube_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*shocktube_ctx), shocktube_ctx));
144e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(shocktube_qfctx, CEED_MEM_HOST, FreeContextPetsc));
145e07531f7SJames Wright   problem->apply_vol_rhs.qfctx = shocktube_qfctx;
146d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
147af8870a9STimothy Aiken }
148