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; 18b7f03f12SJed Brown SetupContext setup_context; 19a515125bSLeila Ghaffari User user = *(User *)ctx; 20a515125bSLeila Ghaffari MPI_Comm comm = PETSC_COMM_WORLD; 21a515125bSLeila Ghaffari 22b7f03f12SJed Brown PetscFunctionBeginUser; 23b7f03f12SJed Brown ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr); 24a515125bSLeila Ghaffari // ------------------------------------------------------ 25a515125bSLeila Ghaffari // SET UP DENSITY_CURRENT 26a515125bSLeila Ghaffari // ------------------------------------------------------ 279785fe93SJed Brown problem->ics.qfunction = ICsDC; 289785fe93SJed Brown problem->ics.qfunction_loc = ICsDC_loc; 29a515125bSLeila Ghaffari problem->bc = Exact_DC; 30b7f03f12SJed Brown setup_context = problem->bc_ctx; 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]; 42bb8a0c61SJames Wright ierr = DMGetBoundingBox(dm, domain_min, domain_max); 43bb8a0c61SJames Wright CHKERRQ(ierr); 44*493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) 45bb8a0c61SJames Wright 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); 51bb8a0c61SJames Wright ierr = PetscOptionsScalar("-theta0", "Reference potential temperature", NULL, 52bb8a0c61SJames Wright theta0, &theta0, NULL); 53bb8a0c61SJames Wright CHKERRQ(ierr); 54bb8a0c61SJames Wright ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature", 55bb8a0c61SJames Wright NULL, thetaC, &thetaC, NULL); 56bb8a0c61SJames Wright CHKERRQ(ierr); 57bb8a0c61SJames Wright ierr = PetscOptionsScalar("-P0", "Atmospheric pressure", NULL, P0, &P0, NULL); 58bb8a0c61SJames Wright CHKERRQ(ierr); 59bb8a0c61SJames Wright ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency", NULL, N, &N, NULL); 60bb8a0c61SJames Wright CHKERRQ(ierr); 61a515125bSLeila Ghaffari ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", 62bb8a0c61SJames Wright NULL, rc, &rc, NULL); 63bb8a0c61SJames Wright CHKERRQ(ierr); 64*493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) 65bb8a0c61SJames Wright center[i] = .5 * domain_size[i]; 66a515125bSLeila Ghaffari PetscInt n = problem->dim; 67bb8a0c61SJames Wright ierr = PetscOptionsRealArray("-center", "Location of bubble center", NULL, 68bb8a0c61SJames Wright center, &n, NULL); 69bb8a0c61SJames Wright CHKERRQ(ierr); 70a515125bSLeila Ghaffari n = problem->dim; 71a515125bSLeila Ghaffari ierr = PetscOptionsRealArray("-dc_axis", 72bb8a0c61SJames Wright "Axis of density current cylindrical anomaly, " 73bb8a0c61SJames Wright "or {0,0,0} for spherically symmetric", 74bb8a0c61SJames Wright NULL, dc_axis, &n, NULL); 75bb8a0c61SJames Wright CHKERRQ(ierr); 76a515125bSLeila Ghaffari { 77a515125bSLeila Ghaffari PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) + PetscSqr(dc_axis[1]) + 78a515125bSLeila Ghaffari PetscSqr(dc_axis[2])); 79a515125bSLeila Ghaffari if (norm > 0) { 80*493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) 81bb8a0c61SJames Wright dc_axis[i] /= norm; 82a515125bSLeila Ghaffari } 83a515125bSLeila Ghaffari } 84a515125bSLeila Ghaffari 851485969bSJeremy L Thompson PetscOptionsEnd(); 86a515125bSLeila Ghaffari 873a8779fbSJames Wright PetscScalar meter = user->units->meter; 88bb8a0c61SJames Wright PetscScalar second = user->units->second; 89bb8a0c61SJames Wright PetscScalar Kelvin = user->units->Kelvin; 90bb8a0c61SJames Wright PetscScalar Pascal = user->units->Pascal; 91a515125bSLeila Ghaffari rc = fabs(rc) * meter; 92bb8a0c61SJames Wright theta0 *= Kelvin; 93bb8a0c61SJames Wright thetaC *= Kelvin; 94bb8a0c61SJames Wright P0 *= Pascal; 95bb8a0c61SJames Wright N *= (1. / second); 96*493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) 97bb8a0c61SJames Wright center[i] *= meter; 98a515125bSLeila Ghaffari 99bb8a0c61SJames Wright setup_context->theta0 = theta0; 100bb8a0c61SJames Wright setup_context->thetaC = thetaC; 101bb8a0c61SJames Wright setup_context->P0 = P0; 102bb8a0c61SJames Wright setup_context->N = N; 103a515125bSLeila Ghaffari setup_context->rc = rc; 104a515125bSLeila Ghaffari setup_context->center[0] = center[0]; 105a515125bSLeila Ghaffari setup_context->center[1] = center[1]; 106a515125bSLeila Ghaffari setup_context->center[2] = center[2]; 107a515125bSLeila Ghaffari setup_context->dc_axis[0] = dc_axis[0]; 108a515125bSLeila Ghaffari setup_context->dc_axis[1] = dc_axis[1]; 109a515125bSLeila Ghaffari setup_context->dc_axis[2] = dc_axis[2]; 110a515125bSLeila Ghaffari 111a515125bSLeila Ghaffari PetscFunctionReturn(0); 112a515125bSLeila Ghaffari } 113