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) { 15*0c373b74SJames Wright Honee honee = *(Honee *)ctx; 16*0c373b74SJames Wright MPI_Comm comm = honee->comm; 17*0c373b74SJames Wright Ceed ceed = honee->ceed; 18e07531f7SJames Wright DensityCurrentContext density_current_ctx; 19e07531f7SJames Wright CeedQFunctionContext density_current_qfctx; 20cbe60e31SLeila Ghaffari NewtonianIdealGasContext newtonian_ig_ctx; 2166d54740SJames Wright PetscInt dim; 22a515125bSLeila Ghaffari 23b7f03f12SJed Brown PetscFunctionBeginUser; 24d1c51a42SJames Wright PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc)); 25e07531f7SJames Wright PetscCall(PetscCalloc1(1, &density_current_ctx)); 26a515125bSLeila Ghaffari // ------------------------------------------------------ 27a515125bSLeila Ghaffari // SET UP DENSITY_CURRENT 28a515125bSLeila Ghaffari // ------------------------------------------------------ 29e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfctx)); 30e07531f7SJames Wright problem->ics.qf_func_ptr = ICsDC; 31e07531f7SJames Wright problem->ics.qf_loc = ICsDC_loc; 32a515125bSLeila Ghaffari 33a515125bSLeila Ghaffari // ------------------------------------------------------ 34a515125bSLeila Ghaffari // Create the libCEED context 35a515125bSLeila Ghaffari // ------------------------------------------------------ 36bb8a0c61SJames Wright CeedScalar theta0 = 300.; // K 37bb8a0c61SJames Wright CeedScalar thetaC = -15.; // K 38bb8a0c61SJames Wright CeedScalar P0 = 1.e5; // Pa 39bb8a0c61SJames Wright CeedScalar N = 0.01; // 1/s 40a515125bSLeila Ghaffari CeedScalar rc = 1000.; // m (Radius of bubble) 41a515125bSLeila Ghaffari PetscReal center[3], dc_axis[3] = {0, 0, 0}; 4205a512bdSLeila Ghaffari PetscReal domain_min[3], domain_max[3], domain_size[3]; 432b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 4466d54740SJames Wright PetscCall(DMGetDimension(dm, &dim)); 4566d54740SJames Wright for (PetscInt i = 0; i < dim; i++) domain_size[i] = domain_max[i] - domain_min[i]; 46a515125bSLeila Ghaffari 47a515125bSLeila Ghaffari // ------------------------------------------------------ 48a515125bSLeila Ghaffari // Command line Options 49a515125bSLeila Ghaffari // ------------------------------------------------------ 501485969bSJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for DENSITY_CURRENT problem", NULL); 512b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-theta0", "Reference potential temperature", NULL, theta0, &theta0, NULL)); 522b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-thetaC", "Perturbation of potential temperature", NULL, thetaC, &thetaC, NULL)); 532b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-P0", "Atmospheric pressure", NULL, P0, &P0, NULL)); 542b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-N", "Brunt-Vaisala frequency", NULL, N, &N, NULL)); 552b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL)); 5666d54740SJames Wright for (PetscInt i = 0; i < dim; i++) center[i] = .5 * domain_size[i]; 5766d54740SJames Wright PetscInt n = dim; 582b916ea7SJeremy L Thompson PetscCall(PetscOptionsRealArray("-center", "Location of bubble center", NULL, center, &n, NULL)); 5966d54740SJames Wright n = dim; 602b916ea7SJeremy L Thompson PetscCall(PetscOptionsRealArray("-dc_axis", 61bb8a0c61SJames Wright "Axis of density current cylindrical anomaly, " 62bb8a0c61SJames Wright "or {0,0,0} for spherically symmetric", 632b916ea7SJeremy L Thompson NULL, dc_axis, &n, NULL)); 64a515125bSLeila Ghaffari { 652b916ea7SJeremy L Thompson PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) + PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2])); 66a515125bSLeila Ghaffari if (norm > 0) { 6766d54740SJames Wright for (PetscInt i = 0; i < dim; i++) dc_axis[i] /= norm; 68a515125bSLeila Ghaffari } 69a515125bSLeila Ghaffari } 70a515125bSLeila Ghaffari 711485969bSJeremy L Thompson PetscOptionsEnd(); 72a515125bSLeila Ghaffari 73*0c373b74SJames Wright PetscScalar meter = honee->units->meter; 74*0c373b74SJames Wright PetscScalar second = honee->units->second; 75*0c373b74SJames Wright PetscScalar Kelvin = honee->units->Kelvin; 76*0c373b74SJames Wright PetscScalar Pascal = honee->units->Pascal; 77a515125bSLeila Ghaffari rc = fabs(rc) * meter; 78bb8a0c61SJames Wright theta0 *= Kelvin; 79bb8a0c61SJames Wright thetaC *= Kelvin; 80bb8a0c61SJames Wright P0 *= Pascal; 81bb8a0c61SJames Wright N *= (1. / second); 8266d54740SJames Wright for (PetscInt i = 0; i < dim; i++) center[i] *= meter; 83a515125bSLeila Ghaffari 84e07531f7SJames Wright density_current_ctx->theta0 = theta0; 85e07531f7SJames Wright density_current_ctx->thetaC = thetaC; 86e07531f7SJames Wright density_current_ctx->P0 = P0; 87e07531f7SJames Wright density_current_ctx->N = N; 88e07531f7SJames Wright density_current_ctx->rc = rc; 8966d54740SJames Wright PetscCall(PetscArraycpy(density_current_ctx->center, center, 3)); 9066d54740SJames Wright PetscCall(PetscArraycpy(density_current_ctx->dc_axis, dc_axis, 3)); 91a515125bSLeila Ghaffari 92e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ig_ctx)); 93e07531f7SJames Wright density_current_ctx->newtonian_ctx = *newtonian_ig_ctx; 94e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ig_ctx)); 95*0c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &density_current_qfctx)); 96e07531f7SJames Wright PetscCallCeed( 97e07531f7SJames Wright ceed, CeedQFunctionContextSetData(density_current_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*density_current_ctx), density_current_ctx)); 98e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(density_current_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 99e07531f7SJames Wright problem->ics.qfctx = density_current_qfctx; 100d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 101a515125bSLeila Ghaffari } 102