xref: /honee/problems/densitycurrent.c (revision e07531f7cad484157c281bb3cb12a29a67a96955)
1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3a515125bSLeila Ghaffari 
4a515125bSLeila Ghaffari /// @file
5a515125bSLeila Ghaffari /// Utility functions for setting up DENSITY_CURRENT
6a515125bSLeila Ghaffari 
7a515125bSLeila Ghaffari #include "../qfunctions/densitycurrent.h"
82b916ea7SJeremy L Thompson 
9e419654dSJeremy L Thompson #include <ceed.h>
10e419654dSJeremy L Thompson #include <petscdm.h>
11e419654dSJeremy L Thompson 
12149fb536SJames Wright #include <navierstokes.h>
13a515125bSLeila Ghaffari 
14991aef52SJames Wright PetscErrorCode NS_DENSITY_CURRENT(ProblemData problem, DM dm, void *ctx, SimpleBC bc) {
15cbe60e31SLeila Ghaffari   User                     user = *(User *)ctx;
16b4c37c5cSJames Wright   MPI_Comm                 comm = user->comm;
17b4c37c5cSJames Wright   Ceed                     ceed = user->ceed;
18*e07531f7SJames Wright   DensityCurrentContext    density_current_ctx;
19*e07531f7SJames Wright   CeedQFunctionContext     density_current_qfctx;
20cbe60e31SLeila Ghaffari   NewtonianIdealGasContext newtonian_ig_ctx;
21a515125bSLeila Ghaffari 
22b7f03f12SJed Brown   PetscFunctionBeginUser;
23d1c51a42SJames Wright   PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc));
24*e07531f7SJames Wright   PetscCall(PetscCalloc1(1, &density_current_ctx));
25a515125bSLeila Ghaffari   // ------------------------------------------------------
26a515125bSLeila Ghaffari   //               SET UP DENSITY_CURRENT
27a515125bSLeila Ghaffari   // ------------------------------------------------------
28*e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfctx));
29*e07531f7SJames Wright   problem->ics.qf_func_ptr = ICsDC;
30*e07531f7SJames Wright   problem->ics.qf_loc      = ICsDC_loc;
31a515125bSLeila Ghaffari 
32a515125bSLeila Ghaffari   // ------------------------------------------------------
33a515125bSLeila Ghaffari   //             Create the libCEED context
34a515125bSLeila Ghaffari   // ------------------------------------------------------
35bb8a0c61SJames Wright   CeedScalar theta0 = 300.;   // K
36bb8a0c61SJames Wright   CeedScalar thetaC = -15.;   // K
37bb8a0c61SJames Wright   CeedScalar P0     = 1.e5;   // Pa
38bb8a0c61SJames Wright   CeedScalar N      = 0.01;   // 1/s
39a515125bSLeila Ghaffari   CeedScalar rc     = 1000.;  // m (Radius of bubble)
40a515125bSLeila Ghaffari   PetscReal  center[3], dc_axis[3] = {0, 0, 0};
4105a512bdSLeila Ghaffari   PetscReal  domain_min[3], domain_max[3], domain_size[3];
422b916ea7SJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
432b916ea7SJeremy L Thompson   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
44a515125bSLeila Ghaffari 
45a515125bSLeila Ghaffari   // ------------------------------------------------------
46a515125bSLeila Ghaffari   //              Command line Options
47a515125bSLeila Ghaffari   // ------------------------------------------------------
481485969bSJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for DENSITY_CURRENT problem", NULL);
492b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-theta0", "Reference potential temperature", NULL, theta0, &theta0, NULL));
502b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-thetaC", "Perturbation of potential temperature", NULL, thetaC, &thetaC, NULL));
512b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-P0", "Atmospheric pressure", NULL, P0, &P0, NULL));
522b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-N", "Brunt-Vaisala frequency", NULL, N, &N, NULL));
532b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL));
542b916ea7SJeremy L Thompson   for (PetscInt i = 0; i < 3; i++) center[i] = .5 * domain_size[i];
55a515125bSLeila Ghaffari   PetscInt n = problem->dim;
562b916ea7SJeremy L Thompson   PetscCall(PetscOptionsRealArray("-center", "Location of bubble center", NULL, center, &n, NULL));
57a515125bSLeila Ghaffari   n = problem->dim;
582b916ea7SJeremy L Thompson   PetscCall(PetscOptionsRealArray("-dc_axis",
59bb8a0c61SJames Wright                                   "Axis of density current cylindrical anomaly, "
60bb8a0c61SJames Wright                                   "or {0,0,0} for spherically symmetric",
612b916ea7SJeremy L Thompson                                   NULL, dc_axis, &n, NULL));
62a515125bSLeila Ghaffari   {
632b916ea7SJeremy L Thompson     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) + PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2]));
64a515125bSLeila Ghaffari     if (norm > 0) {
652b916ea7SJeremy L Thompson       for (PetscInt i = 0; i < 3; i++) dc_axis[i] /= norm;
66a515125bSLeila Ghaffari     }
67a515125bSLeila Ghaffari   }
68a515125bSLeila Ghaffari 
691485969bSJeremy L Thompson   PetscOptionsEnd();
70a515125bSLeila Ghaffari 
713a8779fbSJames Wright   PetscScalar meter  = user->units->meter;
72bb8a0c61SJames Wright   PetscScalar second = user->units->second;
73bb8a0c61SJames Wright   PetscScalar Kelvin = user->units->Kelvin;
74bb8a0c61SJames Wright   PetscScalar Pascal = user->units->Pascal;
75a515125bSLeila Ghaffari   rc                 = fabs(rc) * meter;
76bb8a0c61SJames Wright   theta0 *= Kelvin;
77bb8a0c61SJames Wright   thetaC *= Kelvin;
78bb8a0c61SJames Wright   P0 *= Pascal;
79bb8a0c61SJames Wright   N *= (1. / second);
802b916ea7SJeremy L Thompson   for (PetscInt i = 0; i < 3; i++) center[i] *= meter;
81a515125bSLeila Ghaffari 
82*e07531f7SJames Wright   density_current_ctx->theta0     = theta0;
83*e07531f7SJames Wright   density_current_ctx->thetaC     = thetaC;
84*e07531f7SJames Wright   density_current_ctx->P0         = P0;
85*e07531f7SJames Wright   density_current_ctx->N          = N;
86*e07531f7SJames Wright   density_current_ctx->rc         = rc;
87*e07531f7SJames Wright   density_current_ctx->center[0]  = center[0];
88*e07531f7SJames Wright   density_current_ctx->center[1]  = center[1];
89*e07531f7SJames Wright   density_current_ctx->center[2]  = center[2];
90*e07531f7SJames Wright   density_current_ctx->dc_axis[0] = dc_axis[0];
91*e07531f7SJames Wright   density_current_ctx->dc_axis[1] = dc_axis[1];
92*e07531f7SJames Wright   density_current_ctx->dc_axis[2] = dc_axis[2];
93a515125bSLeila Ghaffari 
94*e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ig_ctx));
95*e07531f7SJames Wright   density_current_ctx->newtonian_ctx = *newtonian_ig_ctx;
96*e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ig_ctx));
97*e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &density_current_qfctx));
98*e07531f7SJames Wright   PetscCallCeed(
99*e07531f7SJames Wright       ceed, CeedQFunctionContextSetData(density_current_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*density_current_ctx), density_current_ctx));
100*e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(density_current_qfctx, CEED_MEM_HOST, FreeContextPetsc));
101*e07531f7SJames Wright   problem->ics.qfctx = density_current_qfctx;
102d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
103a515125bSLeila Ghaffari }
104