1bb8a0c61SJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other 2bb8a0c61SJames Wright // CEED contributors. All Rights Reserved. See the top-level LICENSE and NOTICE 3bb8a0c61SJames Wright // files for details. 4a515125bSLeila Ghaffari // 5727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 6a515125bSLeila Ghaffari // 7727da7e7SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 8a515125bSLeila Ghaffari 9a515125bSLeila Ghaffari /// @file 10a515125bSLeila Ghaffari /// Utility functions for setting up DENSITY_CURRENT 11a515125bSLeila Ghaffari 12a515125bSLeila Ghaffari #include "../qfunctions/densitycurrent.h" 13bb8a0c61SJames Wright #include "../navierstokes.h" 14a515125bSLeila Ghaffari 15b7f03f12SJed Brown PetscErrorCode NS_DENSITY_CURRENT(ProblemData *problem, DM dm, void *ctx) { 163a8779fbSJames Wright 173a8779fbSJames Wright PetscInt ierr; 18a515125bSLeila Ghaffari MPI_Comm comm = PETSC_COMM_WORLD; 19*cbe60e31SLeila Ghaffari User user = *(User *)ctx; 20*cbe60e31SLeila Ghaffari DensityCurrentContext dc_ctx; 21*cbe60e31SLeila Ghaffari CeedQFunctionContext density_current_context; 22*cbe60e31SLeila Ghaffari NewtonianIdealGasContext newtonian_ig_ctx; 23a515125bSLeila Ghaffari 24b7f03f12SJed Brown PetscFunctionBeginUser; 25b7f03f12SJed Brown ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr); 26*cbe60e31SLeila Ghaffari ierr = PetscCalloc1(1, &dc_ctx); CHKERRQ(ierr); 27a515125bSLeila Ghaffari // ------------------------------------------------------ 28a515125bSLeila Ghaffari // SET UP DENSITY_CURRENT 29a515125bSLeila Ghaffari // ------------------------------------------------------ 30*cbe60e31SLeila Ghaffari CeedQFunctionContextDestroy(&problem->ics.qfunction_context); 319785fe93SJed Brown problem->ics.qfunction = ICsDC; 329785fe93SJed Brown problem->ics.qfunction_loc = ICsDC_loc; 33a515125bSLeila Ghaffari 34a515125bSLeila Ghaffari // ------------------------------------------------------ 35a515125bSLeila Ghaffari // Create the libCEED context 36a515125bSLeila Ghaffari // ------------------------------------------------------ 37bb8a0c61SJames Wright CeedScalar theta0 = 300.; // K 38bb8a0c61SJames Wright CeedScalar thetaC = -15.; // K 39bb8a0c61SJames Wright CeedScalar P0 = 1.e5; // Pa 40bb8a0c61SJames Wright CeedScalar N = 0.01; // 1/s 41a515125bSLeila Ghaffari CeedScalar rc = 1000.; // m (Radius of bubble) 42a515125bSLeila Ghaffari PetscReal center[3], dc_axis[3] = {0, 0, 0}; 4305a512bdSLeila Ghaffari PetscReal domain_min[3], domain_max[3], domain_size[3]; 44bb8a0c61SJames Wright ierr = DMGetBoundingBox(dm, domain_min, domain_max); 45bb8a0c61SJames Wright CHKERRQ(ierr); 46493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) 47bb8a0c61SJames Wright 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); 53bb8a0c61SJames Wright ierr = PetscOptionsScalar("-theta0", "Reference potential temperature", NULL, 54bb8a0c61SJames Wright theta0, &theta0, NULL); 55bb8a0c61SJames Wright CHKERRQ(ierr); 56bb8a0c61SJames Wright ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature", 57bb8a0c61SJames Wright NULL, thetaC, &thetaC, NULL); 58bb8a0c61SJames Wright CHKERRQ(ierr); 59bb8a0c61SJames Wright ierr = PetscOptionsScalar("-P0", "Atmospheric pressure", NULL, P0, &P0, NULL); 60bb8a0c61SJames Wright CHKERRQ(ierr); 61bb8a0c61SJames Wright ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency", NULL, N, &N, NULL); 62bb8a0c61SJames Wright CHKERRQ(ierr); 63a515125bSLeila Ghaffari ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", 64bb8a0c61SJames Wright NULL, rc, &rc, NULL); 65bb8a0c61SJames Wright CHKERRQ(ierr); 66493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) 67bb8a0c61SJames Wright center[i] = .5 * domain_size[i]; 68a515125bSLeila Ghaffari PetscInt n = problem->dim; 69bb8a0c61SJames Wright ierr = PetscOptionsRealArray("-center", "Location of bubble center", NULL, 70bb8a0c61SJames Wright center, &n, NULL); 71bb8a0c61SJames Wright CHKERRQ(ierr); 72a515125bSLeila Ghaffari n = problem->dim; 73a515125bSLeila Ghaffari ierr = PetscOptionsRealArray("-dc_axis", 74bb8a0c61SJames Wright "Axis of density current cylindrical anomaly, " 75bb8a0c61SJames Wright "or {0,0,0} for spherically symmetric", 76bb8a0c61SJames Wright NULL, dc_axis, &n, NULL); 77bb8a0c61SJames Wright CHKERRQ(ierr); 78a515125bSLeila Ghaffari { 79a515125bSLeila Ghaffari PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) + PetscSqr(dc_axis[1]) + 80a515125bSLeila Ghaffari PetscSqr(dc_axis[2])); 81a515125bSLeila Ghaffari if (norm > 0) { 82493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) 83bb8a0c61SJames Wright dc_axis[i] /= norm; 84a515125bSLeila Ghaffari } 85a515125bSLeila Ghaffari } 86a515125bSLeila Ghaffari 871485969bSJeremy L Thompson PetscOptionsEnd(); 88a515125bSLeila Ghaffari 893a8779fbSJames Wright PetscScalar meter = user->units->meter; 90bb8a0c61SJames Wright PetscScalar second = user->units->second; 91bb8a0c61SJames Wright PetscScalar Kelvin = user->units->Kelvin; 92bb8a0c61SJames Wright PetscScalar Pascal = user->units->Pascal; 93a515125bSLeila Ghaffari rc = fabs(rc) * meter; 94bb8a0c61SJames Wright theta0 *= Kelvin; 95bb8a0c61SJames Wright thetaC *= Kelvin; 96bb8a0c61SJames Wright P0 *= Pascal; 97bb8a0c61SJames Wright N *= (1. / second); 98493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) 99bb8a0c61SJames Wright center[i] *= meter; 100a515125bSLeila Ghaffari 101*cbe60e31SLeila Ghaffari dc_ctx->theta0 = theta0; 102*cbe60e31SLeila Ghaffari dc_ctx->thetaC = thetaC; 103*cbe60e31SLeila Ghaffari dc_ctx->P0 = P0; 104*cbe60e31SLeila Ghaffari dc_ctx->N = N; 105*cbe60e31SLeila Ghaffari dc_ctx->rc = rc; 106*cbe60e31SLeila Ghaffari dc_ctx->center[0] = center[0]; 107*cbe60e31SLeila Ghaffari dc_ctx->center[1] = center[1]; 108*cbe60e31SLeila Ghaffari dc_ctx->center[2] = center[2]; 109*cbe60e31SLeila Ghaffari dc_ctx->dc_axis[0] = dc_axis[0]; 110*cbe60e31SLeila Ghaffari dc_ctx->dc_axis[1] = dc_axis[1]; 111*cbe60e31SLeila Ghaffari dc_ctx->dc_axis[2] = dc_axis[2]; 112a515125bSLeila Ghaffari 113*cbe60e31SLeila Ghaffari CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, 114*cbe60e31SLeila Ghaffari CEED_MEM_HOST, &newtonian_ig_ctx); 115*cbe60e31SLeila Ghaffari dc_ctx->newtonian_ctx = *newtonian_ig_ctx; 116*cbe60e31SLeila Ghaffari CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, 117*cbe60e31SLeila Ghaffari &newtonian_ig_ctx); 118*cbe60e31SLeila Ghaffari CeedQFunctionContextCreate(user->ceed, &density_current_context); 119*cbe60e31SLeila Ghaffari CeedQFunctionContextSetData(density_current_context, CEED_MEM_HOST, 120*cbe60e31SLeila Ghaffari CEED_USE_POINTER, sizeof(*dc_ctx), dc_ctx); 121*cbe60e31SLeila Ghaffari CeedQFunctionContextSetDataDestroy(density_current_context, CEED_MEM_HOST, 122*cbe60e31SLeila Ghaffari FreeContextPetsc); 123*cbe60e31SLeila Ghaffari problem->ics.qfunction_context = density_current_context; 124270bbb13SJeremy L Thompson 125a515125bSLeila Ghaffari PetscFunctionReturn(0); 126a515125bSLeila Ghaffari } 127