1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other 2 // CEED contributors. All Rights Reserved. See the top-level LICENSE and NOTICE 3 // files for details. 4 // 5 // SPDX-License-Identifier: BSD-2-Clause 6 // 7 // This file is part of CEED: http://github.com/ceed 8 9 /// @file 10 /// Utility functions for setting up DENSITY_CURRENT 11 12 #include "../qfunctions/densitycurrent.h" 13 #include "../navierstokes.h" 14 15 PetscErrorCode NS_DENSITY_CURRENT(ProblemData *problem, DM dm, void *ctx) { 16 17 PetscInt ierr; 18 SetupContext setup_context; 19 User user = *(User *)ctx; 20 MPI_Comm comm = PETSC_COMM_WORLD; 21 22 PetscFunctionBeginUser; 23 ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr); 24 // ------------------------------------------------------ 25 // SET UP DENSITY_CURRENT 26 // ------------------------------------------------------ 27 problem->ics.qfunction = ICsDC; 28 problem->ics.qfunction_loc = ICsDC_loc; 29 problem->bc = Exact_DC; 30 setup_context = problem->bc_ctx; 31 32 // ------------------------------------------------------ 33 // Create the libCEED context 34 // ------------------------------------------------------ 35 CeedScalar theta0 = 300.; // K 36 CeedScalar thetaC = -15.; // K 37 CeedScalar P0 = 1.e5; // Pa 38 CeedScalar N = 0.01; // 1/s 39 CeedScalar rc = 1000.; // m (Radius of bubble) 40 PetscReal center[3], dc_axis[3] = {0, 0, 0}; 41 PetscReal domain_min[3], domain_max[3], domain_size[3]; 42 ierr = DMGetBoundingBox(dm, domain_min, domain_max); 43 CHKERRQ(ierr); 44 for (PetscInt i = 0; i < 3; i++) 45 domain_size[i] = domain_max[i] - domain_min[i]; 46 47 // ------------------------------------------------------ 48 // Command line Options 49 // ------------------------------------------------------ 50 PetscOptionsBegin(comm, NULL, "Options for DENSITY_CURRENT problem", NULL); 51 ierr = PetscOptionsScalar("-theta0", "Reference potential temperature", NULL, 52 theta0, &theta0, NULL); 53 CHKERRQ(ierr); 54 ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature", 55 NULL, thetaC, &thetaC, NULL); 56 CHKERRQ(ierr); 57 ierr = PetscOptionsScalar("-P0", "Atmospheric pressure", NULL, P0, &P0, NULL); 58 CHKERRQ(ierr); 59 ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency", NULL, N, &N, NULL); 60 CHKERRQ(ierr); 61 ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", 62 NULL, rc, &rc, NULL); 63 CHKERRQ(ierr); 64 for (PetscInt i = 0; i < 3; i++) 65 center[i] = .5 * domain_size[i]; 66 PetscInt n = problem->dim; 67 ierr = PetscOptionsRealArray("-center", "Location of bubble center", NULL, 68 center, &n, NULL); 69 CHKERRQ(ierr); 70 n = problem->dim; 71 ierr = PetscOptionsRealArray("-dc_axis", 72 "Axis of density current cylindrical anomaly, " 73 "or {0,0,0} for spherically symmetric", 74 NULL, dc_axis, &n, NULL); 75 CHKERRQ(ierr); 76 { 77 PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) + PetscSqr(dc_axis[1]) + 78 PetscSqr(dc_axis[2])); 79 if (norm > 0) { 80 for (PetscInt i = 0; i < 3; i++) 81 dc_axis[i] /= norm; 82 } 83 } 84 85 PetscOptionsEnd(); 86 87 PetscScalar meter = user->units->meter; 88 PetscScalar second = user->units->second; 89 PetscScalar Kelvin = user->units->Kelvin; 90 PetscScalar Pascal = user->units->Pascal; 91 rc = fabs(rc) * meter; 92 theta0 *= Kelvin; 93 thetaC *= Kelvin; 94 P0 *= Pascal; 95 N *= (1. / second); 96 for (PetscInt i = 0; i < 3; i++) 97 center[i] *= meter; 98 99 setup_context->theta0 = theta0; 100 setup_context->thetaC = thetaC; 101 setup_context->P0 = P0; 102 setup_context->N = N; 103 setup_context->rc = rc; 104 setup_context->center[0] = center[0]; 105 setup_context->center[1] = center[1]; 106 setup_context->center[2] = center[2]; 107 setup_context->dc_axis[0] = dc_axis[0]; 108 setup_context->dc_axis[1] = dc_axis[1]; 109 setup_context->dc_axis[2] = dc_axis[2]; 110 111 PetscFunctionReturn(0); 112 } 113