xref: /libCEED/examples/fluids/problems/densitycurrent.c (revision ee4ca9cbfe2be39196684117442f3ce8d00b6506)
1ea61e9acSJeremy L Thompson // Copyright (c) 2017-2022, 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 
1846de7363SJames Wright PetscErrorCode NS_DENSITY_CURRENT(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) {
1977841947SLeila Ghaffari   MPI_Comm                 comm = PETSC_COMM_WORLD;
20dc805cc4SLeila Ghaffari   User                     user = *(User *)ctx;
21dc805cc4SLeila Ghaffari   DensityCurrentContext    dc_ctx;
22dc805cc4SLeila Ghaffari   CeedQFunctionContext     density_current_context;
23dc805cc4SLeila Ghaffari   NewtonianIdealGasContext newtonian_ig_ctx;
2477841947SLeila Ghaffari 
25a0add3c9SJed Brown   PetscFunctionBeginUser;
2646de7363SJames Wright   PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc));
272b730f8bSJeremy L Thompson   PetscCall(PetscCalloc1(1, &dc_ctx));
2877841947SLeila Ghaffari   // ------------------------------------------------------
2977841947SLeila Ghaffari   //               SET UP DENSITY_CURRENT
3077841947SLeila Ghaffari   // ------------------------------------------------------
31dc805cc4SLeila Ghaffari   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
3291e5af17SJed Brown   problem->ics.qfunction     = ICsDC;
3391e5af17SJed Brown   problem->ics.qfunction_loc = ICsDC_loc;
3477841947SLeila Ghaffari 
3577841947SLeila Ghaffari   // ------------------------------------------------------
3677841947SLeila Ghaffari   //             Create the libCEED context
3777841947SLeila Ghaffari   // ------------------------------------------------------
3888626eedSJames Wright   CeedScalar theta0 = 300.;   // K
3988626eedSJames Wright   CeedScalar thetaC = -15.;   // K
4088626eedSJames Wright   CeedScalar P0     = 1.e5;   // Pa
4188626eedSJames Wright   CeedScalar N      = 0.01;   // 1/s
4277841947SLeila Ghaffari   CeedScalar rc     = 1000.;  // m (Radius of bubble)
4377841947SLeila Ghaffari   PetscReal  center[3], dc_axis[3] = {0, 0, 0};
441864f1c2SLeila Ghaffari   PetscReal  domain_min[3], domain_max[3], domain_size[3];
452b730f8bSJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
462b730f8bSJeremy L Thompson   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
4777841947SLeila Ghaffari 
4877841947SLeila Ghaffari   // ------------------------------------------------------
4977841947SLeila Ghaffari   //              Command line Options
5077841947SLeila Ghaffari   // ------------------------------------------------------
5167490bc6SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for DENSITY_CURRENT problem", NULL);
522b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-theta0", "Reference potential temperature", NULL, theta0, &theta0, NULL));
532b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-thetaC", "Perturbation of potential temperature", NULL, thetaC, &thetaC, NULL));
542b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-P0", "Atmospheric pressure", NULL, P0, &P0, NULL));
552b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-N", "Brunt-Vaisala frequency", NULL, N, &N, NULL));
562b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL));
572b730f8bSJeremy L Thompson   for (PetscInt i = 0; i < 3; i++) center[i] = .5 * domain_size[i];
5877841947SLeila Ghaffari   PetscInt n = problem->dim;
592b730f8bSJeremy L Thompson   PetscCall(PetscOptionsRealArray("-center", "Location of bubble center", NULL, center, &n, NULL));
6077841947SLeila Ghaffari   n = problem->dim;
612b730f8bSJeremy L Thompson   PetscCall(PetscOptionsRealArray("-dc_axis",
6288626eedSJames Wright                                   "Axis of density current cylindrical anomaly, "
6388626eedSJames Wright                                   "or {0,0,0} for spherically symmetric",
642b730f8bSJeremy L Thompson                                   NULL, dc_axis, &n, NULL));
6577841947SLeila Ghaffari   {
662b730f8bSJeremy L Thompson     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) + PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2]));
6777841947SLeila Ghaffari     if (norm > 0) {
682b730f8bSJeremy L Thompson       for (PetscInt i = 0; i < 3; i++) dc_axis[i] /= norm;
6977841947SLeila Ghaffari     }
7077841947SLeila Ghaffari   }
7177841947SLeila Ghaffari 
7267490bc6SJeremy L Thompson   PetscOptionsEnd();
7377841947SLeila Ghaffari 
7488b783a1SJames Wright   PetscScalar meter  = user->units->meter;
7588626eedSJames Wright   PetscScalar second = user->units->second;
7688626eedSJames Wright   PetscScalar Kelvin = user->units->Kelvin;
7788626eedSJames Wright   PetscScalar Pascal = user->units->Pascal;
7877841947SLeila Ghaffari   rc                 = fabs(rc) * meter;
7988626eedSJames Wright   theta0 *= Kelvin;
8088626eedSJames Wright   thetaC *= Kelvin;
8188626eedSJames Wright   P0 *= Pascal;
8288626eedSJames Wright   N *= (1. / second);
832b730f8bSJeremy L Thompson   for (PetscInt i = 0; i < 3; i++) center[i] *= meter;
8477841947SLeila Ghaffari 
85dc805cc4SLeila Ghaffari   dc_ctx->theta0     = theta0;
86dc805cc4SLeila Ghaffari   dc_ctx->thetaC     = thetaC;
87dc805cc4SLeila Ghaffari   dc_ctx->P0         = P0;
88dc805cc4SLeila Ghaffari   dc_ctx->N          = N;
89dc805cc4SLeila Ghaffari   dc_ctx->rc         = rc;
90dc805cc4SLeila Ghaffari   dc_ctx->center[0]  = center[0];
91dc805cc4SLeila Ghaffari   dc_ctx->center[1]  = center[1];
92dc805cc4SLeila Ghaffari   dc_ctx->center[2]  = center[2];
93dc805cc4SLeila Ghaffari   dc_ctx->dc_axis[0] = dc_axis[0];
94dc805cc4SLeila Ghaffari   dc_ctx->dc_axis[1] = dc_axis[1];
95dc805cc4SLeila Ghaffari   dc_ctx->dc_axis[2] = dc_axis[2];
9677841947SLeila Ghaffari 
972b730f8bSJeremy L Thompson   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx);
98dc805cc4SLeila Ghaffari   dc_ctx->newtonian_ctx = *newtonian_ig_ctx;
992b730f8bSJeremy L Thompson   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx);
100dc805cc4SLeila Ghaffari   CeedQFunctionContextCreate(user->ceed, &density_current_context);
1012b730f8bSJeremy L Thompson   CeedQFunctionContextSetData(density_current_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*dc_ctx), dc_ctx);
1022b730f8bSJeremy L Thompson   CeedQFunctionContextSetDataDestroy(density_current_context, CEED_MEM_HOST, FreeContextPetsc);
103dc805cc4SLeila Ghaffari   problem->ics.qfunction_context = density_current_context;
10417ce10faSJeremy L Thompson 
105*ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
10677841947SLeila Ghaffari }
107