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