xref: /libCEED/examples/fluids/problems/densitycurrent.c (revision 9ba83ac0e4b1fca39d6fa6737a318a9f0cbc172d)
1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors.
2ea61e9acSJeremy 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 
1349aac155SJeremy L Thompson #include <ceed.h>
1449aac155SJeremy L Thompson #include <petscdm.h>
1549aac155SJeremy L Thompson 
1688626eedSJames Wright #include "../navierstokes.h"
1777841947SLeila Ghaffari 
18731c13d7SJames Wright PetscErrorCode NS_DENSITY_CURRENT(ProblemData problem, DM dm, void *ctx, SimpleBC bc) {
19dc805cc4SLeila Ghaffari   User                     user = *(User *)ctx;
20a424bcd0SJames Wright   MPI_Comm                 comm = user->comm;
21a424bcd0SJames Wright   Ceed                     ceed = user->ceed;
22dc805cc4SLeila Ghaffari   DensityCurrentContext    dc_ctx;
23dc805cc4SLeila Ghaffari   CeedQFunctionContext     density_current_context;
24dc805cc4SLeila Ghaffari   NewtonianIdealGasContext newtonian_ig_ctx;
2577841947SLeila Ghaffari 
26a0add3c9SJed Brown   PetscFunctionBeginUser;
2746de7363SJames Wright   PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc));
282b730f8bSJeremy L Thompson   PetscCall(PetscCalloc1(1, &dc_ctx));
2977841947SLeila Ghaffari   // ------------------------------------------------------
3077841947SLeila Ghaffari   //               SET UP DENSITY_CURRENT
3177841947SLeila Ghaffari   // ------------------------------------------------------
32a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfunction_context));
3391e5af17SJed Brown   problem->ics.qfunction     = ICsDC;
3491e5af17SJed Brown   problem->ics.qfunction_loc = ICsDC_loc;
3577841947SLeila Ghaffari 
3677841947SLeila Ghaffari   // ------------------------------------------------------
3777841947SLeila Ghaffari   //             Create the libCEED context
3877841947SLeila Ghaffari   // ------------------------------------------------------
3988626eedSJames Wright   CeedScalar theta0 = 300.;   // K
4088626eedSJames Wright   CeedScalar thetaC = -15.;   // K
4188626eedSJames Wright   CeedScalar P0     = 1.e5;   // Pa
4288626eedSJames Wright   CeedScalar N      = 0.01;   // 1/s
4377841947SLeila Ghaffari   CeedScalar rc     = 1000.;  // m (Radius of bubble)
4477841947SLeila Ghaffari   PetscReal  center[3], dc_axis[3] = {0, 0, 0};
451864f1c2SLeila Ghaffari   PetscReal  domain_min[3], domain_max[3], domain_size[3];
462b730f8bSJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
472b730f8bSJeremy L Thompson   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
4877841947SLeila Ghaffari 
4977841947SLeila Ghaffari   // ------------------------------------------------------
5077841947SLeila Ghaffari   //              Command line Options
5177841947SLeila Ghaffari   // ------------------------------------------------------
5267490bc6SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for DENSITY_CURRENT problem", NULL);
532b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-theta0", "Reference potential temperature", NULL, theta0, &theta0, NULL));
542b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-thetaC", "Perturbation of potential temperature", NULL, thetaC, &thetaC, NULL));
552b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-P0", "Atmospheric pressure", NULL, P0, &P0, NULL));
562b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-N", "Brunt-Vaisala frequency", NULL, N, &N, NULL));
572b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL));
582b730f8bSJeremy L Thompson   for (PetscInt i = 0; i < 3; i++) center[i] = .5 * domain_size[i];
5977841947SLeila Ghaffari   PetscInt n = problem->dim;
602b730f8bSJeremy L Thompson   PetscCall(PetscOptionsRealArray("-center", "Location of bubble center", NULL, center, &n, NULL));
6177841947SLeila Ghaffari   n = problem->dim;
622b730f8bSJeremy L Thompson   PetscCall(PetscOptionsRealArray("-dc_axis",
6388626eedSJames Wright                                   "Axis of density current cylindrical anomaly, "
6488626eedSJames Wright                                   "or {0,0,0} for spherically symmetric",
652b730f8bSJeremy L Thompson                                   NULL, dc_axis, &n, NULL));
6677841947SLeila Ghaffari   {
672b730f8bSJeremy L Thompson     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) + PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2]));
6877841947SLeila Ghaffari     if (norm > 0) {
692b730f8bSJeremy L Thompson       for (PetscInt i = 0; i < 3; i++) dc_axis[i] /= norm;
7077841947SLeila Ghaffari     }
7177841947SLeila Ghaffari   }
7277841947SLeila Ghaffari 
7367490bc6SJeremy L Thompson   PetscOptionsEnd();
7477841947SLeila Ghaffari 
7588b783a1SJames Wright   PetscScalar meter  = user->units->meter;
7688626eedSJames Wright   PetscScalar second = user->units->second;
7788626eedSJames Wright   PetscScalar Kelvin = user->units->Kelvin;
7888626eedSJames Wright   PetscScalar Pascal = user->units->Pascal;
7977841947SLeila Ghaffari   rc                 = fabs(rc) * meter;
8088626eedSJames Wright   theta0 *= Kelvin;
8188626eedSJames Wright   thetaC *= Kelvin;
8288626eedSJames Wright   P0 *= Pascal;
8388626eedSJames Wright   N *= (1. / second);
842b730f8bSJeremy L Thompson   for (PetscInt i = 0; i < 3; i++) center[i] *= meter;
8577841947SLeila Ghaffari 
86dc805cc4SLeila Ghaffari   dc_ctx->theta0     = theta0;
87dc805cc4SLeila Ghaffari   dc_ctx->thetaC     = thetaC;
88dc805cc4SLeila Ghaffari   dc_ctx->P0         = P0;
89dc805cc4SLeila Ghaffari   dc_ctx->N          = N;
90dc805cc4SLeila Ghaffari   dc_ctx->rc         = rc;
91dc805cc4SLeila Ghaffari   dc_ctx->center[0]  = center[0];
92dc805cc4SLeila Ghaffari   dc_ctx->center[1]  = center[1];
93dc805cc4SLeila Ghaffari   dc_ctx->center[2]  = center[2];
94dc805cc4SLeila Ghaffari   dc_ctx->dc_axis[0] = dc_axis[0];
95dc805cc4SLeila Ghaffari   dc_ctx->dc_axis[1] = dc_axis[1];
96dc805cc4SLeila Ghaffari   dc_ctx->dc_axis[2] = dc_axis[2];
9777841947SLeila Ghaffari 
98a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx));
99dc805cc4SLeila Ghaffari   dc_ctx->newtonian_ctx = *newtonian_ig_ctx;
100a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx));
101a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &density_current_context));
102a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(density_current_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*dc_ctx), dc_ctx));
103a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(density_current_context, CEED_MEM_HOST, FreeContextPetsc));
104dc805cc4SLeila Ghaffari   problem->ics.qfunction_context = density_current_context;
105ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
10677841947SLeila Ghaffari }
107