xref: /honee/problems/shocktube.c (revision 2d898fa6800fcddb9cd7cb214e287e8576fd592f)
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 
149b05e62eSJames Wright static PetscErrorCode PRINT_SHOCKTUBE(Honee honee, ProblemData problem, AppCtx app_ctx) {
159b05e62eSJames Wright   MPI_Comm comm = honee->comm;
169b05e62eSJames Wright 
179b05e62eSJames Wright   PetscFunctionBeginUser;
189b05e62eSJames Wright   PetscCall(PetscPrintf(comm,
199b05e62eSJames Wright                         "  Problem:\n"
209b05e62eSJames Wright                         "    Problem Name                       : %s\n",
219b05e62eSJames Wright                         app_ctx->problem_name));
229b05e62eSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
239b05e62eSJames Wright }
249b05e62eSJames Wright 
25f27c2204SJames Wright static const char *const component_names[] = {"Density", "MomentumX", "MomentumY", "MomentumZ", "TotalEnergy"};
26f27c2204SJames Wright 
27d3c60affSJames Wright PetscErrorCode NS_SHOCKTUBE(ProblemData problem, DM dm, void *ctx) {
283636f6a4SJames Wright   SetupContextShock    setup_context;
290c373b74SJames Wright   Honee                honee = *(Honee *)ctx;
300c373b74SJames Wright   MPI_Comm             comm  = honee->comm;
310c373b74SJames Wright   Ceed                 ceed  = honee->ceed;
32af8870a9STimothy Aiken   PetscBool            implicit;
33af8870a9STimothy Aiken   PetscBool            yzb;
34af8870a9STimothy Aiken   PetscInt             stab;
3515a3537eSJed Brown   ShockTubeContext     shocktube_ctx;
36e07531f7SJames Wright   CeedQFunctionContext shocktube_qfctx;
37af8870a9STimothy Aiken 
3815a3537eSJed Brown   PetscFunctionBeginUser;
39*2d898fa6SJames Wright   PetscCall(PetscNew(&setup_context));
40*2d898fa6SJames Wright   PetscCall(PetscNew(&shocktube_ctx));
41af8870a9STimothy Aiken 
42af8870a9STimothy Aiken   // ------------------------------------------------------
43af8870a9STimothy Aiken   //               SET UP SHOCKTUBE
44af8870a9STimothy Aiken   // ------------------------------------------------------
45e07531f7SJames Wright   problem->ics.qf_func_ptr                 = ICsShockTube;
46e07531f7SJames Wright   problem->ics.qf_loc                      = ICsShockTube_loc;
47e07531f7SJames Wright   problem->apply_vol_rhs.qf_func_ptr       = EulerShockTube;
48e07531f7SJames Wright   problem->apply_vol_rhs.qf_loc            = EulerShockTube_loc;
49e07531f7SJames Wright   problem->apply_vol_ifunction.qf_func_ptr = NULL;
50e07531f7SJames Wright   problem->apply_vol_ifunction.qf_loc      = NULL;
511abc2837SJames Wright   problem->num_comps_jac_data              = 0;
5258ce1233SJames Wright   problem->compute_exact_solution_error    = PETSC_FALSE;
53af8870a9STimothy Aiken   problem->print_info                      = PRINT_SHOCKTUBE;
54af8870a9STimothy Aiken 
55f27c2204SJames Wright   problem->num_components = 5;
56f27c2204SJames Wright   PetscCall(PetscMalloc1(problem->num_components, &problem->component_names));
57f27c2204SJames Wright   for (PetscInt i = 0; i < 5; i++) PetscCall(PetscStrallocpy(component_names[i], &problem->component_names[i]));
58f27c2204SJames Wright 
59af8870a9STimothy Aiken   // ------------------------------------------------------
60ea615d4cSJames Wright   //             Create the QFunction context
61af8870a9STimothy Aiken   // ------------------------------------------------------
62af8870a9STimothy Aiken   // Driver section initial conditions
63af8870a9STimothy Aiken   CeedScalar P_high   = 1.0;  // Pa
64af8870a9STimothy Aiken   CeedScalar rho_high = 1.0;  // kg/m^3
65af8870a9STimothy Aiken   // Driven section initial conditions
66af8870a9STimothy Aiken   CeedScalar P_low   = 0.1;    // Pa
67af8870a9STimothy Aiken   CeedScalar rho_low = 0.125;  // kg/m^3
68af8870a9STimothy Aiken   // Stabilization parameter
69af8870a9STimothy Aiken   CeedScalar c_tau = 0.5;  // -, based on Hughes et al (2010)
70af8870a9STimothy Aiken   // Tuning parameters for the YZB shock capturing
71af8870a9STimothy Aiken   CeedScalar Cyzb = 0.1;  // -, used in approximation of (Na),x
72af8870a9STimothy Aiken   CeedScalar Byzb = 2.0;  // -, 1 for smooth shocks
73af8870a9STimothy Aiken   //                                          2 for sharp shocks
74af8870a9STimothy Aiken   PetscReal domain_min[3], domain_max[3], domain_size[3];
752b916ea7SJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
76493642f1SJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
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   // -- Warnings
90af8870a9STimothy Aiken   if (stab == STAB_SUPG) {
912b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! -stab supg not implemented for the shocktube problem. \n"));
92af8870a9STimothy Aiken   }
93af8870a9STimothy Aiken   if (yzb && implicit) {
942b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! -yzb only implemented for explicit timestepping. \n"));
95af8870a9STimothy Aiken   }
96af8870a9STimothy Aiken 
97af8870a9STimothy Aiken   PetscOptionsEnd();
98af8870a9STimothy Aiken 
99af8870a9STimothy Aiken   // ------------------------------------------------------
100ea615d4cSJames Wright   //           Set up the QFunction context
101af8870a9STimothy Aiken   // ------------------------------------------------------
102af8870a9STimothy Aiken   // -- Scale variables to desired units
103c9f37605SMohammed Amin   Units units = honee->units;
104493642f1SJames Wright   for (PetscInt i = 0; i < 3; i++) {
105c9f37605SMohammed Amin     domain_size[i] *= units->meter;
106c9f37605SMohammed Amin     domain_min[i] *= units->meter;
107af8870a9STimothy Aiken   }
108af8870a9STimothy Aiken   CeedScalar mid_point = 0.5 * (domain_size[0] + domain_min[0]);
109af8870a9STimothy Aiken 
110af8870a9STimothy Aiken   // -- Setup Context
111af8870a9STimothy Aiken   setup_context->mid_point = mid_point;
112af8870a9STimothy Aiken   setup_context->time      = 0.0;
113af8870a9STimothy Aiken   setup_context->P_high    = P_high;
114af8870a9STimothy Aiken   setup_context->rho_high  = rho_high;
115af8870a9STimothy Aiken   setup_context->P_low     = P_low;
116af8870a9STimothy Aiken   setup_context->rho_low   = rho_low;
117af8870a9STimothy Aiken 
118af8870a9STimothy Aiken   // -- QFunction Context
1190c373b74SJames Wright   honee->phys->implicit        = implicit;
12015a3537eSJed Brown   shocktube_ctx->implicit      = implicit;
12115a3537eSJed Brown   shocktube_ctx->stabilization = stab;
12215a3537eSJed Brown   shocktube_ctx->yzb           = yzb;
12315a3537eSJed Brown   shocktube_ctx->Cyzb          = Cyzb;
12415a3537eSJed Brown   shocktube_ctx->Byzb          = Byzb;
12515a3537eSJed Brown   shocktube_ctx->c_tau         = c_tau;
126af8870a9STimothy Aiken 
1270c373b74SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &problem->ics.qfctx));
128e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(problem->ics.qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context));
129e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfctx, CEED_MEM_HOST, FreeContextPetsc));
130af8870a9STimothy Aiken 
1310c373b74SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &shocktube_qfctx));
132e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(shocktube_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*shocktube_ctx), shocktube_ctx));
133e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(shocktube_qfctx, CEED_MEM_HOST, FreeContextPetsc));
134e07531f7SJames Wright   problem->apply_vol_rhs.qfctx = shocktube_qfctx;
135d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
136af8870a9STimothy Aiken }
137