1dc936754SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 204e40bb6SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3a515125bSLeila Ghaffari // 4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 5a515125bSLeila Ghaffari // 6727da7e7SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 7a515125bSLeila Ghaffari 8a515125bSLeila Ghaffari /// @file 9a515125bSLeila Ghaffari /// Utility functions for setting up DENSITY_CURRENT 10a515125bSLeila Ghaffari 11a515125bSLeila Ghaffari #include "../qfunctions/densitycurrent.h" 122b916ea7SJeremy L Thompson 13e419654dSJeremy L Thompson #include <ceed.h> 14e419654dSJeremy L Thompson #include <petscdm.h> 15e419654dSJeremy L Thompson 16*149fb536SJames Wright #include <navierstokes.h> 17a515125bSLeila Ghaffari 18991aef52SJames Wright PetscErrorCode NS_DENSITY_CURRENT(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 19cbe60e31SLeila Ghaffari User user = *(User *)ctx; 20b4c37c5cSJames Wright MPI_Comm comm = user->comm; 21b4c37c5cSJames Wright Ceed ceed = user->ceed; 22cbe60e31SLeila Ghaffari DensityCurrentContext dc_ctx; 23cbe60e31SLeila Ghaffari CeedQFunctionContext density_current_context; 24cbe60e31SLeila Ghaffari NewtonianIdealGasContext newtonian_ig_ctx; 25a515125bSLeila Ghaffari 26b7f03f12SJed Brown PetscFunctionBeginUser; 27d1c51a42SJames Wright PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc)); 282b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &dc_ctx)); 29a515125bSLeila Ghaffari // ------------------------------------------------------ 30a515125bSLeila Ghaffari // SET UP DENSITY_CURRENT 31a515125bSLeila Ghaffari // ------------------------------------------------------ 32b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfunction_context)); 339785fe93SJed Brown problem->ics.qfunction = ICsDC; 349785fe93SJed Brown problem->ics.qfunction_loc = ICsDC_loc; 35a515125bSLeila Ghaffari 36a515125bSLeila Ghaffari // ------------------------------------------------------ 37a515125bSLeila Ghaffari // Create the libCEED context 38a515125bSLeila Ghaffari // ------------------------------------------------------ 39bb8a0c61SJames Wright CeedScalar theta0 = 300.; // K 40bb8a0c61SJames Wright CeedScalar thetaC = -15.; // K 41bb8a0c61SJames Wright CeedScalar P0 = 1.e5; // Pa 42bb8a0c61SJames Wright CeedScalar N = 0.01; // 1/s 43a515125bSLeila Ghaffari CeedScalar rc = 1000.; // m (Radius of bubble) 44a515125bSLeila Ghaffari PetscReal center[3], dc_axis[3] = {0, 0, 0}; 4505a512bdSLeila Ghaffari PetscReal domain_min[3], domain_max[3], domain_size[3]; 462b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 472b916ea7SJeremy L Thompson for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 48a515125bSLeila Ghaffari 49a515125bSLeila Ghaffari // ------------------------------------------------------ 50a515125bSLeila Ghaffari // Command line Options 51a515125bSLeila Ghaffari // ------------------------------------------------------ 521485969bSJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for DENSITY_CURRENT problem", NULL); 532b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-theta0", "Reference potential temperature", NULL, theta0, &theta0, NULL)); 542b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-thetaC", "Perturbation of potential temperature", NULL, thetaC, &thetaC, NULL)); 552b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-P0", "Atmospheric pressure", NULL, P0, &P0, NULL)); 562b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-N", "Brunt-Vaisala frequency", NULL, N, &N, NULL)); 572b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL)); 582b916ea7SJeremy L Thompson for (PetscInt i = 0; i < 3; i++) center[i] = .5 * domain_size[i]; 59a515125bSLeila Ghaffari PetscInt n = problem->dim; 602b916ea7SJeremy L Thompson PetscCall(PetscOptionsRealArray("-center", "Location of bubble center", NULL, center, &n, NULL)); 61a515125bSLeila Ghaffari n = problem->dim; 622b916ea7SJeremy L Thompson PetscCall(PetscOptionsRealArray("-dc_axis", 63bb8a0c61SJames Wright "Axis of density current cylindrical anomaly, " 64bb8a0c61SJames Wright "or {0,0,0} for spherically symmetric", 652b916ea7SJeremy L Thompson NULL, dc_axis, &n, NULL)); 66a515125bSLeila Ghaffari { 672b916ea7SJeremy L Thompson PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) + PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2])); 68a515125bSLeila Ghaffari if (norm > 0) { 692b916ea7SJeremy L Thompson for (PetscInt i = 0; i < 3; i++) dc_axis[i] /= norm; 70a515125bSLeila Ghaffari } 71a515125bSLeila Ghaffari } 72a515125bSLeila Ghaffari 731485969bSJeremy L Thompson PetscOptionsEnd(); 74a515125bSLeila Ghaffari 753a8779fbSJames Wright PetscScalar meter = user->units->meter; 76bb8a0c61SJames Wright PetscScalar second = user->units->second; 77bb8a0c61SJames Wright PetscScalar Kelvin = user->units->Kelvin; 78bb8a0c61SJames Wright PetscScalar Pascal = user->units->Pascal; 79a515125bSLeila Ghaffari rc = fabs(rc) * meter; 80bb8a0c61SJames Wright theta0 *= Kelvin; 81bb8a0c61SJames Wright thetaC *= Kelvin; 82bb8a0c61SJames Wright P0 *= Pascal; 83bb8a0c61SJames Wright N *= (1. / second); 842b916ea7SJeremy L Thompson for (PetscInt i = 0; i < 3; i++) center[i] *= meter; 85a515125bSLeila Ghaffari 86cbe60e31SLeila Ghaffari dc_ctx->theta0 = theta0; 87cbe60e31SLeila Ghaffari dc_ctx->thetaC = thetaC; 88cbe60e31SLeila Ghaffari dc_ctx->P0 = P0; 89cbe60e31SLeila Ghaffari dc_ctx->N = N; 90cbe60e31SLeila Ghaffari dc_ctx->rc = rc; 91cbe60e31SLeila Ghaffari dc_ctx->center[0] = center[0]; 92cbe60e31SLeila Ghaffari dc_ctx->center[1] = center[1]; 93cbe60e31SLeila Ghaffari dc_ctx->center[2] = center[2]; 94cbe60e31SLeila Ghaffari dc_ctx->dc_axis[0] = dc_axis[0]; 95cbe60e31SLeila Ghaffari dc_ctx->dc_axis[1] = dc_axis[1]; 96cbe60e31SLeila Ghaffari dc_ctx->dc_axis[2] = dc_axis[2]; 97a515125bSLeila Ghaffari 98b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx)); 99cbe60e31SLeila Ghaffari dc_ctx->newtonian_ctx = *newtonian_ig_ctx; 100b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx)); 101b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &density_current_context)); 102b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(density_current_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*dc_ctx), dc_ctx)); 103b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(density_current_context, CEED_MEM_HOST, FreeContextPetsc)); 104cbe60e31SLeila Ghaffari problem->ics.qfunction_context = density_current_context; 105d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 106a515125bSLeila Ghaffari } 107