1*88b783a1SJames Wright // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2*88b783a1SJames Wright // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3*88b783a1SJames Wright // reserved. See files LICENSE and NOTICE for details. 4*88b783a1SJames Wright // 5*88b783a1SJames Wright // This file is part of CEED, a collection of benchmarks, miniapps, software 6*88b783a1SJames Wright // libraries and APIs for efficient high-order finite element and spectral 7*88b783a1SJames Wright // element discretizations for exascale applications. For more information and 8*88b783a1SJames Wright // source code availability see http://github.com/ceed. 9*88b783a1SJames Wright // 10*88b783a1SJames Wright // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11*88b783a1SJames Wright // a collaborative effort of two U.S. Department of Energy organizations (Office 12*88b783a1SJames Wright // of Science and the National Nuclear Security Administration) responsible for 13*88b783a1SJames Wright // the planning and preparation of a capable exascale ecosystem, including 14*88b783a1SJames Wright // software, applications, hardware, advanced system engineering and early 15*88b783a1SJames Wright // testbed platforms, in support of the nation's exascale computing imperative. 16*88b783a1SJames Wright 17*88b783a1SJames Wright /// @file 18*88b783a1SJames Wright /// Utility functions for setting up problems using the Newtonian Qfunction 19*88b783a1SJames Wright 20*88b783a1SJames Wright #include "../navierstokes.h" 21*88b783a1SJames Wright #include "../qfunctions/setupgeo.h" 22*88b783a1SJames Wright #include "../qfunctions/newtonian.h" 23*88b783a1SJames Wright 24*88b783a1SJames Wright PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *setup_ctx, 25*88b783a1SJames Wright void *ctx) { 26*88b783a1SJames Wright SetupContext setup_context = *(SetupContext *)setup_ctx; 27*88b783a1SJames Wright User user = *(User *)ctx; 28*88b783a1SJames Wright StabilizationType stab; 29*88b783a1SJames Wright MPI_Comm comm = PETSC_COMM_WORLD; 30*88b783a1SJames Wright PetscBool implicit; 31*88b783a1SJames Wright PetscBool has_curr_time = PETSC_FALSE; 32*88b783a1SJames Wright PetscInt ierr; 33*88b783a1SJames Wright PetscFunctionBeginUser; 34*88b783a1SJames Wright 35*88b783a1SJames Wright ierr = PetscCalloc1(1, &user->phys->newtonian_ig_ctx); CHKERRQ(ierr); 36*88b783a1SJames Wright 37*88b783a1SJames Wright // ------------------------------------------------------ 38*88b783a1SJames Wright // Setup Generic Newtonian IG Problem 39*88b783a1SJames Wright // ------------------------------------------------------ 40*88b783a1SJames Wright problem->dim = 3; 41*88b783a1SJames Wright problem->q_data_size_vol = 10; 42*88b783a1SJames Wright problem->q_data_size_sur = 4; 43*88b783a1SJames Wright problem->setup_vol = Setup; 44*88b783a1SJames Wright problem->setup_vol_loc = Setup_loc; 45*88b783a1SJames Wright problem->ics = ICsNewtonianIG; 46*88b783a1SJames Wright problem->ics_loc = ICsNewtonianIG_loc; 47*88b783a1SJames Wright problem->setup_sur = SetupBoundary; 48*88b783a1SJames Wright problem->setup_sur_loc = SetupBoundary_loc; 49*88b783a1SJames Wright problem->apply_vol_rhs = Newtonian; 50*88b783a1SJames Wright problem->apply_vol_rhs_loc = Newtonian_loc; 51*88b783a1SJames Wright problem->apply_vol_ifunction = IFunction_Newtonian; 52*88b783a1SJames Wright problem->apply_vol_ifunction_loc = IFunction_Newtonian_loc; 53*88b783a1SJames Wright problem->setup_ctx = SetupContext_DENSITY_CURRENT; 54*88b783a1SJames Wright problem->non_zero_time = PETSC_FALSE; 55*88b783a1SJames Wright problem->print_info = PRINT_DENSITY_CURRENT; 56*88b783a1SJames Wright 57*88b783a1SJames Wright // ------------------------------------------------------ 58*88b783a1SJames Wright // Create the libCEED context 59*88b783a1SJames Wright // ------------------------------------------------------ 60*88b783a1SJames Wright CeedScalar theta0 = 300.; // K 61*88b783a1SJames Wright CeedScalar thetaC = -15.; // K 62*88b783a1SJames Wright CeedScalar P0 = 1.e5; // Pa 63*88b783a1SJames Wright CeedScalar N = 0.01; // 1/s 64*88b783a1SJames Wright CeedScalar cv = 717.; // J/(kg K) 65*88b783a1SJames Wright CeedScalar cp = 1004.; // J/(kg K) 66*88b783a1SJames Wright CeedScalar g = 9.81; // m/s^2 67*88b783a1SJames Wright CeedScalar lambda = -2./3.; // - 68*88b783a1SJames Wright CeedScalar mu = 75.; // Pa s, dynamic viscosity 69*88b783a1SJames Wright // mu = 75 is not physical for air, but is good for numerical stability 70*88b783a1SJames Wright CeedScalar k = 0.02638; // W/(m K) 71*88b783a1SJames Wright CeedScalar c_tau = 0.5; // - 72*88b783a1SJames Wright // c_tau = 0.5 is reported as "optimal" in Hughes et al 2010 73*88b783a1SJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 74*88b783a1SJames Wright ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 75*88b783a1SJames Wright for (int i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 76*88b783a1SJames Wright 77*88b783a1SJames Wright // ------------------------------------------------------ 78*88b783a1SJames Wright // Create the PETSc context 79*88b783a1SJames Wright // ------------------------------------------------------ 80*88b783a1SJames Wright PetscScalar meter = 1e-2; // 1 meter in scaled length units 81*88b783a1SJames Wright PetscScalar kilogram = 1e-6; // 1 kilogram in scaled mass units 82*88b783a1SJames Wright PetscScalar second = 1e-2; // 1 second in scaled time units 83*88b783a1SJames Wright PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 84*88b783a1SJames Wright PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s; 85*88b783a1SJames Wright 86*88b783a1SJames Wright // ------------------------------------------------------ 87*88b783a1SJames Wright // Command line Options 88*88b783a1SJames Wright // ------------------------------------------------------ 89*88b783a1SJames Wright ierr = PetscOptionsBegin(comm, NULL, 90*88b783a1SJames Wright "Options for Newtonian Ideal Gas based problem", 91*88b783a1SJames Wright NULL); CHKERRQ(ierr); 92*88b783a1SJames Wright // -- Physics 93*88b783a1SJames Wright ierr = PetscOptionsScalar("-theta0", "Reference potential temperature", 94*88b783a1SJames Wright NULL, theta0, &theta0, NULL); CHKERRQ(ierr); 95*88b783a1SJames Wright ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature", 96*88b783a1SJames Wright NULL, thetaC, &thetaC, NULL); CHKERRQ(ierr); 97*88b783a1SJames Wright ierr = PetscOptionsScalar("-P0", "Atmospheric pressure", 98*88b783a1SJames Wright NULL, P0, &P0, NULL); CHKERRQ(ierr); 99*88b783a1SJames Wright ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency", 100*88b783a1SJames Wright NULL, N, &N, NULL); CHKERRQ(ierr); 101*88b783a1SJames Wright ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume", 102*88b783a1SJames Wright NULL, cv, &cv, NULL); CHKERRQ(ierr); 103*88b783a1SJames Wright ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure", 104*88b783a1SJames Wright NULL, cp, &cp, NULL); CHKERRQ(ierr); 105*88b783a1SJames Wright ierr = PetscOptionsScalar("-g", "Gravitational acceleration", 106*88b783a1SJames Wright NULL, g, &g, NULL); CHKERRQ(ierr); 107*88b783a1SJames Wright ierr = PetscOptionsScalar("-lambda", 108*88b783a1SJames Wright "Stokes hypothesis second viscosity coefficient", 109*88b783a1SJames Wright NULL, lambda, &lambda, NULL); CHKERRQ(ierr); 110*88b783a1SJames Wright ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", 111*88b783a1SJames Wright NULL, mu, &mu, NULL); CHKERRQ(ierr); 112*88b783a1SJames Wright ierr = PetscOptionsScalar("-k", "Thermal conductivity", 113*88b783a1SJames Wright NULL, k, &k, NULL); CHKERRQ(ierr); 114*88b783a1SJames Wright 115*88b783a1SJames Wright ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL, 116*88b783a1SJames Wright StabilizationTypes, (PetscEnum)(stab = STAB_NONE), 117*88b783a1SJames Wright (PetscEnum *)&stab, NULL); CHKERRQ(ierr); 118*88b783a1SJames Wright ierr = PetscOptionsScalar("-c_tau", "Stabilization constant", 119*88b783a1SJames Wright NULL, c_tau, &c_tau, NULL); CHKERRQ(ierr); 120*88b783a1SJames Wright ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", 121*88b783a1SJames Wright NULL, implicit=PETSC_FALSE, &implicit, NULL); 122*88b783a1SJames Wright CHKERRQ(ierr); 123*88b783a1SJames Wright 124*88b783a1SJames Wright // -- Units 125*88b783a1SJames Wright ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units", 126*88b783a1SJames Wright NULL, meter, &meter, NULL); CHKERRQ(ierr); 127*88b783a1SJames Wright meter = fabs(meter); 128*88b783a1SJames Wright ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units", 129*88b783a1SJames Wright NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr); 130*88b783a1SJames Wright kilogram = fabs(kilogram); 131*88b783a1SJames Wright ierr = PetscOptionsScalar("-units_second","1 second in scaled time units", 132*88b783a1SJames Wright NULL, second, &second, NULL); CHKERRQ(ierr); 133*88b783a1SJames Wright second = fabs(second); 134*88b783a1SJames Wright ierr = PetscOptionsScalar("-units_Kelvin", 135*88b783a1SJames Wright "1 Kelvin in scaled temperature units", 136*88b783a1SJames Wright NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr); 137*88b783a1SJames Wright Kelvin = fabs(Kelvin); 138*88b783a1SJames Wright 139*88b783a1SJames Wright // -- Warnings 140*88b783a1SJames Wright if (stab == STAB_SUPG && !implicit) { 141*88b783a1SJames Wright ierr = PetscPrintf(comm, 142*88b783a1SJames Wright "Warning! Use -stab supg only with -implicit\n"); 143*88b783a1SJames Wright CHKERRQ(ierr); 144*88b783a1SJames Wright } 145*88b783a1SJames Wright ierr = PetscOptionsEnd(); CHKERRQ(ierr); 146*88b783a1SJames Wright 147*88b783a1SJames Wright // ------------------------------------------------------ 148*88b783a1SJames Wright // Set up the PETSc context 149*88b783a1SJames Wright // ------------------------------------------------------ 150*88b783a1SJames Wright // -- Define derived units 151*88b783a1SJames Wright Pascal = kilogram / (meter * PetscSqr(second)); 152*88b783a1SJames Wright J_per_kg_K = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 153*88b783a1SJames Wright m_per_squared_s = meter / PetscSqr(second); 154*88b783a1SJames Wright W_per_m_K = kilogram * meter / (pow(second,3) * Kelvin); 155*88b783a1SJames Wright 156*88b783a1SJames Wright user->units->meter = meter; 157*88b783a1SJames Wright user->units->kilogram = kilogram; 158*88b783a1SJames Wright user->units->second = second; 159*88b783a1SJames Wright user->units->Kelvin = Kelvin; 160*88b783a1SJames Wright user->units->Pascal = Pascal; 161*88b783a1SJames Wright user->units->J_per_kg_K = J_per_kg_K; 162*88b783a1SJames Wright user->units->m_per_squared_s = m_per_squared_s; 163*88b783a1SJames Wright user->units->W_per_m_K = W_per_m_K; 164*88b783a1SJames Wright 165*88b783a1SJames Wright // ------------------------------------------------------ 166*88b783a1SJames Wright // Set up the libCEED context 167*88b783a1SJames Wright // ------------------------------------------------------ 168*88b783a1SJames Wright // -- Scale variables to desired units 169*88b783a1SJames Wright theta0 *= Kelvin; 170*88b783a1SJames Wright thetaC *= Kelvin; 171*88b783a1SJames Wright P0 *= Pascal; 172*88b783a1SJames Wright N *= (1./second); 173*88b783a1SJames Wright cv *= J_per_kg_K; 174*88b783a1SJames Wright cp *= J_per_kg_K; 175*88b783a1SJames Wright g *= m_per_squared_s; 176*88b783a1SJames Wright mu *= Pascal * second; 177*88b783a1SJames Wright k *= W_per_m_K; 178*88b783a1SJames Wright for (int i=0; i<3; i++) domain_size[i] *= meter; 179*88b783a1SJames Wright problem->dm_scale = meter; 180*88b783a1SJames Wright 181*88b783a1SJames Wright // -- Setup Context 182*88b783a1SJames Wright setup_context->theta0 = theta0; 183*88b783a1SJames Wright setup_context->thetaC = thetaC; 184*88b783a1SJames Wright setup_context->P0 = P0; 185*88b783a1SJames Wright setup_context->N = N; 186*88b783a1SJames Wright setup_context->cv = cv; 187*88b783a1SJames Wright setup_context->cp = cp; 188*88b783a1SJames Wright setup_context->g = g; 189*88b783a1SJames Wright setup_context->lx = domain_size[0]; 190*88b783a1SJames Wright setup_context->ly = domain_size[1]; 191*88b783a1SJames Wright setup_context->lz = domain_size[2]; 192*88b783a1SJames Wright setup_context->time = 0; 193*88b783a1SJames Wright 194*88b783a1SJames Wright // -- Solver Settings 195*88b783a1SJames Wright user->phys->stab = stab; 196*88b783a1SJames Wright user->phys->implicit = implicit; 197*88b783a1SJames Wright user->phys->has_curr_time = has_curr_time; 198*88b783a1SJames Wright 199*88b783a1SJames Wright // -- QFunction Context 200*88b783a1SJames Wright user->phys->newtonian_ig_ctx->lambda = lambda; 201*88b783a1SJames Wright user->phys->newtonian_ig_ctx->mu = mu; 202*88b783a1SJames Wright user->phys->newtonian_ig_ctx->k = k; 203*88b783a1SJames Wright user->phys->newtonian_ig_ctx->cv = cv; 204*88b783a1SJames Wright user->phys->newtonian_ig_ctx->cp = cp; 205*88b783a1SJames Wright user->phys->newtonian_ig_ctx->g = g; 206*88b783a1SJames Wright user->phys->newtonian_ig_ctx->c_tau = c_tau; 207*88b783a1SJames Wright user->phys->newtonian_ig_ctx->stabilization = stab; 208*88b783a1SJames Wright 209*88b783a1SJames Wright PetscFunctionReturn(0); 210*88b783a1SJames Wright } 211*88b783a1SJames Wright 212*88b783a1SJames Wright PetscErrorCode SetupContext_NEWTONIAN_IG(Ceed ceed, CeedData ceed_data, 213*88b783a1SJames Wright AppCtx app_ctx, SetupContext setup_ctx, Physics phys) { 214*88b783a1SJames Wright PetscFunctionBeginUser; 215*88b783a1SJames Wright CeedQFunctionContextCreate(ceed, &ceed_data->setup_context); 216*88b783a1SJames Wright CeedQFunctionContextSetData(ceed_data->setup_context, CEED_MEM_HOST, 217*88b783a1SJames Wright CEED_USE_POINTER, sizeof(*setup_ctx), setup_ctx); 218*88b783a1SJames Wright CeedQFunctionSetContext(ceed_data->qf_ics, ceed_data->setup_context); 219*88b783a1SJames Wright CeedQFunctionContextCreate(ceed, &ceed_data->newt_ig_context); 220*88b783a1SJames Wright CeedQFunctionContextSetData(ceed_data->newt_ig_context, CEED_MEM_HOST, 221*88b783a1SJames Wright CEED_USE_POINTER, 222*88b783a1SJames Wright sizeof(*phys->newtonian_ig_ctx), phys->newtonian_ig_ctx); 223*88b783a1SJames Wright if (ceed_data->qf_rhs_vol) 224*88b783a1SJames Wright CeedQFunctionSetContext(ceed_data->qf_rhs_vol, ceed_data->newt_ig_context); 225*88b783a1SJames Wright if (ceed_data->qf_ifunction_vol) 226*88b783a1SJames Wright CeedQFunctionSetContext(ceed_data->qf_ifunction_vol, 227*88b783a1SJames Wright ceed_data->newt_ig_context); 228*88b783a1SJames Wright PetscFunctionReturn(0); 229*88b783a1SJames Wright } 230