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 Channel flow
6bb8a0c61SJames Wright
7bb8a0c61SJames Wright #include "../qfunctions/channel.h"
8bb8a0c61SJames Wright
9e419654dSJeremy L Thompson #include <ceed.h>
10e419654dSJeremy L Thompson #include <petscdm.h>
11e419654dSJeremy L Thompson
12149fb536SJames Wright #include <navierstokes.h>
13bb8a0c61SJames Wright
14b3b24828SJames Wright static PetscErrorCode DivDiffFluxVerifyMesh(DM dm);
15b3b24828SJames Wright
ChannelOutflowBCSetup_CreateIFunctionQF(BCDefinition bc_def,CeedQFunction * qf)16f978755dSJames Wright static PetscErrorCode ChannelOutflowBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) {
17f978755dSJames Wright HoneeBCStruct honee_bc;
18f978755dSJames Wright
19f978755dSJames Wright PetscFunctionBeginUser;
20f978755dSJames Wright PetscCall(BCDefinitionGetContext(bc_def, &honee_bc));
21f978755dSJames Wright PetscCheck(honee_bc->honee->phys->state_var == STATEVAR_CONSERVATIVE, PETSC_COMM_WORLD, PETSC_ERR_SUP,
22f978755dSJames Wright "Channel outflow only valid for Conservative variables, recieved %s", StateVariables[honee_bc->honee->phys->state_var]);
23f978755dSJames Wright PetscCall(HoneeBCCreateIFunctionQF(bc_def, Channel_Outflow, Channel_Outflow_loc, honee_bc->qfctx, qf));
24f978755dSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
25f978755dSJames Wright }
26f978755dSJames Wright
ChannelInflowBCSetup_CreateIFunctionQF(BCDefinition bc_def,CeedQFunction * qf)27d6cac220SJames Wright static PetscErrorCode ChannelInflowBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) {
28d6cac220SJames Wright HoneeBCStruct honee_bc;
29d6cac220SJames Wright
30d6cac220SJames Wright PetscFunctionBeginUser;
31d6cac220SJames Wright PetscCall(BCDefinitionGetContext(bc_def, &honee_bc));
32d6cac220SJames Wright PetscCheck(honee_bc->honee->phys->state_var == STATEVAR_CONSERVATIVE, PETSC_COMM_WORLD, PETSC_ERR_SUP,
33d6cac220SJames Wright "Channel inflow only valid for Conservative variables, recieved %s", StateVariables[honee_bc->honee->phys->state_var]);
34d6cac220SJames Wright PetscCall(HoneeBCCreateIFunctionQF(bc_def, Channel_Inflow, Channel_Inflow_loc, honee_bc->qfctx, qf));
35d6cac220SJames Wright PetscFunctionReturn(PETSC_SUCCESS);
36d6cac220SJames Wright }
37d6cac220SJames Wright
NS_CHANNEL(ProblemData problem,DM dm,void * ctx)38d3c60affSJames Wright PetscErrorCode NS_CHANNEL(ProblemData problem, DM dm, void *ctx) {
390c373b74SJames Wright Honee honee = *(Honee *)ctx;
400c373b74SJames Wright MPI_Comm comm = honee->comm;
410c373b74SJames Wright Ceed ceed = honee->ceed;
4215a3537eSJed Brown ChannelContext channel_ctx;
4315a3537eSJed Brown NewtonianIdealGasContext newtonian_ig_ctx;
44e07531f7SJames Wright CeedQFunctionContext channel_qfctx;
45b3b24828SJames Wright PetscBool use_divdiff_verify_mesh = PETSC_FALSE;
4615a3537eSJed Brown
47bb8a0c61SJames Wright PetscFunctionBeginUser;
48d3c60affSJames Wright PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx));
492d898fa6SJames Wright PetscCall(PetscNew(&channel_ctx));
50bb8a0c61SJames Wright
51b3b24828SJames Wright PetscCall(PetscOptionsGetBool(NULL, NULL, "-mesh_transform_channel_div_diff_projection_verify", &use_divdiff_verify_mesh, NULL));
52b3b24828SJames Wright if (use_divdiff_verify_mesh) PetscCall(DivDiffFluxVerifyMesh(dm));
53b3b24828SJames Wright
54bb8a0c61SJames Wright // -- Command Line Options
55bb8a0c61SJames Wright CeedScalar umax = 10.; // m/s
56bb8a0c61SJames Wright CeedScalar theta0 = 300.; // K
57bb8a0c61SJames Wright CeedScalar P0 = 1.e5; // Pa
5810786903SJed Brown PetscReal body_force_scale = 1.;
59bb8a0c61SJames Wright PetscOptionsBegin(comm, NULL, "Options for CHANNEL problem", NULL);
602b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-umax", "Centerline velocity of the Channel", NULL, umax, &umax, NULL));
612b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-theta0", "Wall temperature", NULL, theta0, &theta0, NULL));
622b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-P0", "Pressure at outflow", NULL, P0, &P0, NULL));
632b916ea7SJeremy L Thompson PetscCall(PetscOptionsReal("-body_force_scale", "Multiplier for body force", NULL, body_force_scale = 1, &body_force_scale, NULL));
64bb8a0c61SJames Wright PetscOptionsEnd();
65bb8a0c61SJames Wright
66c9f37605SMohammed Amin Units units = honee->units;
67bb8a0c61SJames Wright
68c9f37605SMohammed Amin theta0 *= units->Kelvin;
69c9f37605SMohammed Amin P0 *= units->Pascal;
70c9f37605SMohammed Amin umax *= units->meter / units->second;
71bb8a0c61SJames Wright
72bb8a0c61SJames Wright //-- Setup Problem information
73bb8a0c61SJames Wright CeedScalar H, center;
74bb8a0c61SJames Wright {
75bb8a0c61SJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3];
762b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
77493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
78bb8a0c61SJames Wright
797e3656bdSJames Wright H = 0.5 * domain_size[1];
807e3656bdSJames Wright center = H + domain_min[1];
81bb8a0c61SJames Wright }
82bb8a0c61SJames Wright
8315a3537eSJed Brown // Some properties depend on parameters from NewtonianIdealGas
84e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ig_ctx));
8515a3537eSJed Brown
8615a3537eSJed Brown channel_ctx->center = center;
8715a3537eSJed Brown channel_ctx->H = H;
8815a3537eSJed Brown channel_ctx->theta0 = theta0;
8915a3537eSJed Brown channel_ctx->P0 = P0;
9015a3537eSJed Brown channel_ctx->umax = umax;
910c373b74SJames Wright channel_ctx->implicit = honee->phys->implicit;
92cde3d787SJames Wright channel_ctx->B = body_force_scale * 2 * umax * newtonian_ig_ctx->gas.mu / (H * H);
93bb8a0c61SJames Wright
94bb8a0c61SJames Wright {
95bb8a0c61SJames Wright // Calculate Body force
96cde3d787SJames Wright CeedScalar cv = newtonian_ig_ctx->gas.cv, cp = newtonian_ig_ctx->gas.cp;
97bb8a0c61SJames Wright CeedScalar Rd = cp - cv;
98bb8a0c61SJames Wright CeedScalar rho = P0 / (Rd * theta0);
9915a3537eSJed Brown CeedScalar g[] = {channel_ctx->B / rho, 0., 0.};
1002b916ea7SJeremy L Thompson PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3));
101bb8a0c61SJames Wright }
102cde3d787SJames Wright channel_ctx->newt_ctx = *newtonian_ig_ctx;
103e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ig_ctx));
104bb8a0c61SJames Wright
1050c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &channel_qfctx));
106e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(channel_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*channel_ctx), channel_ctx));
107e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(channel_qfctx, CEED_MEM_HOST, FreeContextPetsc));
108f5dc303cSJames Wright
109f5dc303cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfctx));
110f5dc303cSJames Wright problem->ics = (HoneeQFSpec){.qf_func_ptr = ICsChannel, .qf_loc = ICsChannel_loc, .qfctx = channel_qfctx};
111f978755dSJames Wright
112f978755dSJames Wright for (PetscCount b = 0; b < problem->num_bc_defs; b++) {
113f978755dSJames Wright BCDefinition bc_def = problem->bc_defs[b];
114f978755dSJames Wright const char *name;
115f978755dSJames Wright
116f978755dSJames Wright PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL));
117f978755dSJames Wright if (honee->phys->state_var == STATEVAR_CONSERVATIVE && !strcmp(name, "outflow")) {
118f978755dSJames Wright HoneeBCStruct honee_bc;
119f978755dSJames Wright
120f978755dSJames Wright PetscCall(PetscPrintf(comm, "WARNING! Channel flow with Inflow and Outflow is currently broken.\n"));
121f978755dSJames Wright PetscCall(PetscNew(&honee_bc));
122f978755dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(channel_qfctx, &honee_bc->qfctx));
123f978755dSJames Wright honee_bc->honee = honee;
1241abc2837SJames Wright honee_bc->num_comps_jac_data = honee->phys->implicit ? 11 : 0;
125*26d401f3SJames Wright PetscCall(BCDefinitionSetContext(bc_def, (PetscCtxDestroyFn *)HoneeBCDestroy, honee_bc));
126f978755dSJames Wright
127f978755dSJames Wright PetscCall(BCDefinitionSetIFunction(bc_def, ChannelOutflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp));
128f978755dSJames Wright PetscCall(BCDefinitionSetIJacobian(bc_def, NULL, NULL));
129d6cac220SJames Wright } else if (honee->phys->state_var == STATEVAR_CONSERVATIVE && !strcmp(name, "inflow")) {
130d6cac220SJames Wright HoneeBCStruct honee_bc;
131d6cac220SJames Wright
132d6cac220SJames Wright PetscCall(PetscPrintf(comm, "WARNING! Channel flow with Inflow and Outflow is currently broken.\n"));
133d6cac220SJames Wright PetscCall(PetscNew(&honee_bc));
134d6cac220SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(channel_qfctx, &honee_bc->qfctx));
135d6cac220SJames Wright honee_bc->honee = honee;
1361abc2837SJames Wright honee_bc->num_comps_jac_data = honee->phys->implicit ? 11 : 0;
137*26d401f3SJames Wright PetscCall(BCDefinitionSetContext(bc_def, (PetscCtxDestroyFn *)HoneeBCDestroy, honee_bc));
138d6cac220SJames Wright
139d6cac220SJames Wright PetscCall(BCDefinitionSetIFunction(bc_def, ChannelInflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp));
140d6cac220SJames Wright PetscCall(BCDefinitionSetIJacobian(bc_def, NULL, NULL));
141f978755dSJames Wright }
142f978755dSJames Wright }
143d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
144bb8a0c61SJames Wright }
145b3b24828SJames Wright
146b3b24828SJames Wright // This function transforms the mesh coordinates to mimic the mesh used in
147b3b24828SJames Wright // *A better consistency for low-order stabilized finite element methods* Jansen et. al. 1999
148b3b24828SJames Wright // which is used to verify the projection of divergence of diffusive flux. See !27 and !94 for more details.
DivDiffFluxVerifyMesh(DM dm)149b3b24828SJames Wright static PetscErrorCode DivDiffFluxVerifyMesh(DM dm) {
150b3b24828SJames Wright PetscInt narr, ncoords, dim;
151b3b24828SJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3];
152b3b24828SJames Wright PetscScalar *arr_coords;
153b3b24828SJames Wright Vec vec_coords;
154b3b24828SJames Wright
155b3b24828SJames Wright PetscFunctionBeginUser;
156b3b24828SJames Wright PetscCall(DMGetDimension(dm, &dim));
157b3b24828SJames Wright // Get domain boundary information
158b3b24828SJames Wright PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
159b3b24828SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
160b3b24828SJames Wright
161b3b24828SJames Wright // Get coords array from DM
162b3b24828SJames Wright PetscCall(DMGetCoordinatesLocal(dm, &vec_coords));
163b3b24828SJames Wright PetscCall(VecGetLocalSize(vec_coords, &narr));
164b3b24828SJames Wright PetscCall(VecGetArray(vec_coords, &arr_coords));
165b3b24828SJames Wright
166b3b24828SJames Wright PetscScalar(*coords)[dim] = (PetscScalar(*)[dim])arr_coords;
167b3b24828SJames Wright ncoords = narr / dim;
168b3b24828SJames Wright
169b3b24828SJames Wright // Get mesh information
170b3b24828SJames Wright PetscInt nmax = 3, faces[3];
171b3b24828SJames Wright PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL));
172b3b24828SJames Wright // Get element size of the box mesh, for indexing each node
173b3b24828SJames Wright const PetscReal dxbox = domain_size[0] / (faces[0]);
174b3b24828SJames Wright
175b3b24828SJames Wright for (PetscInt i = 0; i < ncoords; i++) {
176b3b24828SJames Wright PetscInt x_box_index = round(coords[i][0] / dxbox);
177b3b24828SJames Wright if (x_box_index % 2) {
178b3b24828SJames Wright coords[i][0] = (x_box_index - 1) * dxbox + 0.5 * dxbox;
179b3b24828SJames Wright }
180b3b24828SJames Wright }
181b3b24828SJames Wright
182b3b24828SJames Wright PetscCall(VecRestoreArray(vec_coords, &arr_coords));
183b3b24828SJames Wright PetscCall(DMSetCoordinatesLocal(dm, vec_coords));
184b3b24828SJames Wright PetscFunctionReturn(0);
185b3b24828SJames Wright }
186