xref: /honee/problems/densitycurrent.c (revision d949ddfc92e3b3e8aad90700e3fb60ca7f4585db)
104e40bb6SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
204e40bb6SJeremy 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 
13e419654dSJeremy L Thompson #include <ceed.h>
14e419654dSJeremy L Thompson #include <petscdm.h>
15e419654dSJeremy L Thompson 
16bb8a0c61SJames Wright #include "../navierstokes.h"
17a515125bSLeila Ghaffari 
18d1c51a42SJames Wright PetscErrorCode NS_DENSITY_CURRENT(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) {
19a515125bSLeila Ghaffari   MPI_Comm                 comm = PETSC_COMM_WORLD;
20cbe60e31SLeila Ghaffari   User                     user = *(User *)ctx;
21cbe60e31SLeila Ghaffari   DensityCurrentContext    dc_ctx;
22cbe60e31SLeila Ghaffari   CeedQFunctionContext     density_current_context;
23cbe60e31SLeila Ghaffari   NewtonianIdealGasContext newtonian_ig_ctx;
24a515125bSLeila Ghaffari 
25b7f03f12SJed Brown   PetscFunctionBeginUser;
26d1c51a42SJames Wright   PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc));
272b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &dc_ctx));
28a515125bSLeila Ghaffari   // ------------------------------------------------------
29a515125bSLeila Ghaffari   //               SET UP DENSITY_CURRENT
30a515125bSLeila Ghaffari   // ------------------------------------------------------
31cbe60e31SLeila Ghaffari   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
329785fe93SJed Brown   problem->ics.qfunction     = ICsDC;
339785fe93SJed Brown   problem->ics.qfunction_loc = ICsDC_loc;
34a515125bSLeila Ghaffari 
35a515125bSLeila Ghaffari   // ------------------------------------------------------
36a515125bSLeila Ghaffari   //             Create the libCEED context
37a515125bSLeila Ghaffari   // ------------------------------------------------------
38bb8a0c61SJames Wright   CeedScalar theta0 = 300.;   // K
39bb8a0c61SJames Wright   CeedScalar thetaC = -15.;   // K
40bb8a0c61SJames Wright   CeedScalar P0     = 1.e5;   // Pa
41bb8a0c61SJames Wright   CeedScalar N      = 0.01;   // 1/s
42a515125bSLeila Ghaffari   CeedScalar rc     = 1000.;  // m (Radius of bubble)
43a515125bSLeila Ghaffari   PetscReal  center[3], dc_axis[3] = {0, 0, 0};
4405a512bdSLeila Ghaffari   PetscReal  domain_min[3], domain_max[3], domain_size[3];
452b916ea7SJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
462b916ea7SJeremy L Thompson   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
47a515125bSLeila Ghaffari 
48a515125bSLeila Ghaffari   // ------------------------------------------------------
49a515125bSLeila Ghaffari   //              Command line Options
50a515125bSLeila Ghaffari   // ------------------------------------------------------
511485969bSJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for DENSITY_CURRENT problem", NULL);
522b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-theta0", "Reference potential temperature", NULL, theta0, &theta0, NULL));
532b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-thetaC", "Perturbation of potential temperature", NULL, thetaC, &thetaC, NULL));
542b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-P0", "Atmospheric pressure", NULL, P0, &P0, NULL));
552b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-N", "Brunt-Vaisala frequency", NULL, N, &N, NULL));
562b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL));
572b916ea7SJeremy L Thompson   for (PetscInt i = 0; i < 3; i++) center[i] = .5 * domain_size[i];
58a515125bSLeila Ghaffari   PetscInt n = problem->dim;
592b916ea7SJeremy L Thompson   PetscCall(PetscOptionsRealArray("-center", "Location of bubble center", NULL, center, &n, NULL));
60a515125bSLeila Ghaffari   n = problem->dim;
612b916ea7SJeremy L Thompson   PetscCall(PetscOptionsRealArray("-dc_axis",
62bb8a0c61SJames Wright                                   "Axis of density current cylindrical anomaly, "
63bb8a0c61SJames Wright                                   "or {0,0,0} for spherically symmetric",
642b916ea7SJeremy L Thompson                                   NULL, dc_axis, &n, NULL));
65a515125bSLeila Ghaffari   {
662b916ea7SJeremy L Thompson     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) + PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2]));
67a515125bSLeila Ghaffari     if (norm > 0) {
682b916ea7SJeremy L Thompson       for (PetscInt i = 0; i < 3; i++) dc_axis[i] /= norm;
69a515125bSLeila Ghaffari     }
70a515125bSLeila Ghaffari   }
71a515125bSLeila Ghaffari 
721485969bSJeremy L Thompson   PetscOptionsEnd();
73a515125bSLeila Ghaffari 
743a8779fbSJames Wright   PetscScalar meter  = user->units->meter;
75bb8a0c61SJames Wright   PetscScalar second = user->units->second;
76bb8a0c61SJames Wright   PetscScalar Kelvin = user->units->Kelvin;
77bb8a0c61SJames Wright   PetscScalar Pascal = user->units->Pascal;
78a515125bSLeila Ghaffari   rc                 = fabs(rc) * meter;
79bb8a0c61SJames Wright   theta0 *= Kelvin;
80bb8a0c61SJames Wright   thetaC *= Kelvin;
81bb8a0c61SJames Wright   P0 *= Pascal;
82bb8a0c61SJames Wright   N *= (1. / second);
832b916ea7SJeremy L Thompson   for (PetscInt i = 0; i < 3; i++) center[i] *= meter;
84a515125bSLeila Ghaffari 
85cbe60e31SLeila Ghaffari   dc_ctx->theta0     = theta0;
86cbe60e31SLeila Ghaffari   dc_ctx->thetaC     = thetaC;
87cbe60e31SLeila Ghaffari   dc_ctx->P0         = P0;
88cbe60e31SLeila Ghaffari   dc_ctx->N          = N;
89cbe60e31SLeila Ghaffari   dc_ctx->rc         = rc;
90cbe60e31SLeila Ghaffari   dc_ctx->center[0]  = center[0];
91cbe60e31SLeila Ghaffari   dc_ctx->center[1]  = center[1];
92cbe60e31SLeila Ghaffari   dc_ctx->center[2]  = center[2];
93cbe60e31SLeila Ghaffari   dc_ctx->dc_axis[0] = dc_axis[0];
94cbe60e31SLeila Ghaffari   dc_ctx->dc_axis[1] = dc_axis[1];
95cbe60e31SLeila Ghaffari   dc_ctx->dc_axis[2] = dc_axis[2];
96a515125bSLeila Ghaffari 
972b916ea7SJeremy L Thompson   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx);
98cbe60e31SLeila Ghaffari   dc_ctx->newtonian_ctx = *newtonian_ig_ctx;
992b916ea7SJeremy L Thompson   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx);
100cbe60e31SLeila Ghaffari   CeedQFunctionContextCreate(user->ceed, &density_current_context);
1012b916ea7SJeremy L Thompson   CeedQFunctionContextSetData(density_current_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*dc_ctx), dc_ctx);
1022b916ea7SJeremy L Thompson   CeedQFunctionContextSetDataDestroy(density_current_context, CEED_MEM_HOST, FreeContextPetsc);
103cbe60e31SLeila Ghaffari   problem->ics.qfunction_context = density_current_context;
104270bbb13SJeremy L Thompson 
105*d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
106a515125bSLeila Ghaffari }
107