xref: /libCEED/examples/fluids/problems/densitycurrent.c (revision 2b730f8b5a9c809740a0b3b302db43a719c636b1)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other
2 // CEED contributors. All Rights Reserved. See the top-level LICENSE and NOTICE
3 // files for details.
4 //
5 // SPDX-License-Identifier: BSD-2-Clause
6 //
7 // This file is part of CEED:  http://github.com/ceed
8 
9 /// @file
10 /// Utility functions for setting up DENSITY_CURRENT
11 
12 #include "../qfunctions/densitycurrent.h"
13 
14 #include "../navierstokes.h"
15 
16 PetscErrorCode NS_DENSITY_CURRENT(ProblemData *problem, DM dm, void *ctx) {
17   MPI_Comm                 comm = PETSC_COMM_WORLD;
18   User                     user = *(User *)ctx;
19   DensityCurrentContext    dc_ctx;
20   CeedQFunctionContext     density_current_context;
21   NewtonianIdealGasContext newtonian_ig_ctx;
22 
23   PetscFunctionBeginUser;
24   PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx));
25   PetscCall(PetscCalloc1(1, &dc_ctx));
26   // ------------------------------------------------------
27   //               SET UP DENSITY_CURRENT
28   // ------------------------------------------------------
29   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
30   problem->ics.qfunction     = ICsDC;
31   problem->ics.qfunction_loc = ICsDC_loc;
32 
33   // ------------------------------------------------------
34   //             Create the libCEED context
35   // ------------------------------------------------------
36   CeedScalar theta0 = 300.;   // K
37   CeedScalar thetaC = -15.;   // K
38   CeedScalar P0     = 1.e5;   // Pa
39   CeedScalar N      = 0.01;   // 1/s
40   CeedScalar rc     = 1000.;  // m (Radius of bubble)
41   PetscReal  center[3], dc_axis[3] = {0, 0, 0};
42   PetscReal  domain_min[3], domain_max[3], domain_size[3];
43   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
44   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
45 
46   // ------------------------------------------------------
47   //              Command line Options
48   // ------------------------------------------------------
49   PetscOptionsBegin(comm, NULL, "Options for DENSITY_CURRENT problem", NULL);
50   PetscCall(PetscOptionsScalar("-theta0", "Reference potential temperature", NULL, theta0, &theta0, NULL));
51   PetscCall(PetscOptionsScalar("-thetaC", "Perturbation of potential temperature", NULL, thetaC, &thetaC, NULL));
52   PetscCall(PetscOptionsScalar("-P0", "Atmospheric pressure", NULL, P0, &P0, NULL));
53   PetscCall(PetscOptionsScalar("-N", "Brunt-Vaisala frequency", NULL, N, &N, NULL));
54   PetscCall(PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL));
55   for (PetscInt i = 0; i < 3; i++) center[i] = .5 * domain_size[i];
56   PetscInt n = problem->dim;
57   PetscCall(PetscOptionsRealArray("-center", "Location of bubble center", NULL, center, &n, NULL));
58   n = problem->dim;
59   PetscCall(PetscOptionsRealArray("-dc_axis",
60                                   "Axis of density current cylindrical anomaly, "
61                                   "or {0,0,0} for spherically symmetric",
62                                   NULL, dc_axis, &n, NULL));
63   {
64     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) + PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2]));
65     if (norm > 0) {
66       for (PetscInt i = 0; i < 3; i++) dc_axis[i] /= norm;
67     }
68   }
69 
70   PetscOptionsEnd();
71 
72   PetscScalar meter  = user->units->meter;
73   PetscScalar second = user->units->second;
74   PetscScalar Kelvin = user->units->Kelvin;
75   PetscScalar Pascal = user->units->Pascal;
76   rc                 = fabs(rc) * meter;
77   theta0 *= Kelvin;
78   thetaC *= Kelvin;
79   P0 *= Pascal;
80   N *= (1. / second);
81   for (PetscInt i = 0; i < 3; i++) center[i] *= meter;
82 
83   dc_ctx->theta0     = theta0;
84   dc_ctx->thetaC     = thetaC;
85   dc_ctx->P0         = P0;
86   dc_ctx->N          = N;
87   dc_ctx->rc         = rc;
88   dc_ctx->center[0]  = center[0];
89   dc_ctx->center[1]  = center[1];
90   dc_ctx->center[2]  = center[2];
91   dc_ctx->dc_axis[0] = dc_axis[0];
92   dc_ctx->dc_axis[1] = dc_axis[1];
93   dc_ctx->dc_axis[2] = dc_axis[2];
94 
95   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx);
96   dc_ctx->newtonian_ctx = *newtonian_ig_ctx;
97   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx);
98   CeedQFunctionContextCreate(user->ceed, &density_current_context);
99   CeedQFunctionContextSetData(density_current_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*dc_ctx), dc_ctx);
100   CeedQFunctionContextSetDataDestroy(density_current_context, CEED_MEM_HOST, FreeContextPetsc);
101   problem->ics.qfunction_context = density_current_context;
102 
103   PetscFunctionReturn(0);
104 }
105