xref: /honee/problems/blasius.c (revision 337840fc9f8b699b98947f09e30443f8c0c7f962)
1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3bb8a0c61SJames Wright 
4bb8a0c61SJames Wright /// @file
5bb8a0c61SJames Wright /// Utility functions for setting up Blasius Boundary Layer
6bb8a0c61SJames Wright 
7bb8a0c61SJames Wright #include "../qfunctions/blasius.h"
82b916ea7SJeremy L Thompson 
9e419654dSJeremy L Thompson #include <ceed.h>
10e419654dSJeremy L Thompson #include <petscdm.h>
11e419654dSJeremy L Thompson #include <petscdt.h>
12e419654dSJeremy L Thompson 
138d78d7c8SJames Wright #include <differential_filter.h>
14149fb536SJames Wright #include <navierstokes.h>
15493642f1SJames Wright #include "stg_shur14.h"
16bb8a0c61SJames Wright 
CompressibleBlasiusResidual(SNES snes,Vec X,Vec R,void * ctx)17e0d1a4dfSLeila Ghaffari PetscErrorCode CompressibleBlasiusResidual(SNES snes, Vec X, Vec R, void *ctx) {
18e0d1a4dfSLeila Ghaffari   const BlasiusContext  blasius = (BlasiusContext)ctx;
19e0d1a4dfSLeila Ghaffari   const PetscScalar    *Tf, *Th;  // Chebyshev coefficients
20e0d1a4dfSLeila Ghaffari   PetscScalar          *r, f[4], h[4];
21e0d1a4dfSLeila Ghaffari   PetscInt              N       = blasius->n_cheb;
22fcb2c22aSJames Wright   State                 S_infty = blasius->S_infty;
2364667825SJames Wright   CeedScalar            U_infty = Norm3(S_infty.Y.velocity);
24cde3d787SJames Wright   NewtonianIGProperties gas     = blasius->newt_ctx.gas;
2506f41313SJames Wright 
2606f41313SJames Wright   PetscFunctionBeginUser;
27cde3d787SJames Wright   PetscScalar Ma = Mach(gas, S_infty.Y.temperature, U_infty), Pr = Prandtl(gas), gamma = HeatCapacityRatio(gas);
2806f41313SJames Wright 
29e0d1a4dfSLeila Ghaffari   PetscCall(VecGetArrayRead(X, &Tf));
30e0d1a4dfSLeila Ghaffari   Th = Tf + N;
31e0d1a4dfSLeila Ghaffari   PetscCall(VecGetArray(R, &r));
32e0d1a4dfSLeila Ghaffari 
33e0d1a4dfSLeila Ghaffari   // Left boundary conditions f = f' = 0
34e0d1a4dfSLeila Ghaffari   ChebyshevEval(N, Tf, -1., blasius->eta_max, f);
35e0d1a4dfSLeila Ghaffari   r[0] = f[0];
36e0d1a4dfSLeila Ghaffari   r[1] = f[1];
37e0d1a4dfSLeila Ghaffari 
38e0d1a4dfSLeila Ghaffari   // f - right end boundary condition
39e0d1a4dfSLeila Ghaffari   ChebyshevEval(N, Tf, 1., blasius->eta_max, f);
40e0d1a4dfSLeila Ghaffari   r[2] = f[1] - 1.;
41e0d1a4dfSLeila Ghaffari 
42e0d1a4dfSLeila Ghaffari   for (int i = 0; i < N - 3; i++) {
43e0d1a4dfSLeila Ghaffari     ChebyshevEval(N, Tf, blasius->X[i], blasius->eta_max, f);
44e0d1a4dfSLeila Ghaffari     ChebyshevEval(N - 1, Th, blasius->X[i], blasius->eta_max, h);
4504e40bb6SJeremy L Thompson     // mu and rho generally depend on h.
4604e40bb6SJeremy L Thompson     // We naively assume constant mu.
470d850f2eSLeila Ghaffari     // For an ideal gas at constant pressure, density is inversely proportional to enthalpy.
480d850f2eSLeila Ghaffari     // The *_tilde values are *relative* to their freestream values, and we proved first derivatives here.
490d850f2eSLeila Ghaffari     const PetscScalar mu_tilde[2]     = {1, 0};
500d850f2eSLeila Ghaffari     const PetscScalar rho_tilde[2]    = {1 / h[0], -h[1] / PetscSqr(h[0])};
510d850f2eSLeila Ghaffari     const PetscScalar mu_rho_tilde[2] = {
520d850f2eSLeila Ghaffari         mu_tilde[0] * rho_tilde[0],
530d850f2eSLeila Ghaffari         mu_tilde[1] * rho_tilde[0] + mu_tilde[0] * rho_tilde[1],
540d850f2eSLeila Ghaffari     };
550d850f2eSLeila Ghaffari     r[3 + i]     = 2 * (mu_rho_tilde[0] * f[3] + mu_rho_tilde[1] * f[2]) + f[2] * f[0];
562b916ea7SJeremy L Thompson     r[N + 2 + i] = (mu_rho_tilde[0] * h[2] + mu_rho_tilde[1] * h[1]) + Pr * f[0] * h[1] + Pr * (gamma - 1) * mu_rho_tilde[0] * PetscSqr(Ma * f[2]);
57e0d1a4dfSLeila Ghaffari   }
58e0d1a4dfSLeila Ghaffari 
59e0d1a4dfSLeila Ghaffari   // h - left end boundary condition
60e0d1a4dfSLeila Ghaffari   ChebyshevEval(N - 1, Th, -1., blasius->eta_max, h);
61fcb2c22aSJames Wright   r[N] = h[0] - blasius->T_wall / S_infty.Y.temperature;
62e0d1a4dfSLeila Ghaffari 
63e0d1a4dfSLeila Ghaffari   // h - right end boundary condition
64e0d1a4dfSLeila Ghaffari   ChebyshevEval(N - 1, Th, 1., blasius->eta_max, h);
65e0d1a4dfSLeila Ghaffari   r[N + 1] = h[0] - 1.;
66e0d1a4dfSLeila Ghaffari 
67e0d1a4dfSLeila Ghaffari   // Restore vectors
68e0d1a4dfSLeila Ghaffari   PetscCall(VecRestoreArrayRead(X, &Tf));
69e0d1a4dfSLeila Ghaffari   PetscCall(VecRestoreArray(R, &r));
70d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
71e0d1a4dfSLeila Ghaffari }
72e0d1a4dfSLeila Ghaffari 
ComputeChebyshevCoefficients(BlasiusContext blasius)73e0d1a4dfSLeila Ghaffari PetscErrorCode ComputeChebyshevCoefficients(BlasiusContext blasius) {
74e0d1a4dfSLeila Ghaffari   SNES                snes;
75e0d1a4dfSLeila Ghaffari   Vec                 sol, res;
76e0d1a4dfSLeila Ghaffari   PetscReal          *w;
77e0d1a4dfSLeila Ghaffari   PetscInt            N = blasius->n_cheb;
780d850f2eSLeila Ghaffari   SNESConvergedReason reason;
79e0d1a4dfSLeila Ghaffari   const PetscScalar  *cheb_coefs;
800d850f2eSLeila Ghaffari 
8106f41313SJames Wright   PetscFunctionBeginUser;
820d850f2eSLeila Ghaffari   // Allocate memory
83e0d1a4dfSLeila Ghaffari   PetscCall(PetscMalloc2(N - 3, &blasius->X, N - 3, &w));
84e0d1a4dfSLeila Ghaffari   PetscCall(PetscDTGaussQuadrature(N - 3, -1., 1., blasius->X, w));
850d850f2eSLeila Ghaffari 
860d850f2eSLeila Ghaffari   // Snes solve
87e0d1a4dfSLeila Ghaffari   PetscCall(SNESCreate(PETSC_COMM_SELF, &snes));
88e0d1a4dfSLeila Ghaffari   PetscCall(VecCreate(PETSC_COMM_SELF, &sol));
89e0d1a4dfSLeila Ghaffari   PetscCall(VecSetSizes(sol, PETSC_DECIDE, 2 * N - 1));
90e0d1a4dfSLeila Ghaffari   PetscCall(VecSetFromOptions(sol));
910d850f2eSLeila Ghaffari   // Constant relative enthalpy 1 as initial guess
920d850f2eSLeila Ghaffari   PetscCall(VecSetValue(sol, N, 1., INSERT_VALUES));
93e0d1a4dfSLeila Ghaffari   PetscCall(VecDuplicate(sol, &res));
94e0d1a4dfSLeila Ghaffari   PetscCall(SNESSetFunction(snes, res, CompressibleBlasiusResidual, blasius));
950d850f2eSLeila Ghaffari   PetscCall(SNESSetOptionsPrefix(snes, "chebyshev_"));
96e0d1a4dfSLeila Ghaffari   PetscCall(SNESSetFromOptions(snes));
97e0d1a4dfSLeila Ghaffari   PetscCall(SNESSolve(snes, NULL, sol));
980d850f2eSLeila Ghaffari   PetscCall(SNESGetConvergedReason(snes, &reason));
9914bd2a07SJames Wright   PetscCheck(reason >= 0, PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "The Chebyshev solve failed.");
1000d850f2eSLeila Ghaffari 
1010d850f2eSLeila Ghaffari   // Assign Chebyshev coefficients
102e0d1a4dfSLeila Ghaffari   PetscCall(VecGetArrayRead(sol, &cheb_coefs));
103e0d1a4dfSLeila Ghaffari   for (int i = 0; i < N; i++) blasius->Tf_cheb[i] = cheb_coefs[i];
104e0d1a4dfSLeila Ghaffari   for (int i = 0; i < N - 1; i++) blasius->Th_cheb[i] = cheb_coefs[i + N];
1050d850f2eSLeila Ghaffari 
1060d850f2eSLeila Ghaffari   // Destroy objects
107e0d1a4dfSLeila Ghaffari   PetscCall(PetscFree2(blasius->X, w));
108e0d1a4dfSLeila Ghaffari   PetscCall(VecDestroy(&sol));
109e0d1a4dfSLeila Ghaffari   PetscCall(VecDestroy(&res));
110e0d1a4dfSLeila Ghaffari   PetscCall(SNESDestroy(&snes));
111d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
112e0d1a4dfSLeila Ghaffari }
113e0d1a4dfSLeila Ghaffari 
BlasiusInflowBCSetup_CreateIFunctionQF(BCDefinition bc_def,CeedQFunction * qf)114d6cac220SJames Wright static PetscErrorCode BlasiusInflowBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) {
115d6cac220SJames Wright   HoneeBCStruct honee_bc;
116d6cac220SJames Wright 
117d6cac220SJames Wright   PetscFunctionBeginUser;
118d6cac220SJames Wright   PetscCall(BCDefinitionGetContext(bc_def, &honee_bc));
119d6cac220SJames Wright   PetscCall(HoneeBCCreateIFunctionQF(bc_def, Blasius_Inflow, Blasius_Inflow_loc, honee_bc->qfctx, qf));
120d6cac220SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
121d6cac220SJames Wright }
122d6cac220SJames Wright 
BlasiusInflowBCSetup_CreateIJacobianQF(BCDefinition bc_def,CeedQFunction * qf)123d6cac220SJames Wright static PetscErrorCode BlasiusInflowBCSetup_CreateIJacobianQF(BCDefinition bc_def, CeedQFunction *qf) {
124d6cac220SJames Wright   HoneeBCStruct honee_bc;
125d6cac220SJames Wright 
126d6cac220SJames Wright   PetscFunctionBeginUser;
127d6cac220SJames Wright   PetscCall(BCDefinitionGetContext(bc_def, &honee_bc));
128d6cac220SJames Wright   PetscCall(HoneeBCCreateIJacobianQF(bc_def, Blasius_Inflow_Jacobian, Blasius_Inflow_Jacobian_loc, honee_bc->qfctx, qf));
129d6cac220SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
130d6cac220SJames Wright }
131d6cac220SJames Wright 
NS_BLASIUS(ProblemData problem,DM dm,void * ctx)132d3c60affSJames Wright PetscErrorCode NS_BLASIUS(ProblemData problem, DM dm, void *ctx) {
1330c373b74SJames Wright   Honee                    honee   = *(Honee *)ctx;
1340c373b74SJames Wright   MPI_Comm                 comm    = honee->comm;
1350c373b74SJames Wright   Ceed                     ceed    = honee->ceed;
136493642f1SJames Wright   PetscBool                use_stg = PETSC_FALSE;
13715a3537eSJed Brown   BlasiusContext           blasius_ctx;
13815a3537eSJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
139e07531f7SJames Wright   CeedQFunctionContext     blasius_qfctx;
14015a3537eSJed Brown 
141bb8a0c61SJames Wright   PetscFunctionBeginUser;
142d3c60affSJames Wright   PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx));
143bb8a0c61SJames Wright 
144bb8a0c61SJames Wright   // ------------------------------------------------------
145bb8a0c61SJames Wright   //               SET UP Blasius
146bb8a0c61SJames Wright   // ------------------------------------------------------
147f5dc303cSJames Wright   problem->ics = (HoneeQFSpec){.qf_func_ptr = ICsBlasius, .qf_loc = ICsBlasius_loc, .qfctx = problem->ics.qfctx};
148bb8a0c61SJames Wright 
149aef1eb53SLeila Ghaffari   CeedScalar U_inf  = 40;           // m/s
150aef1eb53SLeila Ghaffari   CeedScalar T_inf  = 288.;         // K
1510d850f2eSLeila Ghaffari   CeedScalar T_wall = 288.;         // K
152e0d1a4dfSLeila Ghaffari   CeedScalar delta0 = 4.2e-3;       // m
153fcb2c22aSJames Wright   CeedScalar P_inf  = 1.01e5;       // Pa
154defe8520SJames Wright   PetscInt   N      = 20;           // Number of Chebyshev terms
1552acc7cbcSKenneth E. Jansen   PetscBool  weakT  = PETSC_FALSE;  // weak density or temperature
156a213b8aaSJames Wright   PetscBool  P0_set;
157bb8a0c61SJames Wright 
1583fd71269SJames Wright   PetscOptionsBegin(comm, NULL, "Options for BLASIUS problem", NULL);
1592b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", NULL, weakT, &weakT, NULL));
1602b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-velocity_infinity", "Velocity at boundary layer edge", NULL, U_inf, &U_inf, NULL));
1612b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-temperature_infinity", "Temperature at boundary layer edge", NULL, T_inf, &T_inf, NULL));
162a213b8aaSJames Wright   PetscCall(PetscOptionsHasName(NULL, NULL, "-P0", &P0_set));  // For maintaining behavior of -P0 flag (which is deprecated)
1639eadbee4SJames Wright   PetscCall(PetscOptionsDeprecated("-P0", "-pressure_infinity", "libCEED 0.12.0",
1649eadbee4SJames Wright                                    "Use -pressure_infinity to set pressure at boundary layer edge and -idl_pressure to set the IDL reference "
1659eadbee4SJames Wright                                    "pressure"));
166a213b8aaSJames Wright   PetscCall(PetscOptionsScalar("-pressure_infinity", "Pressure at boundary layer edge", NULL, P_inf, &P_inf, NULL));
1672b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-temperature_wall", "Temperature at wall", NULL, T_wall, &T_wall, NULL));
1682b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-delta0", "Boundary layer height at inflow", NULL, delta0, &delta0, NULL));
1692b916ea7SJeremy L Thompson   PetscCall(PetscOptionsInt("-n_chebyshev", "Number of Chebyshev terms", NULL, N, &N, NULL));
1702b916ea7SJeremy L Thompson   PetscCheck(3 <= N && N <= BLASIUS_MAX_N_CHEBYSHEV, comm, PETSC_ERR_ARG_OUTOFRANGE, "-n_chebyshev %" PetscInt_FMT " must be in range [3, %d]", N,
1710d850f2eSLeila Ghaffari              BLASIUS_MAX_N_CHEBYSHEV);
1722b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", NULL, use_stg, &use_stg, NULL));
173bb8a0c61SJames Wright   PetscOptionsEnd();
174bb8a0c61SJames Wright 
175c9f37605SMohammed Amin   Units units = honee->units;
176bb8a0c61SJames Wright 
177c9f37605SMohammed Amin   T_inf *= units->Kelvin;
178c9f37605SMohammed Amin   T_wall *= units->Kelvin;
179c9f37605SMohammed Amin   P_inf *= units->Pascal;
180c9f37605SMohammed Amin   U_inf *= units->meter / units->second;
181c9f37605SMohammed Amin   delta0 *= units->meter;
182bb8a0c61SJames Wright 
18315a3537eSJed Brown   // Some properties depend on parameters from NewtonianIdealGas
184e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ig_ctx));
185bb8a0c61SJames Wright 
186fcb2c22aSJames Wright   StatePrimitive Y_inf = {
187fcb2c22aSJames Wright       .pressure = P_inf, .velocity = {U_inf, 0, 0},
188fcb2c22aSJames Wright            .temperature = T_inf
189fcb2c22aSJames Wright   };
190cde3d787SJames Wright   State S_infty = StateFromPrimitive(newtonian_ig_ctx->gas, Y_inf);
191fcb2c22aSJames Wright 
1927e3656bdSJames Wright   PetscCall(PetscNew(&blasius_ctx));
193493642f1SJames Wright   blasius_ctx->weakT    = weakT;
194e0d1a4dfSLeila Ghaffari   blasius_ctx->T_wall   = T_wall;
19515a3537eSJed Brown   blasius_ctx->delta0   = delta0;
196fcb2c22aSJames Wright   blasius_ctx->S_infty  = S_infty;
197e0d1a4dfSLeila Ghaffari   blasius_ctx->n_cheb   = N;
1980c373b74SJames Wright   blasius_ctx->implicit = honee->phys->implicit;
199a213b8aaSJames Wright   if (P0_set) newtonian_ig_ctx->idl_pressure = P_inf;  // For maintaining behavior of -P0 flag (which is deprecated)
200cde3d787SJames Wright   blasius_ctx->newt_ctx = *newtonian_ig_ctx;
201493642f1SJames Wright 
202ef2c71fdSJames Wright   {
203e0d1a4dfSLeila Ghaffari     PetscReal domain_min[3], domain_max[3];
2042b916ea7SJeremy L Thompson     PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
205ef2c71fdSJames Wright     blasius_ctx->x_inflow = domain_min[0];
206e0d1a4dfSLeila Ghaffari     blasius_ctx->eta_max  = 5 * domain_max[1] / blasius_ctx->delta0;
207ef2c71fdSJames Wright   }
20825c92e8fSJames Wright   PetscBool diff_filter_mms = PETSC_FALSE;
20925c92e8fSJames Wright   PetscCall(PetscOptionsGetBool(NULL, NULL, "-diff_filter_mms", &diff_filter_mms, NULL));
21025c92e8fSJames Wright   if (!use_stg && !diff_filter_mms) PetscCall(ComputeChebyshevCoefficients(blasius_ctx));
211ef2c71fdSJames Wright 
212e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ig_ctx));
213bb8a0c61SJames Wright 
2140c373b74SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &blasius_qfctx));
215e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(blasius_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*blasius_ctx), blasius_ctx));
216e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(blasius_qfctx, CEED_MEM_HOST, FreeContextPetsc));
217bb8a0c61SJames Wright 
218e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfctx));
219e07531f7SJames Wright   problem->ics.qfctx = blasius_qfctx;
220493642f1SJames Wright   if (use_stg) {
2210c373b74SJames Wright     PetscCall(SetupStg(comm, dm, problem, honee, weakT, S_infty.Y.temperature, S_infty.Y.pressure));
22225c92e8fSJames Wright   } else if (diff_filter_mms) {
2238d78d7c8SJames Wright     PetscCall(DifferentialFilterMmsICSetup(honee));
2248085925cSJames Wright   } else {
2250c373b74SJames Wright     PetscCheck((honee->phys->state_var == STATEVAR_CONSERVATIVE) || (honee->app_ctx->test_type == TESTTYPE_DIFF_FILTER), honee->comm,
226727ec889SJames Wright                PETSC_ERR_ARG_INCOMP, "Can only use conservative variables with Blasius and weak inflow");
227d6cac220SJames Wright     for (PetscCount b = 0; b < problem->num_bc_defs; b++) {
228d6cac220SJames Wright       BCDefinition bc_def = problem->bc_defs[b];
229d6cac220SJames Wright       const char  *name;
230d6cac220SJames Wright 
231d6cac220SJames Wright       PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL));
232d6cac220SJames Wright       if (!strcmp(name, "inflow")) {
233d6cac220SJames Wright         HoneeBCStruct honee_bc;
234d6cac220SJames Wright 
235d6cac220SJames Wright         PetscCall(PetscNew(&honee_bc));
236d6cac220SJames Wright         PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(blasius_qfctx, &honee_bc->qfctx));
237d6cac220SJames Wright         honee_bc->honee              = honee;
2381abc2837SJames Wright         honee_bc->num_comps_jac_data = 0;
239*26d401f3SJames Wright         PetscCall(BCDefinitionSetContext(bc_def, (PetscCtxDestroyFn *)HoneeBCDestroy, honee_bc));
240d6cac220SJames Wright 
241d6cac220SJames Wright         PetscCall(BCDefinitionSetIFunction(bc_def, BlasiusInflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp));
242d6cac220SJames Wright         PetscCall(BCDefinitionSetIJacobian(bc_def, BlasiusInflowBCSetup_CreateIJacobianQF, HoneeBCAddIJacobianOp));
243d6cac220SJames Wright       }
244d6cac220SJames Wright     }
245493642f1SJames Wright   }
246d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
247bb8a0c61SJames Wright }
248