xref: /libCEED/examples/fluids/problems/densitycurrent.c (revision ea61e9ac44808524e4667c1525a05976f536c19c)
1*ea61e9acSJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2*ea61e9acSJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
377841947SLeila Ghaffari //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
577841947SLeila Ghaffari //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
777841947SLeila Ghaffari 
877841947SLeila Ghaffari /// @file
977841947SLeila Ghaffari /// Utility functions for setting up DENSITY_CURRENT
1077841947SLeila Ghaffari 
1177841947SLeila Ghaffari #include "../qfunctions/densitycurrent.h"
122b730f8bSJeremy L Thompson 
1388626eedSJames Wright #include "../navierstokes.h"
1477841947SLeila Ghaffari 
1546de7363SJames Wright PetscErrorCode NS_DENSITY_CURRENT(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) {
1677841947SLeila Ghaffari   MPI_Comm                 comm = PETSC_COMM_WORLD;
17dc805cc4SLeila Ghaffari   User                     user = *(User *)ctx;
18dc805cc4SLeila Ghaffari   DensityCurrentContext    dc_ctx;
19dc805cc4SLeila Ghaffari   CeedQFunctionContext     density_current_context;
20dc805cc4SLeila Ghaffari   NewtonianIdealGasContext newtonian_ig_ctx;
2177841947SLeila Ghaffari 
22a0add3c9SJed Brown   PetscFunctionBeginUser;
2346de7363SJames Wright   PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc));
242b730f8bSJeremy L Thompson   PetscCall(PetscCalloc1(1, &dc_ctx));
2577841947SLeila Ghaffari   // ------------------------------------------------------
2677841947SLeila Ghaffari   //               SET UP DENSITY_CURRENT
2777841947SLeila Ghaffari   // ------------------------------------------------------
28dc805cc4SLeila Ghaffari   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
2991e5af17SJed Brown   problem->ics.qfunction     = ICsDC;
3091e5af17SJed Brown   problem->ics.qfunction_loc = ICsDC_loc;
3177841947SLeila Ghaffari 
3277841947SLeila Ghaffari   // ------------------------------------------------------
3377841947SLeila Ghaffari   //             Create the libCEED context
3477841947SLeila Ghaffari   // ------------------------------------------------------
3588626eedSJames Wright   CeedScalar theta0 = 300.;   // K
3688626eedSJames Wright   CeedScalar thetaC = -15.;   // K
3788626eedSJames Wright   CeedScalar P0     = 1.e5;   // Pa
3888626eedSJames Wright   CeedScalar N      = 0.01;   // 1/s
3977841947SLeila Ghaffari   CeedScalar rc     = 1000.;  // m (Radius of bubble)
4077841947SLeila Ghaffari   PetscReal  center[3], dc_axis[3] = {0, 0, 0};
411864f1c2SLeila Ghaffari   PetscReal  domain_min[3], domain_max[3], domain_size[3];
422b730f8bSJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
432b730f8bSJeremy L Thompson   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
4477841947SLeila Ghaffari 
4577841947SLeila Ghaffari   // ------------------------------------------------------
4677841947SLeila Ghaffari   //              Command line Options
4777841947SLeila Ghaffari   // ------------------------------------------------------
4867490bc6SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for DENSITY_CURRENT problem", NULL);
492b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-theta0", "Reference potential temperature", NULL, theta0, &theta0, NULL));
502b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-thetaC", "Perturbation of potential temperature", NULL, thetaC, &thetaC, NULL));
512b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-P0", "Atmospheric pressure", NULL, P0, &P0, NULL));
522b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-N", "Brunt-Vaisala frequency", NULL, N, &N, NULL));
532b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL));
542b730f8bSJeremy L Thompson   for (PetscInt i = 0; i < 3; i++) center[i] = .5 * domain_size[i];
5577841947SLeila Ghaffari   PetscInt n = problem->dim;
562b730f8bSJeremy L Thompson   PetscCall(PetscOptionsRealArray("-center", "Location of bubble center", NULL, center, &n, NULL));
5777841947SLeila Ghaffari   n = problem->dim;
582b730f8bSJeremy L Thompson   PetscCall(PetscOptionsRealArray("-dc_axis",
5988626eedSJames Wright                                   "Axis of density current cylindrical anomaly, "
6088626eedSJames Wright                                   "or {0,0,0} for spherically symmetric",
612b730f8bSJeremy L Thompson                                   NULL, dc_axis, &n, NULL));
6277841947SLeila Ghaffari   {
632b730f8bSJeremy L Thompson     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) + PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2]));
6477841947SLeila Ghaffari     if (norm > 0) {
652b730f8bSJeremy L Thompson       for (PetscInt i = 0; i < 3; i++) dc_axis[i] /= norm;
6677841947SLeila Ghaffari     }
6777841947SLeila Ghaffari   }
6877841947SLeila Ghaffari 
6967490bc6SJeremy L Thompson   PetscOptionsEnd();
7077841947SLeila Ghaffari 
7188b783a1SJames Wright   PetscScalar meter  = user->units->meter;
7288626eedSJames Wright   PetscScalar second = user->units->second;
7388626eedSJames Wright   PetscScalar Kelvin = user->units->Kelvin;
7488626eedSJames Wright   PetscScalar Pascal = user->units->Pascal;
7577841947SLeila Ghaffari   rc                 = fabs(rc) * meter;
7688626eedSJames Wright   theta0 *= Kelvin;
7788626eedSJames Wright   thetaC *= Kelvin;
7888626eedSJames Wright   P0 *= Pascal;
7988626eedSJames Wright   N *= (1. / second);
802b730f8bSJeremy L Thompson   for (PetscInt i = 0; i < 3; i++) center[i] *= meter;
8177841947SLeila Ghaffari 
82dc805cc4SLeila Ghaffari   dc_ctx->theta0     = theta0;
83dc805cc4SLeila Ghaffari   dc_ctx->thetaC     = thetaC;
84dc805cc4SLeila Ghaffari   dc_ctx->P0         = P0;
85dc805cc4SLeila Ghaffari   dc_ctx->N          = N;
86dc805cc4SLeila Ghaffari   dc_ctx->rc         = rc;
87dc805cc4SLeila Ghaffari   dc_ctx->center[0]  = center[0];
88dc805cc4SLeila Ghaffari   dc_ctx->center[1]  = center[1];
89dc805cc4SLeila Ghaffari   dc_ctx->center[2]  = center[2];
90dc805cc4SLeila Ghaffari   dc_ctx->dc_axis[0] = dc_axis[0];
91dc805cc4SLeila Ghaffari   dc_ctx->dc_axis[1] = dc_axis[1];
92dc805cc4SLeila Ghaffari   dc_ctx->dc_axis[2] = dc_axis[2];
9377841947SLeila Ghaffari 
942b730f8bSJeremy L Thompson   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx);
95dc805cc4SLeila Ghaffari   dc_ctx->newtonian_ctx = *newtonian_ig_ctx;
962b730f8bSJeremy L Thompson   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx);
97dc805cc4SLeila Ghaffari   CeedQFunctionContextCreate(user->ceed, &density_current_context);
982b730f8bSJeremy L Thompson   CeedQFunctionContextSetData(density_current_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*dc_ctx), dc_ctx);
992b730f8bSJeremy L Thompson   CeedQFunctionContextSetDataDestroy(density_current_context, CEED_MEM_HOST, FreeContextPetsc);
100dc805cc4SLeila Ghaffari   problem->ics.qfunction_context = density_current_context;
10117ce10faSJeremy L Thompson 
10277841947SLeila Ghaffari   PetscFunctionReturn(0);
10377841947SLeila Ghaffari }
104