xref: /libCEED/examples/fluids/problems/densitycurrent.c (revision dc805cc4a09d29f27b3febd084feb659e74b9d08)
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 
15a0add3c9SJed Brown PetscErrorCode NS_DENSITY_CURRENT(ProblemData *problem, DM dm, void *ctx) {
1688b783a1SJames Wright 
1788b783a1SJames Wright   PetscInt                 ierr;
1877841947SLeila Ghaffari   MPI_Comm                 comm = PETSC_COMM_WORLD;
19*dc805cc4SLeila Ghaffari   User                     user = *(User *)ctx;
20*dc805cc4SLeila Ghaffari   DensityCurrentContext    dc_ctx;
21*dc805cc4SLeila Ghaffari   CeedQFunctionContext     density_current_context;
22*dc805cc4SLeila Ghaffari   NewtonianIdealGasContext newtonian_ig_ctx;
2377841947SLeila Ghaffari 
24a0add3c9SJed Brown   PetscFunctionBeginUser;
25a0add3c9SJed Brown   ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr);
26*dc805cc4SLeila Ghaffari   ierr = PetscCalloc1(1, &dc_ctx); CHKERRQ(ierr);
2777841947SLeila Ghaffari   // ------------------------------------------------------
2877841947SLeila Ghaffari   //               SET UP DENSITY_CURRENT
2977841947SLeila Ghaffari   // ------------------------------------------------------
30*dc805cc4SLeila Ghaffari   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
3191e5af17SJed Brown   problem->ics.qfunction     = ICsDC;
3291e5af17SJed Brown   problem->ics.qfunction_loc = ICsDC_loc;
3377841947SLeila Ghaffari 
3477841947SLeila Ghaffari   // ------------------------------------------------------
3577841947SLeila Ghaffari   //             Create the libCEED context
3677841947SLeila Ghaffari   // ------------------------------------------------------
3788626eedSJames Wright   CeedScalar theta0 = 300.; // K
3888626eedSJames Wright   CeedScalar thetaC = -15.; // K
3988626eedSJames Wright   CeedScalar P0 = 1.e5;     // Pa
4088626eedSJames Wright   CeedScalar N = 0.01;      // 1/s
4177841947SLeila Ghaffari   CeedScalar rc = 1000.;    // m (Radius of bubble)
4277841947SLeila Ghaffari   PetscReal center[3], dc_axis[3] = {0, 0, 0};
431864f1c2SLeila Ghaffari   PetscReal domain_min[3], domain_max[3], domain_size[3];
4488626eedSJames Wright   ierr = DMGetBoundingBox(dm, domain_min, domain_max);
4588626eedSJames Wright   CHKERRQ(ierr);
46ba6664aeSJames Wright   for (PetscInt i = 0; i < 3; i++)
4788626eedSJames Wright     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);
5388626eedSJames Wright   ierr = PetscOptionsScalar("-theta0", "Reference potential temperature", NULL,
5488626eedSJames Wright                             theta0, &theta0, NULL);
5588626eedSJames Wright   CHKERRQ(ierr);
5688626eedSJames Wright   ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature",
5788626eedSJames Wright                             NULL, thetaC, &thetaC, NULL);
5888626eedSJames Wright   CHKERRQ(ierr);
5988626eedSJames Wright   ierr = PetscOptionsScalar("-P0", "Atmospheric pressure", NULL, P0, &P0, NULL);
6088626eedSJames Wright   CHKERRQ(ierr);
6188626eedSJames Wright   ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency", NULL, N, &N, NULL);
6288626eedSJames Wright   CHKERRQ(ierr);
6377841947SLeila Ghaffari   ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble",
6488626eedSJames Wright                             NULL, rc, &rc, NULL);
6588626eedSJames Wright   CHKERRQ(ierr);
66ba6664aeSJames Wright   for (PetscInt i = 0; i < 3; i++)
6788626eedSJames Wright     center[i] = .5 * domain_size[i];
6877841947SLeila Ghaffari   PetscInt n = problem->dim;
6988626eedSJames Wright   ierr = PetscOptionsRealArray("-center", "Location of bubble center", NULL,
7088626eedSJames Wright                                center, &n, NULL);
7188626eedSJames Wright   CHKERRQ(ierr);
7277841947SLeila Ghaffari   n = problem->dim;
7377841947SLeila Ghaffari   ierr = PetscOptionsRealArray("-dc_axis",
7488626eedSJames Wright                                "Axis of density current cylindrical anomaly, "
7588626eedSJames Wright                                "or {0,0,0} for spherically symmetric",
7688626eedSJames Wright                                NULL, dc_axis, &n, NULL);
7788626eedSJames Wright   CHKERRQ(ierr);
7877841947SLeila Ghaffari   {
7977841947SLeila Ghaffari     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) + PetscSqr(dc_axis[1]) +
8077841947SLeila Ghaffari                                    PetscSqr(dc_axis[2]));
8177841947SLeila Ghaffari     if (norm > 0) {
82ba6664aeSJames Wright       for (PetscInt i = 0; i < 3; i++)
8388626eedSJames Wright         dc_axis[i] /= norm;
8477841947SLeila Ghaffari     }
8577841947SLeila Ghaffari   }
8677841947SLeila Ghaffari 
8767490bc6SJeremy L Thompson   PetscOptionsEnd();
8877841947SLeila Ghaffari 
8988b783a1SJames Wright   PetscScalar meter  = user->units->meter;
9088626eedSJames Wright   PetscScalar second = user->units->second;
9188626eedSJames Wright   PetscScalar Kelvin = user->units->Kelvin;
9288626eedSJames Wright   PetscScalar Pascal = user->units->Pascal;
9377841947SLeila Ghaffari   rc = fabs(rc) * meter;
9488626eedSJames Wright   theta0 *= Kelvin;
9588626eedSJames Wright   thetaC *= Kelvin;
9688626eedSJames Wright   P0 *= Pascal;
9788626eedSJames Wright   N *= (1. / second);
98ba6664aeSJames Wright   for (PetscInt i = 0; i < 3; i++)
9988626eedSJames Wright     center[i] *= meter;
10077841947SLeila Ghaffari 
101*dc805cc4SLeila Ghaffari   dc_ctx->theta0 = theta0;
102*dc805cc4SLeila Ghaffari   dc_ctx->thetaC = thetaC;
103*dc805cc4SLeila Ghaffari   dc_ctx->P0 = P0;
104*dc805cc4SLeila Ghaffari   dc_ctx->N = N;
105*dc805cc4SLeila Ghaffari   dc_ctx->rc = rc;
106*dc805cc4SLeila Ghaffari   dc_ctx->center[0] = center[0];
107*dc805cc4SLeila Ghaffari   dc_ctx->center[1] = center[1];
108*dc805cc4SLeila Ghaffari   dc_ctx->center[2] = center[2];
109*dc805cc4SLeila Ghaffari   dc_ctx->dc_axis[0] = dc_axis[0];
110*dc805cc4SLeila Ghaffari   dc_ctx->dc_axis[1] = dc_axis[1];
111*dc805cc4SLeila Ghaffari   dc_ctx->dc_axis[2] = dc_axis[2];
11277841947SLeila Ghaffari 
113*dc805cc4SLeila Ghaffari   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context,
114*dc805cc4SLeila Ghaffari                               CEED_MEM_HOST, &newtonian_ig_ctx);
115*dc805cc4SLeila Ghaffari   dc_ctx->newtonian_ctx = *newtonian_ig_ctx;
116*dc805cc4SLeila Ghaffari   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context,
117*dc805cc4SLeila Ghaffari                                   &newtonian_ig_ctx);
118*dc805cc4SLeila Ghaffari   CeedQFunctionContextCreate(user->ceed, &density_current_context);
119*dc805cc4SLeila Ghaffari   CeedQFunctionContextSetData(density_current_context, CEED_MEM_HOST,
120*dc805cc4SLeila Ghaffari                               CEED_USE_POINTER, sizeof(*dc_ctx), dc_ctx);
121*dc805cc4SLeila Ghaffari   CeedQFunctionContextSetDataDestroy(density_current_context, CEED_MEM_HOST,
122*dc805cc4SLeila Ghaffari                                      FreeContextPetsc);
123*dc805cc4SLeila Ghaffari   problem->ics.qfunction_context = density_current_context;
12417ce10faSJeremy L Thompson 
12577841947SLeila Ghaffari   PetscFunctionReturn(0);
12677841947SLeila Ghaffari }
127