xref: /honee/problems/channel.c (revision d949ddfc92e3b3e8aad90700e3fb60ca7f4585db)
1bb8a0c61SJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2bb8a0c61SJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3bb8a0c61SJames Wright //
4bb8a0c61SJames Wright // SPDX-License-Identifier: BSD-2-Clause
5bb8a0c61SJames Wright //
6bb8a0c61SJames Wright // This file is part of CEED:  http://github.com/ceed
7bb8a0c61SJames Wright 
8bb8a0c61SJames Wright /// @file
9bb8a0c61SJames Wright /// Utility functions for setting up Channel flow
10bb8a0c61SJames Wright 
11bb8a0c61SJames Wright #include "../qfunctions/channel.h"
12bb8a0c61SJames Wright 
13e419654dSJeremy L Thompson #include <ceed.h>
14e419654dSJeremy L Thompson #include <petscdm.h>
15e419654dSJeremy L Thompson 
162b916ea7SJeremy L Thompson #include "../navierstokes.h"
17bb8a0c61SJames Wright 
18d1c51a42SJames Wright PetscErrorCode NS_CHANNEL(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) {
19bb8a0c61SJames Wright   User                     user = *(User *)ctx;
20bb8a0c61SJames Wright   MPI_Comm                 comm = PETSC_COMM_WORLD;
2115a3537eSJed Brown   ChannelContext           channel_ctx;
2215a3537eSJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
2315a3537eSJed Brown   CeedQFunctionContext     channel_context;
2415a3537eSJed Brown 
25bb8a0c61SJames Wright   PetscFunctionBeginUser;
26d1c51a42SJames Wright   PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc));
272b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &channel_ctx));
28bb8a0c61SJames Wright 
29bb8a0c61SJames Wright   // ------------------------------------------------------
30bb8a0c61SJames Wright   //               SET UP Channel
31bb8a0c61SJames Wright   // ------------------------------------------------------
3215a3537eSJed Brown   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
339785fe93SJed Brown   problem->ics.qfunction     = ICsChannel;
349785fe93SJed Brown   problem->ics.qfunction_loc = ICsChannel_loc;
353636f6a4SJames Wright   if (user->phys->state_var == STATEVAR_CONSERVATIVE) {
369785fe93SJed Brown     problem->apply_inflow.qfunction      = Channel_Inflow;
379785fe93SJed Brown     problem->apply_inflow.qfunction_loc  = Channel_Inflow_loc;
389785fe93SJed Brown     problem->apply_outflow.qfunction     = Channel_Outflow;
399785fe93SJed Brown     problem->apply_outflow.qfunction_loc = Channel_Outflow_loc;
40cbe60e31SLeila Ghaffari   }
41bb8a0c61SJames Wright 
42bb8a0c61SJames Wright   // -- Command Line Options
43bb8a0c61SJames Wright   CeedScalar umax             = 10.;   // m/s
44bb8a0c61SJames Wright   CeedScalar theta0           = 300.;  // K
45bb8a0c61SJames Wright   CeedScalar P0               = 1.e5;  // Pa
4610786903SJed Brown   PetscReal  body_force_scale = 1.;
47bb8a0c61SJames Wright   PetscOptionsBegin(comm, NULL, "Options for CHANNEL problem", NULL);
482b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-umax", "Centerline velocity of the Channel", NULL, umax, &umax, NULL));
492b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-theta0", "Wall temperature", NULL, theta0, &theta0, NULL));
502b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-P0", "Pressure at outflow", NULL, P0, &P0, NULL));
512b916ea7SJeremy L Thompson   PetscCall(PetscOptionsReal("-body_force_scale", "Multiplier for body force", NULL, body_force_scale = 1, &body_force_scale, NULL));
52bb8a0c61SJames Wright   PetscOptionsEnd();
53bb8a0c61SJames Wright 
54bb8a0c61SJames Wright   PetscScalar meter  = user->units->meter;
55bb8a0c61SJames Wright   PetscScalar second = user->units->second;
56bb8a0c61SJames Wright   PetscScalar Kelvin = user->units->Kelvin;
57bb8a0c61SJames Wright   PetscScalar Pascal = user->units->Pascal;
58bb8a0c61SJames Wright 
59bb8a0c61SJames Wright   theta0 *= Kelvin;
60bb8a0c61SJames Wright   P0 *= Pascal;
61bb8a0c61SJames Wright   umax *= meter / second;
62bb8a0c61SJames Wright 
63bb8a0c61SJames Wright   //-- Setup Problem information
64bb8a0c61SJames Wright   CeedScalar H, center;
65bb8a0c61SJames Wright   {
66bb8a0c61SJames Wright     PetscReal domain_min[3], domain_max[3], domain_size[3];
672b916ea7SJeremy L Thompson     PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
68493642f1SJames Wright     for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
69bb8a0c61SJames Wright 
70bb8a0c61SJames Wright     H      = 0.5 * domain_size[1] * meter;
71bb8a0c61SJames Wright     center = H + domain_min[1] * meter;
72bb8a0c61SJames Wright   }
73bb8a0c61SJames Wright 
7415a3537eSJed Brown   // Some properties depend on parameters from NewtonianIdealGas
752b916ea7SJeremy L Thompson   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx);
7615a3537eSJed Brown 
7715a3537eSJed Brown   channel_ctx->center   = center;
7815a3537eSJed Brown   channel_ctx->H        = H;
7915a3537eSJed Brown   channel_ctx->theta0   = theta0;
8015a3537eSJed Brown   channel_ctx->P0       = P0;
8115a3537eSJed Brown   channel_ctx->umax     = umax;
8215a3537eSJed Brown   channel_ctx->implicit = user->phys->implicit;
8310786903SJed Brown   channel_ctx->B        = body_force_scale * 2 * umax * newtonian_ig_ctx->mu / (H * H);
84bb8a0c61SJames Wright 
85bb8a0c61SJames Wright   {
86bb8a0c61SJames Wright     // Calculate Body force
872b916ea7SJeremy L Thompson     CeedScalar cv = newtonian_ig_ctx->cv, cp = newtonian_ig_ctx->cp;
88bb8a0c61SJames Wright     CeedScalar Rd  = cp - cv;
89bb8a0c61SJames Wright     CeedScalar rho = P0 / (Rd * theta0);
9015a3537eSJed Brown     CeedScalar g[] = {channel_ctx->B / rho, 0., 0.};
912b916ea7SJeremy L Thompson     PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3));
92bb8a0c61SJames Wright   }
9315a3537eSJed Brown   channel_ctx->newtonian_ctx = *newtonian_ig_ctx;
942b916ea7SJeremy L Thompson   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx);
95bb8a0c61SJames Wright 
9615a3537eSJed Brown   CeedQFunctionContextCreate(user->ceed, &channel_context);
972b916ea7SJeremy L Thompson   CeedQFunctionContextSetData(channel_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*channel_ctx), channel_ctx);
982b916ea7SJeremy L Thompson   CeedQFunctionContextSetDataDestroy(channel_context, CEED_MEM_HOST, FreeContextPetsc);
9915a3537eSJed Brown 
10015a3537eSJed Brown   problem->ics.qfunction_context = channel_context;
1012b916ea7SJeremy L Thompson   CeedQFunctionContextReferenceCopy(channel_context, &problem->apply_inflow.qfunction_context);
1022b916ea7SJeremy L Thompson   CeedQFunctionContextReferenceCopy(channel_context, &problem->apply_outflow.qfunction_context);
103*d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
104bb8a0c61SJames Wright }
105