xref: /honee/problems/densitycurrent.c (revision 04e40bb60650195adcc92556a3eb81ec7887ccc8)
1*04e40bb6SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2*04e40bb6SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3a515125bSLeila Ghaffari //
4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5a515125bSLeila Ghaffari //
6727da7e7SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7a515125bSLeila Ghaffari 
8a515125bSLeila Ghaffari /// @file
9a515125bSLeila Ghaffari /// Utility functions for setting up DENSITY_CURRENT
10a515125bSLeila Ghaffari 
11a515125bSLeila Ghaffari #include "../qfunctions/densitycurrent.h"
122b916ea7SJeremy L Thompson 
13bb8a0c61SJames Wright #include "../navierstokes.h"
14a515125bSLeila Ghaffari 
15d1c51a42SJames Wright PetscErrorCode NS_DENSITY_CURRENT(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) {
16a515125bSLeila Ghaffari   MPI_Comm                 comm = PETSC_COMM_WORLD;
17cbe60e31SLeila Ghaffari   User                     user = *(User *)ctx;
18cbe60e31SLeila Ghaffari   DensityCurrentContext    dc_ctx;
19cbe60e31SLeila Ghaffari   CeedQFunctionContext     density_current_context;
20cbe60e31SLeila Ghaffari   NewtonianIdealGasContext newtonian_ig_ctx;
21a515125bSLeila Ghaffari 
22b7f03f12SJed Brown   PetscFunctionBeginUser;
23d1c51a42SJames Wright   PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc));
242b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &dc_ctx));
25a515125bSLeila Ghaffari   // ------------------------------------------------------
26a515125bSLeila Ghaffari   //               SET UP DENSITY_CURRENT
27a515125bSLeila Ghaffari   // ------------------------------------------------------
28cbe60e31SLeila Ghaffari   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
299785fe93SJed Brown   problem->ics.qfunction     = ICsDC;
309785fe93SJed Brown   problem->ics.qfunction_loc = ICsDC_loc;
31a515125bSLeila Ghaffari 
32a515125bSLeila Ghaffari   // ------------------------------------------------------
33a515125bSLeila Ghaffari   //             Create the libCEED context
34a515125bSLeila Ghaffari   // ------------------------------------------------------
35bb8a0c61SJames Wright   CeedScalar theta0 = 300.;   // K
36bb8a0c61SJames Wright   CeedScalar thetaC = -15.;   // K
37bb8a0c61SJames Wright   CeedScalar P0     = 1.e5;   // Pa
38bb8a0c61SJames Wright   CeedScalar N      = 0.01;   // 1/s
39a515125bSLeila Ghaffari   CeedScalar rc     = 1000.;  // m (Radius of bubble)
40a515125bSLeila Ghaffari   PetscReal  center[3], dc_axis[3] = {0, 0, 0};
4105a512bdSLeila Ghaffari   PetscReal  domain_min[3], domain_max[3], domain_size[3];
422b916ea7SJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
432b916ea7SJeremy L Thompson   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
44a515125bSLeila Ghaffari 
45a515125bSLeila Ghaffari   // ------------------------------------------------------
46a515125bSLeila Ghaffari   //              Command line Options
47a515125bSLeila Ghaffari   // ------------------------------------------------------
481485969bSJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for DENSITY_CURRENT problem", NULL);
492b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-theta0", "Reference potential temperature", NULL, theta0, &theta0, NULL));
502b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-thetaC", "Perturbation of potential temperature", NULL, thetaC, &thetaC, NULL));
512b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-P0", "Atmospheric pressure", NULL, P0, &P0, NULL));
522b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-N", "Brunt-Vaisala frequency", NULL, N, &N, NULL));
532b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL));
542b916ea7SJeremy L Thompson   for (PetscInt i = 0; i < 3; i++) center[i] = .5 * domain_size[i];
55a515125bSLeila Ghaffari   PetscInt n = problem->dim;
562b916ea7SJeremy L Thompson   PetscCall(PetscOptionsRealArray("-center", "Location of bubble center", NULL, center, &n, NULL));
57a515125bSLeila Ghaffari   n = problem->dim;
582b916ea7SJeremy L Thompson   PetscCall(PetscOptionsRealArray("-dc_axis",
59bb8a0c61SJames Wright                                   "Axis of density current cylindrical anomaly, "
60bb8a0c61SJames Wright                                   "or {0,0,0} for spherically symmetric",
612b916ea7SJeremy L Thompson                                   NULL, dc_axis, &n, NULL));
62a515125bSLeila Ghaffari   {
632b916ea7SJeremy L Thompson     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) + PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2]));
64a515125bSLeila Ghaffari     if (norm > 0) {
652b916ea7SJeremy L Thompson       for (PetscInt i = 0; i < 3; i++) dc_axis[i] /= norm;
66a515125bSLeila Ghaffari     }
67a515125bSLeila Ghaffari   }
68a515125bSLeila Ghaffari 
691485969bSJeremy L Thompson   PetscOptionsEnd();
70a515125bSLeila Ghaffari 
713a8779fbSJames Wright   PetscScalar meter  = user->units->meter;
72bb8a0c61SJames Wright   PetscScalar second = user->units->second;
73bb8a0c61SJames Wright   PetscScalar Kelvin = user->units->Kelvin;
74bb8a0c61SJames Wright   PetscScalar Pascal = user->units->Pascal;
75a515125bSLeila Ghaffari   rc                 = fabs(rc) * meter;
76bb8a0c61SJames Wright   theta0 *= Kelvin;
77bb8a0c61SJames Wright   thetaC *= Kelvin;
78bb8a0c61SJames Wright   P0 *= Pascal;
79bb8a0c61SJames Wright   N *= (1. / second);
802b916ea7SJeremy L Thompson   for (PetscInt i = 0; i < 3; i++) center[i] *= meter;
81a515125bSLeila Ghaffari 
82cbe60e31SLeila Ghaffari   dc_ctx->theta0     = theta0;
83cbe60e31SLeila Ghaffari   dc_ctx->thetaC     = thetaC;
84cbe60e31SLeila Ghaffari   dc_ctx->P0         = P0;
85cbe60e31SLeila Ghaffari   dc_ctx->N          = N;
86cbe60e31SLeila Ghaffari   dc_ctx->rc         = rc;
87cbe60e31SLeila Ghaffari   dc_ctx->center[0]  = center[0];
88cbe60e31SLeila Ghaffari   dc_ctx->center[1]  = center[1];
89cbe60e31SLeila Ghaffari   dc_ctx->center[2]  = center[2];
90cbe60e31SLeila Ghaffari   dc_ctx->dc_axis[0] = dc_axis[0];
91cbe60e31SLeila Ghaffari   dc_ctx->dc_axis[1] = dc_axis[1];
92cbe60e31SLeila Ghaffari   dc_ctx->dc_axis[2] = dc_axis[2];
93a515125bSLeila Ghaffari 
942b916ea7SJeremy L Thompson   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx);
95cbe60e31SLeila Ghaffari   dc_ctx->newtonian_ctx = *newtonian_ig_ctx;
962b916ea7SJeremy L Thompson   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx);
97cbe60e31SLeila Ghaffari   CeedQFunctionContextCreate(user->ceed, &density_current_context);
982b916ea7SJeremy L Thompson   CeedQFunctionContextSetData(density_current_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*dc_ctx), dc_ctx);
992b916ea7SJeremy L Thompson   CeedQFunctionContextSetDataDestroy(density_current_context, CEED_MEM_HOST, FreeContextPetsc);
100cbe60e31SLeila Ghaffari   problem->ics.qfunction_context = density_current_context;
101270bbb13SJeremy L Thompson 
102a515125bSLeila Ghaffari   PetscFunctionReturn(0);
103a515125bSLeila Ghaffari }
104