xref: /honee/problems/shocktube.c (revision ea615d4cc464aa6ad650c06fae6d120cc2465bc4)
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 
14991aef52SJames Wright PetscErrorCode NS_SHOCKTUBE(ProblemData problem, DM dm, void *ctx, SimpleBC bc) {
153636f6a4SJames Wright   SetupContextShock    setup_context;
160c373b74SJames Wright   Honee                honee = *(Honee *)ctx;
170c373b74SJames Wright   MPI_Comm             comm  = honee->comm;
180c373b74SJames Wright   Ceed                 ceed  = honee->ceed;
19af8870a9STimothy Aiken   PetscBool            implicit;
20af8870a9STimothy Aiken   PetscBool            yzb;
21af8870a9STimothy Aiken   PetscInt             stab;
2215a3537eSJed Brown   ShockTubeContext     shocktube_ctx;
23e07531f7SJames Wright   CeedQFunctionContext shocktube_qfctx;
24af8870a9STimothy Aiken 
2515a3537eSJed Brown   PetscFunctionBeginUser;
262b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &setup_context));
272b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &shocktube_ctx));
28af8870a9STimothy Aiken 
29af8870a9STimothy Aiken   // ------------------------------------------------------
30af8870a9STimothy Aiken   //               SET UP SHOCKTUBE
31af8870a9STimothy Aiken   // ------------------------------------------------------
32e07531f7SJames Wright   problem->ics.qf_func_ptr                 = ICsShockTube;
33e07531f7SJames Wright   problem->ics.qf_loc                      = ICsShockTube_loc;
34e07531f7SJames Wright   problem->apply_vol_rhs.qf_func_ptr       = EulerShockTube;
35e07531f7SJames Wright   problem->apply_vol_rhs.qf_loc            = EulerShockTube_loc;
36e07531f7SJames Wright   problem->apply_vol_ifunction.qf_func_ptr = NULL;
37e07531f7SJames Wright   problem->apply_vol_ifunction.qf_loc      = NULL;
3828160fc2SJames Wright   problem->jac_data_size_vol               = 0;
3928160fc2SJames Wright   problem->jac_data_size_sur               = 0;
4058ce1233SJames Wright   problem->compute_exact_solution_error    = PETSC_FALSE;
41af8870a9STimothy Aiken   problem->print_info                      = PRINT_SHOCKTUBE;
42af8870a9STimothy Aiken 
43af8870a9STimothy Aiken   // ------------------------------------------------------
44*ea615d4cSJames Wright   //             Create the QFunction context
45af8870a9STimothy Aiken   // ------------------------------------------------------
46af8870a9STimothy Aiken   // Driver section initial conditions
47af8870a9STimothy Aiken   CeedScalar P_high   = 1.0;  // Pa
48af8870a9STimothy Aiken   CeedScalar rho_high = 1.0;  // kg/m^3
49af8870a9STimothy Aiken   // Driven section initial conditions
50af8870a9STimothy Aiken   CeedScalar P_low   = 0.1;    // Pa
51af8870a9STimothy Aiken   CeedScalar rho_low = 0.125;  // kg/m^3
52af8870a9STimothy Aiken   // Stabilization parameter
53af8870a9STimothy Aiken   CeedScalar c_tau = 0.5;  // -, based on Hughes et al (2010)
54af8870a9STimothy Aiken   // Tuning parameters for the YZB shock capturing
55af8870a9STimothy Aiken   CeedScalar Cyzb = 0.1;  // -, used in approximation of (Na),x
56af8870a9STimothy Aiken   CeedScalar Byzb = 2.0;  // -, 1 for smooth shocks
57af8870a9STimothy Aiken   //                                          2 for sharp shocks
58af8870a9STimothy Aiken   PetscReal domain_min[3], domain_max[3], domain_size[3];
592b916ea7SJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
60493642f1SJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
61af8870a9STimothy Aiken 
62af8870a9STimothy Aiken   // ------------------------------------------------------
63af8870a9STimothy Aiken   //             Create the PETSc context
64af8870a9STimothy Aiken   // ------------------------------------------------------
65af8870a9STimothy Aiken   PetscScalar meter  = 1e-2;  // 1 meter in scaled length units
66af8870a9STimothy Aiken   PetscScalar second = 1e-2;  // 1 second in scaled time units
67af8870a9STimothy Aiken 
68af8870a9STimothy Aiken   // ------------------------------------------------------
69af8870a9STimothy Aiken   //              Command line Options
70af8870a9STimothy Aiken   // ------------------------------------------------------
71af8870a9STimothy Aiken   PetscOptionsBegin(comm, NULL, "Options for SHOCKTUBE problem", NULL);
72af8870a9STimothy Aiken 
73af8870a9STimothy Aiken   // -- Numerical formulation options
742b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
752b916ea7SJeremy L Thompson   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
762b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL));
772b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-yzb", "Use YZB discontinuity capturing", NULL, yzb = PETSC_FALSE, &yzb, NULL));
78af8870a9STimothy Aiken 
79af8870a9STimothy Aiken   // -- Units
802b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
81af8870a9STimothy Aiken   meter = fabs(meter);
822b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
83af8870a9STimothy Aiken   second = fabs(second);
84af8870a9STimothy Aiken 
85af8870a9STimothy Aiken   // -- Warnings
86af8870a9STimothy Aiken   if (stab == STAB_SUPG) {
872b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! -stab supg not implemented for the shocktube problem. \n"));
88af8870a9STimothy Aiken   }
89af8870a9STimothy Aiken   if (yzb && implicit) {
902b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! -yzb only implemented for explicit timestepping. \n"));
91af8870a9STimothy Aiken   }
92af8870a9STimothy Aiken 
93af8870a9STimothy Aiken   PetscOptionsEnd();
94af8870a9STimothy Aiken 
95af8870a9STimothy Aiken   // ------------------------------------------------------
96af8870a9STimothy Aiken   //           Set up the PETSc context
97af8870a9STimothy Aiken   // ------------------------------------------------------
980c373b74SJames Wright   honee->units->meter  = meter;
990c373b74SJames Wright   honee->units->second = second;
100af8870a9STimothy Aiken 
101af8870a9STimothy Aiken   // ------------------------------------------------------
102*ea615d4cSJames Wright   //           Set up the QFunction context
103af8870a9STimothy Aiken   // ------------------------------------------------------
104af8870a9STimothy Aiken   // -- Scale variables to desired units
105493642f1SJames Wright   for (PetscInt i = 0; i < 3; i++) {
106af8870a9STimothy Aiken     domain_size[i] *= meter;
107af8870a9STimothy Aiken     domain_min[i] *= meter;
108af8870a9STimothy Aiken   }
109af8870a9STimothy Aiken   CeedScalar mid_point = 0.5 * (domain_size[0] + domain_min[0]);
110af8870a9STimothy Aiken 
111af8870a9STimothy Aiken   // -- Setup Context
112af8870a9STimothy Aiken   setup_context->mid_point = mid_point;
113af8870a9STimothy Aiken   setup_context->time      = 0.0;
114af8870a9STimothy Aiken   setup_context->P_high    = P_high;
115af8870a9STimothy Aiken   setup_context->rho_high  = rho_high;
116af8870a9STimothy Aiken   setup_context->P_low     = P_low;
117af8870a9STimothy Aiken   setup_context->rho_low   = rho_low;
118af8870a9STimothy Aiken 
119af8870a9STimothy Aiken   // -- QFunction Context
1200c373b74SJames Wright   honee->phys->implicit        = implicit;
12115a3537eSJed Brown   shocktube_ctx->implicit      = implicit;
12215a3537eSJed Brown   shocktube_ctx->stabilization = stab;
12315a3537eSJed Brown   shocktube_ctx->yzb           = yzb;
12415a3537eSJed Brown   shocktube_ctx->Cyzb          = Cyzb;
12515a3537eSJed Brown   shocktube_ctx->Byzb          = Byzb;
12615a3537eSJed Brown   shocktube_ctx->c_tau         = c_tau;
127af8870a9STimothy Aiken 
1280c373b74SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &problem->ics.qfctx));
129e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(problem->ics.qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context));
130e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfctx, CEED_MEM_HOST, FreeContextPetsc));
131af8870a9STimothy Aiken 
1320c373b74SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &shocktube_qfctx));
133e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(shocktube_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*shocktube_ctx), shocktube_ctx));
134e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(shocktube_qfctx, CEED_MEM_HOST, FreeContextPetsc));
135e07531f7SJames Wright   problem->apply_vol_rhs.qfctx = shocktube_qfctx;
136d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
137af8870a9STimothy Aiken }
138af8870a9STimothy Aiken 
1390c373b74SJames Wright PetscErrorCode PRINT_SHOCKTUBE(Honee honee, ProblemData problem, AppCtx app_ctx) {
1400c373b74SJames Wright   MPI_Comm comm = honee->comm;
141af8870a9STimothy Aiken 
14206f41313SJames Wright   PetscFunctionBeginUser;
1432b916ea7SJeremy L Thompson   PetscCall(PetscPrintf(comm,
144af8870a9STimothy Aiken                         "  Problem:\n"
145af8870a9STimothy Aiken                         "    Problem Name                       : %s\n",
1462b916ea7SJeremy L Thompson                         app_ctx->problem_name));
147d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
148af8870a9STimothy Aiken }
149