13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 388b783a1SJames Wright // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 588b783a1SJames Wright // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 788b783a1SJames Wright 888b783a1SJames Wright /// @file 988b783a1SJames Wright /// Utility functions for setting up problems using the Newtonian Qfunction 1088b783a1SJames Wright 1188b783a1SJames Wright #include "../navierstokes.h" 1288b783a1SJames Wright #include "../qfunctions/setupgeo.h" 1388b783a1SJames Wright #include "../qfunctions/newtonian.h" 1488b783a1SJames Wright 1588b783a1SJames Wright PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *setup_ctx, 1688b783a1SJames Wright void *ctx) { 1788b783a1SJames Wright SetupContext setup_context = *(SetupContext *)setup_ctx; 1888b783a1SJames Wright User user = *(User *)ctx; 1988b783a1SJames Wright StabilizationType stab; 2088b783a1SJames Wright MPI_Comm comm = PETSC_COMM_WORLD; 2188b783a1SJames Wright PetscBool implicit; 2288b783a1SJames Wright PetscBool has_curr_time = PETSC_FALSE; 2388b783a1SJames Wright PetscInt ierr; 2488b783a1SJames Wright PetscFunctionBeginUser; 2588b783a1SJames Wright 2688b783a1SJames Wright ierr = PetscCalloc1(1, &user->phys->newtonian_ig_ctx); CHKERRQ(ierr); 2788b783a1SJames Wright 2888b783a1SJames Wright // ------------------------------------------------------ 2988b783a1SJames Wright // Setup Generic Newtonian IG Problem 3088b783a1SJames Wright // ------------------------------------------------------ 3188b783a1SJames Wright problem->dim = 3; 3288b783a1SJames Wright problem->q_data_size_vol = 10; 3388b783a1SJames Wright problem->q_data_size_sur = 4; 3488b783a1SJames Wright problem->setup_vol = Setup; 3588b783a1SJames Wright problem->setup_vol_loc = Setup_loc; 3688b783a1SJames Wright problem->ics = ICsNewtonianIG; 3788b783a1SJames Wright problem->ics_loc = ICsNewtonianIG_loc; 3888b783a1SJames Wright problem->setup_sur = SetupBoundary; 3988b783a1SJames Wright problem->setup_sur_loc = SetupBoundary_loc; 4088b783a1SJames Wright problem->apply_vol_rhs = Newtonian; 4188b783a1SJames Wright problem->apply_vol_rhs_loc = Newtonian_loc; 4288b783a1SJames Wright problem->apply_vol_ifunction = IFunction_Newtonian; 4388b783a1SJames Wright problem->apply_vol_ifunction_loc = IFunction_Newtonian_loc; 4488b783a1SJames Wright problem->setup_ctx = SetupContext_DENSITY_CURRENT; 4588b783a1SJames Wright problem->non_zero_time = PETSC_FALSE; 4688b783a1SJames Wright problem->print_info = PRINT_DENSITY_CURRENT; 4788b783a1SJames Wright 4888b783a1SJames Wright // ------------------------------------------------------ 4988b783a1SJames Wright // Create the libCEED context 5088b783a1SJames Wright // ------------------------------------------------------ 5188b783a1SJames Wright CeedScalar theta0 = 300.; // K 5288b783a1SJames Wright CeedScalar thetaC = -15.; // K 5388b783a1SJames Wright CeedScalar P0 = 1.e5; // Pa 5488b783a1SJames Wright CeedScalar N = 0.01; // 1/s 5588b783a1SJames Wright CeedScalar cv = 717.; // J/(kg K) 5688b783a1SJames Wright CeedScalar cp = 1004.; // J/(kg K) 5788b783a1SJames Wright CeedScalar g = 9.81; // m/s^2 5888b783a1SJames Wright CeedScalar lambda = -2./3.; // - 5988b783a1SJames Wright CeedScalar mu = 75.; // Pa s, dynamic viscosity 6088b783a1SJames Wright // mu = 75 is not physical for air, but is good for numerical stability 6188b783a1SJames Wright CeedScalar k = 0.02638; // W/(m K) 6288b783a1SJames Wright CeedScalar c_tau = 0.5; // - 6388b783a1SJames Wright // c_tau = 0.5 is reported as "optimal" in Hughes et al 2010 6488b783a1SJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 6588b783a1SJames Wright ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 6688b783a1SJames Wright for (int i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 6788b783a1SJames Wright 6888b783a1SJames Wright // ------------------------------------------------------ 6988b783a1SJames Wright // Create the PETSc context 7088b783a1SJames Wright // ------------------------------------------------------ 7188b783a1SJames Wright PetscScalar meter = 1e-2; // 1 meter in scaled length units 7288b783a1SJames Wright PetscScalar kilogram = 1e-6; // 1 kilogram in scaled mass units 7388b783a1SJames Wright PetscScalar second = 1e-2; // 1 second in scaled time units 7488b783a1SJames Wright PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 7588b783a1SJames Wright PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s; 7688b783a1SJames Wright 7788b783a1SJames Wright // ------------------------------------------------------ 7888b783a1SJames Wright // Command line Options 7988b783a1SJames Wright // ------------------------------------------------------ 80*67490bc6SJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", 81*67490bc6SJeremy L Thompson NULL); 82*67490bc6SJeremy L Thompson 8388b783a1SJames Wright // -- Physics 8488b783a1SJames Wright ierr = PetscOptionsScalar("-theta0", "Reference potential temperature", 8588b783a1SJames Wright NULL, theta0, &theta0, NULL); CHKERRQ(ierr); 8688b783a1SJames Wright ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature", 8788b783a1SJames Wright NULL, thetaC, &thetaC, NULL); CHKERRQ(ierr); 8888b783a1SJames Wright ierr = PetscOptionsScalar("-P0", "Atmospheric pressure", 8988b783a1SJames Wright NULL, P0, &P0, NULL); CHKERRQ(ierr); 9088b783a1SJames Wright ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency", 9188b783a1SJames Wright NULL, N, &N, NULL); CHKERRQ(ierr); 9288b783a1SJames Wright ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume", 9388b783a1SJames Wright NULL, cv, &cv, NULL); CHKERRQ(ierr); 9488b783a1SJames Wright ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure", 9588b783a1SJames Wright NULL, cp, &cp, NULL); CHKERRQ(ierr); 9688b783a1SJames Wright ierr = PetscOptionsScalar("-g", "Gravitational acceleration", 9788b783a1SJames Wright NULL, g, &g, NULL); CHKERRQ(ierr); 9888b783a1SJames Wright ierr = PetscOptionsScalar("-lambda", 9988b783a1SJames Wright "Stokes hypothesis second viscosity coefficient", 10088b783a1SJames Wright NULL, lambda, &lambda, NULL); CHKERRQ(ierr); 10188b783a1SJames Wright ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", 10288b783a1SJames Wright NULL, mu, &mu, NULL); CHKERRQ(ierr); 10388b783a1SJames Wright ierr = PetscOptionsScalar("-k", "Thermal conductivity", 10488b783a1SJames Wright NULL, k, &k, NULL); CHKERRQ(ierr); 10588b783a1SJames Wright 10688b783a1SJames Wright ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL, 10788b783a1SJames Wright StabilizationTypes, (PetscEnum)(stab = STAB_NONE), 10888b783a1SJames Wright (PetscEnum *)&stab, NULL); CHKERRQ(ierr); 10988b783a1SJames Wright ierr = PetscOptionsScalar("-c_tau", "Stabilization constant", 11088b783a1SJames Wright NULL, c_tau, &c_tau, NULL); CHKERRQ(ierr); 11188b783a1SJames Wright ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", 11288b783a1SJames Wright NULL, implicit=PETSC_FALSE, &implicit, NULL); 11388b783a1SJames Wright CHKERRQ(ierr); 11488b783a1SJames Wright 11588b783a1SJames Wright // -- Units 11688b783a1SJames Wright ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units", 11788b783a1SJames Wright NULL, meter, &meter, NULL); CHKERRQ(ierr); 11888b783a1SJames Wright meter = fabs(meter); 11988b783a1SJames Wright ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units", 12088b783a1SJames Wright NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr); 12188b783a1SJames Wright kilogram = fabs(kilogram); 12288b783a1SJames Wright ierr = PetscOptionsScalar("-units_second","1 second in scaled time units", 12388b783a1SJames Wright NULL, second, &second, NULL); CHKERRQ(ierr); 12488b783a1SJames Wright second = fabs(second); 12588b783a1SJames Wright ierr = PetscOptionsScalar("-units_Kelvin", 12688b783a1SJames Wright "1 Kelvin in scaled temperature units", 12788b783a1SJames Wright NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr); 12888b783a1SJames Wright Kelvin = fabs(Kelvin); 12988b783a1SJames Wright 13088b783a1SJames Wright // -- Warnings 13188b783a1SJames Wright if (stab == STAB_SUPG && !implicit) { 13288b783a1SJames Wright ierr = PetscPrintf(comm, 13388b783a1SJames Wright "Warning! Use -stab supg only with -implicit\n"); 13488b783a1SJames Wright CHKERRQ(ierr); 13588b783a1SJames Wright } 136*67490bc6SJeremy L Thompson PetscOptionsEnd(); 13788b783a1SJames Wright 13888b783a1SJames Wright // ------------------------------------------------------ 13988b783a1SJames Wright // Set up the PETSc context 14088b783a1SJames Wright // ------------------------------------------------------ 14188b783a1SJames Wright // -- Define derived units 14288b783a1SJames Wright Pascal = kilogram / (meter * PetscSqr(second)); 14388b783a1SJames Wright J_per_kg_K = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 14488b783a1SJames Wright m_per_squared_s = meter / PetscSqr(second); 14588b783a1SJames Wright W_per_m_K = kilogram * meter / (pow(second,3) * Kelvin); 14688b783a1SJames Wright 14788b783a1SJames Wright user->units->meter = meter; 14888b783a1SJames Wright user->units->kilogram = kilogram; 14988b783a1SJames Wright user->units->second = second; 15088b783a1SJames Wright user->units->Kelvin = Kelvin; 15188b783a1SJames Wright user->units->Pascal = Pascal; 15288b783a1SJames Wright user->units->J_per_kg_K = J_per_kg_K; 15388b783a1SJames Wright user->units->m_per_squared_s = m_per_squared_s; 15488b783a1SJames Wright user->units->W_per_m_K = W_per_m_K; 15588b783a1SJames Wright 15688b783a1SJames Wright // ------------------------------------------------------ 15788b783a1SJames Wright // Set up the libCEED context 15888b783a1SJames Wright // ------------------------------------------------------ 15988b783a1SJames Wright // -- Scale variables to desired units 16088b783a1SJames Wright theta0 *= Kelvin; 16188b783a1SJames Wright thetaC *= Kelvin; 16288b783a1SJames Wright P0 *= Pascal; 16388b783a1SJames Wright N *= (1./second); 16488b783a1SJames Wright cv *= J_per_kg_K; 16588b783a1SJames Wright cp *= J_per_kg_K; 16688b783a1SJames Wright g *= m_per_squared_s; 16788b783a1SJames Wright mu *= Pascal * second; 16888b783a1SJames Wright k *= W_per_m_K; 16988b783a1SJames Wright for (int i=0; i<3; i++) domain_size[i] *= meter; 17088b783a1SJames Wright problem->dm_scale = meter; 17188b783a1SJames Wright 17288b783a1SJames Wright // -- Setup Context 17388b783a1SJames Wright setup_context->theta0 = theta0; 17488b783a1SJames Wright setup_context->thetaC = thetaC; 17588b783a1SJames Wright setup_context->P0 = P0; 17688b783a1SJames Wright setup_context->N = N; 17788b783a1SJames Wright setup_context->cv = cv; 17888b783a1SJames Wright setup_context->cp = cp; 17988b783a1SJames Wright setup_context->g = g; 18088b783a1SJames Wright setup_context->lx = domain_size[0]; 18188b783a1SJames Wright setup_context->ly = domain_size[1]; 18288b783a1SJames Wright setup_context->lz = domain_size[2]; 18388b783a1SJames Wright setup_context->time = 0; 18488b783a1SJames Wright 18588b783a1SJames Wright // -- Solver Settings 18688b783a1SJames Wright user->phys->stab = stab; 18788b783a1SJames Wright user->phys->implicit = implicit; 18888b783a1SJames Wright user->phys->has_curr_time = has_curr_time; 18988b783a1SJames Wright 19088b783a1SJames Wright // -- QFunction Context 19188b783a1SJames Wright user->phys->newtonian_ig_ctx->lambda = lambda; 19288b783a1SJames Wright user->phys->newtonian_ig_ctx->mu = mu; 19388b783a1SJames Wright user->phys->newtonian_ig_ctx->k = k; 19488b783a1SJames Wright user->phys->newtonian_ig_ctx->cv = cv; 19588b783a1SJames Wright user->phys->newtonian_ig_ctx->cp = cp; 19688b783a1SJames Wright user->phys->newtonian_ig_ctx->g = g; 19788b783a1SJames Wright user->phys->newtonian_ig_ctx->c_tau = c_tau; 19888b783a1SJames Wright user->phys->newtonian_ig_ctx->stabilization = stab; 19988b783a1SJames Wright 20088b783a1SJames Wright PetscFunctionReturn(0); 20188b783a1SJames Wright } 20288b783a1SJames Wright 20388b783a1SJames Wright PetscErrorCode SetupContext_NEWTONIAN_IG(Ceed ceed, CeedData ceed_data, 20488b783a1SJames Wright AppCtx app_ctx, SetupContext setup_ctx, Physics phys) { 20588b783a1SJames Wright PetscFunctionBeginUser; 20688b783a1SJames Wright CeedQFunctionContextCreate(ceed, &ceed_data->setup_context); 20788b783a1SJames Wright CeedQFunctionContextSetData(ceed_data->setup_context, CEED_MEM_HOST, 20888b783a1SJames Wright CEED_USE_POINTER, sizeof(*setup_ctx), setup_ctx); 20988b783a1SJames Wright CeedQFunctionSetContext(ceed_data->qf_ics, ceed_data->setup_context); 21088b783a1SJames Wright CeedQFunctionContextCreate(ceed, &ceed_data->newt_ig_context); 21188b783a1SJames Wright CeedQFunctionContextSetData(ceed_data->newt_ig_context, CEED_MEM_HOST, 21288b783a1SJames Wright CEED_USE_POINTER, 21388b783a1SJames Wright sizeof(*phys->newtonian_ig_ctx), phys->newtonian_ig_ctx); 21488b783a1SJames Wright if (ceed_data->qf_rhs_vol) 21588b783a1SJames Wright CeedQFunctionSetContext(ceed_data->qf_rhs_vol, ceed_data->newt_ig_context); 21688b783a1SJames Wright if (ceed_data->qf_ifunction_vol) 21788b783a1SJames Wright CeedQFunctionSetContext(ceed_data->qf_ifunction_vol, 21888b783a1SJames Wright ceed_data->newt_ig_context); 21988b783a1SJames Wright PetscFunctionReturn(0); 22088b783a1SJames Wright } 221