xref: /libCEED/examples/fluids/problems/densitycurrent.c (revision a0add3c91c6ddcd77c0d376840911ac920dd4230)
188626eedSJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other
288626eedSJames Wright // CEED contributors. All Rights Reserved. See the top-level LICENSE and NOTICE
388626eedSJames Wright // files for details.
477841947SLeila Ghaffari //
53d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
677841947SLeila Ghaffari //
73d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
877841947SLeila Ghaffari 
977841947SLeila Ghaffari /// @file
1077841947SLeila Ghaffari /// Utility functions for setting up DENSITY_CURRENT
1177841947SLeila Ghaffari 
1277841947SLeila Ghaffari #include "../qfunctions/densitycurrent.h"
1388626eedSJames Wright #include "../navierstokes.h"
1477841947SLeila Ghaffari 
15*a0add3c9SJed Brown PetscErrorCode NS_DENSITY_CURRENT(ProblemData *problem, DM dm, void *ctx) {
1688b783a1SJames Wright 
1788b783a1SJames Wright   PetscInt ierr;
18*a0add3c9SJed Brown   SetupContext setup_context;
1977841947SLeila Ghaffari   User user = *(User *)ctx;
2077841947SLeila Ghaffari   MPI_Comm comm = PETSC_COMM_WORLD;
2177841947SLeila Ghaffari 
22*a0add3c9SJed Brown   PetscFunctionBeginUser;
23*a0add3c9SJed Brown   ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr);
2477841947SLeila Ghaffari   // ------------------------------------------------------
2577841947SLeila Ghaffari   //               SET UP DENSITY_CURRENT
2677841947SLeila Ghaffari   // ------------------------------------------------------
2791e5af17SJed Brown   problem->ics.qfunction = ICsDC;
2891e5af17SJed Brown   problem->ics.qfunction_loc = ICsDC_loc;
2977841947SLeila Ghaffari   problem->bc = Exact_DC;
30*a0add3c9SJed Brown   setup_context = problem->bc_ctx;
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];
4288626eedSJames Wright   ierr = DMGetBoundingBox(dm, domain_min, domain_max);
4388626eedSJames Wright   CHKERRQ(ierr);
4488626eedSJames Wright   for (int i = 0; i < 3; i++)
4588626eedSJames Wright     domain_size[i] = domain_max[i] - domain_min[i];
4677841947SLeila Ghaffari 
4777841947SLeila Ghaffari   // ------------------------------------------------------
4877841947SLeila Ghaffari   //              Command line Options
4977841947SLeila Ghaffari   // ------------------------------------------------------
5067490bc6SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for DENSITY_CURRENT problem", NULL);
5188626eedSJames Wright   ierr = PetscOptionsScalar("-theta0", "Reference potential temperature", NULL,
5288626eedSJames Wright                             theta0, &theta0, NULL);
5388626eedSJames Wright   CHKERRQ(ierr);
5488626eedSJames Wright   ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature",
5588626eedSJames Wright                             NULL, thetaC, &thetaC, NULL);
5688626eedSJames Wright   CHKERRQ(ierr);
5788626eedSJames Wright   ierr = PetscOptionsScalar("-P0", "Atmospheric pressure", NULL, P0, &P0, NULL);
5888626eedSJames Wright   CHKERRQ(ierr);
5988626eedSJames Wright   ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency", NULL, N, &N, NULL);
6088626eedSJames Wright   CHKERRQ(ierr);
6177841947SLeila Ghaffari   ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble",
6288626eedSJames Wright                             NULL, rc, &rc, NULL);
6388626eedSJames Wright   CHKERRQ(ierr);
6488626eedSJames Wright   for (int i = 0; i < 3; i++)
6588626eedSJames Wright     center[i] = .5 * domain_size[i];
6677841947SLeila Ghaffari   PetscInt n = problem->dim;
6788626eedSJames Wright   ierr = PetscOptionsRealArray("-center", "Location of bubble center", NULL,
6888626eedSJames Wright                                center, &n, NULL);
6988626eedSJames Wright   CHKERRQ(ierr);
7077841947SLeila Ghaffari   n = problem->dim;
7177841947SLeila Ghaffari   ierr = PetscOptionsRealArray("-dc_axis",
7288626eedSJames Wright                                "Axis of density current cylindrical anomaly, "
7388626eedSJames Wright                                "or {0,0,0} for spherically symmetric",
7488626eedSJames Wright                                NULL, dc_axis, &n, NULL);
7588626eedSJames Wright   CHKERRQ(ierr);
7677841947SLeila Ghaffari   {
7777841947SLeila Ghaffari     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) + PetscSqr(dc_axis[1]) +
7877841947SLeila Ghaffari                                    PetscSqr(dc_axis[2]));
7977841947SLeila Ghaffari     if (norm > 0) {
8088626eedSJames Wright       for (int i = 0; i < 3; i++)
8188626eedSJames Wright         dc_axis[i] /= norm;
8277841947SLeila Ghaffari     }
8377841947SLeila Ghaffari   }
8477841947SLeila Ghaffari 
8567490bc6SJeremy L Thompson   PetscOptionsEnd();
8677841947SLeila Ghaffari 
8788b783a1SJames Wright   PetscScalar meter           = user->units->meter;
8888626eedSJames Wright   PetscScalar second          = user->units->second;
8988626eedSJames Wright   PetscScalar Kelvin          = user->units->Kelvin;
9088626eedSJames Wright   PetscScalar Pascal          = user->units->Pascal;
9177841947SLeila Ghaffari   rc = fabs(rc) * meter;
9288626eedSJames Wright   theta0 *= Kelvin;
9388626eedSJames Wright   thetaC *= Kelvin;
9488626eedSJames Wright   P0 *= Pascal;
9588626eedSJames Wright   N *= (1. / second);
9688626eedSJames Wright   for (int i = 0; i < 3; i++)
9788626eedSJames Wright     center[i] *= meter;
9877841947SLeila Ghaffari 
9988626eedSJames Wright   setup_context->theta0 = theta0;
10088626eedSJames Wright   setup_context->thetaC = thetaC;
10188626eedSJames Wright   setup_context->P0 = P0;
10288626eedSJames Wright   setup_context->N = N;
10377841947SLeila Ghaffari   setup_context->rc = rc;
10477841947SLeila Ghaffari   setup_context->center[0] = center[0];
10577841947SLeila Ghaffari   setup_context->center[1] = center[1];
10677841947SLeila Ghaffari   setup_context->center[2] = center[2];
10777841947SLeila Ghaffari   setup_context->dc_axis[0] = dc_axis[0];
10877841947SLeila Ghaffari   setup_context->dc_axis[1] = dc_axis[1];
10977841947SLeila Ghaffari   setup_context->dc_axis[2] = dc_axis[2];
11077841947SLeila Ghaffari 
11177841947SLeila Ghaffari   PetscFunctionReturn(0);
11277841947SLeila Ghaffari }
113