xref: /honee/problems/channel.c (revision f978755d57fb6574cfe1ff974b5424124ae3c75e)
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 
16*f978755dSJames Wright static PetscErrorCode ChannelOutflowBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) {
17*f978755dSJames Wright   HoneeBCStruct honee_bc;
18*f978755dSJames Wright 
19*f978755dSJames Wright   PetscFunctionBeginUser;
20*f978755dSJames Wright   PetscCall(BCDefinitionGetContext(bc_def, &honee_bc));
21*f978755dSJames Wright   PetscCheck(honee_bc->honee->phys->state_var == STATEVAR_CONSERVATIVE, PETSC_COMM_WORLD, PETSC_ERR_SUP,
22*f978755dSJames Wright              "Channel outflow only valid for Conservative variables, recieved %s", StateVariables[honee_bc->honee->phys->state_var]);
23*f978755dSJames Wright   PetscCall(HoneeBCCreateIFunctionQF(bc_def, Channel_Outflow, Channel_Outflow_loc, honee_bc->qfctx, qf));
24*f978755dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
25*f978755dSJames Wright }
26*f978755dSJames Wright 
27991aef52SJames Wright PetscErrorCode NS_CHANNEL(ProblemData problem, DM dm, void *ctx, SimpleBC bc) {
280c373b74SJames Wright   Honee                    honee = *(Honee *)ctx;
290c373b74SJames Wright   MPI_Comm                 comm  = honee->comm;
300c373b74SJames Wright   Ceed                     ceed  = honee->ceed;
3115a3537eSJed Brown   ChannelContext           channel_ctx;
3215a3537eSJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
33e07531f7SJames Wright   CeedQFunctionContext     channel_qfctx;
34b3b24828SJames Wright   PetscBool                use_divdiff_verify_mesh = PETSC_FALSE;
3515a3537eSJed Brown 
36bb8a0c61SJames Wright   PetscFunctionBeginUser;
37d1c51a42SJames Wright   PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc));
382b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &channel_ctx));
39bb8a0c61SJames Wright 
40b3b24828SJames Wright   PetscCall(PetscOptionsGetBool(NULL, NULL, "-mesh_transform_channel_div_diff_projection_verify", &use_divdiff_verify_mesh, NULL));
41b3b24828SJames Wright   if (use_divdiff_verify_mesh) PetscCall(DivDiffFluxVerifyMesh(dm));
42b3b24828SJames Wright 
43bb8a0c61SJames Wright   // ------------------------------------------------------
44bb8a0c61SJames Wright   //               SET UP Channel
45bb8a0c61SJames Wright   // ------------------------------------------------------
46e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfctx));
47e07531f7SJames Wright   problem->ics.qf_func_ptr = ICsChannel;
48e07531f7SJames Wright   problem->ics.qf_loc      = ICsChannel_loc;
490c373b74SJames Wright   if (honee->phys->state_var == STATEVAR_CONSERVATIVE) {
50e07531f7SJames Wright     problem->apply_inflow.qf_func_ptr = Channel_Inflow;
51e07531f7SJames Wright     problem->apply_inflow.qf_loc      = Channel_Inflow_loc;
52cbe60e31SLeila Ghaffari   }
53bb8a0c61SJames 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 
660c373b74SJames Wright   PetscScalar meter  = honee->units->meter;
670c373b74SJames Wright   PetscScalar second = honee->units->second;
680c373b74SJames Wright   PetscScalar Kelvin = honee->units->Kelvin;
690c373b74SJames Wright   PetscScalar Pascal = honee->units->Pascal;
70bb8a0c61SJames Wright 
71bb8a0c61SJames Wright   theta0 *= Kelvin;
72bb8a0c61SJames Wright   P0 *= Pascal;
73bb8a0c61SJames Wright   umax *= meter / second;
74bb8a0c61SJames Wright 
75bb8a0c61SJames Wright   //-- Setup Problem information
76bb8a0c61SJames Wright   CeedScalar H, center;
77bb8a0c61SJames Wright   {
78bb8a0c61SJames Wright     PetscReal domain_min[3], domain_max[3], domain_size[3];
792b916ea7SJeremy L Thompson     PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
80493642f1SJames Wright     for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
81bb8a0c61SJames Wright 
82bb8a0c61SJames Wright     H      = 0.5 * domain_size[1] * meter;
83bb8a0c61SJames Wright     center = H + domain_min[1] * meter;
84bb8a0c61SJames Wright   }
85bb8a0c61SJames Wright 
8615a3537eSJed Brown   // Some properties depend on parameters from NewtonianIdealGas
87e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ig_ctx));
8815a3537eSJed Brown 
8915a3537eSJed Brown   channel_ctx->center   = center;
9015a3537eSJed Brown   channel_ctx->H        = H;
9115a3537eSJed Brown   channel_ctx->theta0   = theta0;
9215a3537eSJed Brown   channel_ctx->P0       = P0;
9315a3537eSJed Brown   channel_ctx->umax     = umax;
940c373b74SJames Wright   channel_ctx->implicit = honee->phys->implicit;
9510786903SJed Brown   channel_ctx->B        = body_force_scale * 2 * umax * newtonian_ig_ctx->mu / (H * H);
96bb8a0c61SJames Wright 
97bb8a0c61SJames Wright   {
98bb8a0c61SJames Wright     // Calculate Body force
992b916ea7SJeremy L Thompson     CeedScalar cv = newtonian_ig_ctx->cv, cp = newtonian_ig_ctx->cp;
100bb8a0c61SJames Wright     CeedScalar Rd  = cp - cv;
101bb8a0c61SJames Wright     CeedScalar rho = P0 / (Rd * theta0);
10215a3537eSJed Brown     CeedScalar g[] = {channel_ctx->B / rho, 0., 0.};
1032b916ea7SJeremy L Thompson     PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3));
104bb8a0c61SJames Wright   }
10515a3537eSJed Brown   channel_ctx->newtonian_ctx = *newtonian_ig_ctx;
106e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ig_ctx));
107bb8a0c61SJames Wright 
1080c373b74SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &channel_qfctx));
109e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(channel_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*channel_ctx), channel_ctx));
110e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(channel_qfctx, CEED_MEM_HOST, FreeContextPetsc));
11115a3537eSJed Brown 
112e07531f7SJames Wright   problem->ics.qfctx = channel_qfctx;
113e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(channel_qfctx, &problem->apply_inflow.qfctx));
114*f978755dSJames Wright 
115*f978755dSJames Wright   for (PetscCount b = 0; b < problem->num_bc_defs; b++) {
116*f978755dSJames Wright     BCDefinition bc_def = problem->bc_defs[b];
117*f978755dSJames Wright     const char  *name;
118*f978755dSJames Wright 
119*f978755dSJames Wright     PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL));
120*f978755dSJames Wright     if (honee->phys->state_var == STATEVAR_CONSERVATIVE && !strcmp(name, "outflow")) {
121*f978755dSJames Wright       HoneeBCStruct honee_bc;
122*f978755dSJames Wright 
123*f978755dSJames Wright       PetscCall(PetscPrintf(comm, "WARNING! Channel flow with Inflow and Outflow is currently broken.\n"));
124*f978755dSJames Wright       PetscCall(PetscNew(&honee_bc));
125*f978755dSJames Wright       PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(channel_qfctx, &honee_bc->qfctx));
126*f978755dSJames Wright       honee_bc->honee             = honee;
127*f978755dSJames Wright       honee_bc->jac_data_size_sur = honee->phys->implicit ? problem->jac_data_size_sur : 0;
128*f978755dSJames Wright       PetscCall(BCDefinitionSetContext(bc_def, HoneeBCDestroy, honee_bc));
129*f978755dSJames Wright 
130*f978755dSJames Wright       PetscCall(BCDefinitionSetIFunction(bc_def, ChannelOutflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp));
131*f978755dSJames Wright       PetscCall(BCDefinitionSetIJacobian(bc_def, NULL, NULL));
132*f978755dSJames Wright     }
133*f978755dSJames Wright   }
134d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
135bb8a0c61SJames Wright }
136b3b24828SJames Wright 
137b3b24828SJames Wright // This function transforms the mesh coordinates to mimic the mesh used in
138b3b24828SJames Wright // *A better consistency for low-order stabilized finite element methods* Jansen et. al. 1999
139b3b24828SJames Wright // which is used to verify the projection of divergence of diffusive flux. See !27 and !94 for more details.
140b3b24828SJames Wright static PetscErrorCode DivDiffFluxVerifyMesh(DM dm) {
141b3b24828SJames Wright   PetscInt     narr, ncoords, dim;
142b3b24828SJames Wright   PetscReal    domain_min[3], domain_max[3], domain_size[3];
143b3b24828SJames Wright   PetscScalar *arr_coords;
144b3b24828SJames Wright   Vec          vec_coords;
145b3b24828SJames Wright 
146b3b24828SJames Wright   PetscFunctionBeginUser;
147b3b24828SJames Wright   PetscCall(DMGetDimension(dm, &dim));
148b3b24828SJames Wright   // Get domain boundary information
149b3b24828SJames Wright   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
150b3b24828SJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
151b3b24828SJames Wright 
152b3b24828SJames Wright   // Get coords array from DM
153b3b24828SJames Wright   PetscCall(DMGetCoordinatesLocal(dm, &vec_coords));
154b3b24828SJames Wright   PetscCall(VecGetLocalSize(vec_coords, &narr));
155b3b24828SJames Wright   PetscCall(VecGetArray(vec_coords, &arr_coords));
156b3b24828SJames Wright 
157b3b24828SJames Wright   PetscScalar(*coords)[dim] = (PetscScalar(*)[dim])arr_coords;
158b3b24828SJames Wright   ncoords                   = narr / dim;
159b3b24828SJames Wright 
160b3b24828SJames Wright   // Get mesh information
161b3b24828SJames Wright   PetscInt nmax = 3, faces[3];
162b3b24828SJames Wright   PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL));
163b3b24828SJames Wright   // Get element size of the box mesh, for indexing each node
164b3b24828SJames Wright   const PetscReal dxbox = domain_size[0] / (faces[0]);
165b3b24828SJames Wright 
166b3b24828SJames Wright   for (PetscInt i = 0; i < ncoords; i++) {
167b3b24828SJames Wright     PetscInt x_box_index = round(coords[i][0] / dxbox);
168b3b24828SJames Wright     if (x_box_index % 2) {
169b3b24828SJames Wright       coords[i][0] = (x_box_index - 1) * dxbox + 0.5 * dxbox;
170b3b24828SJames Wright     }
171b3b24828SJames Wright   }
172b3b24828SJames Wright 
173b3b24828SJames Wright   PetscCall(VecRestoreArray(vec_coords, &arr_coords));
174b3b24828SJames Wright   PetscCall(DMSetCoordinatesLocal(dm, vec_coords));
175b3b24828SJames Wright   PetscFunctionReturn(0);
176b3b24828SJames Wright }
177