1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 33a8779fbSJames Wright // 4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 53a8779fbSJames Wright // 6727da7e7SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 73a8779fbSJames Wright 83a8779fbSJames Wright /// @file 93a8779fbSJames Wright /// Utility functions for setting up problems using the Newtonian Qfunction 103a8779fbSJames Wright 113a8779fbSJames Wright #include "../navierstokes.h" 123a8779fbSJames Wright #include "../qfunctions/setupgeo.h" 133a8779fbSJames Wright #include "../qfunctions/newtonian.h" 143a8779fbSJames Wright 15bb8a0c61SJames Wright 16bb8a0c61SJames Wright #ifndef newtonian_context_struct 17bb8a0c61SJames Wright #define newtonian_context_struct 18bb8a0c61SJames Wright typedef struct NewtonianIdealGasContext_ *NewtonianIdealGasContext; 19bb8a0c61SJames Wright struct NewtonianIdealGasContext_ { 20bb8a0c61SJames Wright CeedScalar lambda; 21bb8a0c61SJames Wright CeedScalar mu; 22bb8a0c61SJames Wright CeedScalar k; 23bb8a0c61SJames Wright CeedScalar cv; 24bb8a0c61SJames Wright CeedScalar cp; 25bb8a0c61SJames Wright CeedScalar g[3]; 26bb8a0c61SJames Wright CeedScalar c_tau; 27bb8a0c61SJames Wright CeedScalar Ctau_t; 28bb8a0c61SJames Wright CeedScalar Ctau_v; 29bb8a0c61SJames Wright CeedScalar Ctau_C; 30bb8a0c61SJames Wright CeedScalar Ctau_M; 31bb8a0c61SJames Wright CeedScalar Ctau_E; 32bb8a0c61SJames Wright CeedScalar dt; 33bb8a0c61SJames Wright StabilizationType stabilization; 34bb8a0c61SJames Wright }; 35bb8a0c61SJames Wright #endif 36bb8a0c61SJames Wright 373a8779fbSJames Wright PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *setup_ctx, 383a8779fbSJames Wright void *ctx) { 393a8779fbSJames Wright SetupContext setup_context = *(SetupContext *)setup_ctx; 403a8779fbSJames Wright User user = *(User *)ctx; 413a8779fbSJames Wright StabilizationType stab; 423a8779fbSJames Wright MPI_Comm comm = PETSC_COMM_WORLD; 433a8779fbSJames Wright PetscBool implicit; 443a8779fbSJames Wright PetscBool has_curr_time = PETSC_FALSE; 453a8779fbSJames Wright PetscInt ierr; 463a8779fbSJames Wright PetscFunctionBeginUser; 473a8779fbSJames Wright 483a8779fbSJames Wright ierr = PetscCalloc1(1, &user->phys->newtonian_ig_ctx); CHKERRQ(ierr); 493a8779fbSJames Wright 503a8779fbSJames Wright // ------------------------------------------------------ 513a8779fbSJames Wright // Setup Generic Newtonian IG Problem 523a8779fbSJames Wright // ------------------------------------------------------ 533a8779fbSJames Wright problem->dim = 3; 543a8779fbSJames Wright problem->q_data_size_vol = 10; 553a8779fbSJames Wright problem->q_data_size_sur = 4; 56*9785fe93SJed Brown problem->setup_vol.qfunction = Setup; 57*9785fe93SJed Brown problem->setup_vol.qfunction_loc = Setup_loc; 58*9785fe93SJed Brown problem->ics.qfunction = ICsNewtonianIG; 59*9785fe93SJed Brown problem->ics.qfunction_loc = ICsNewtonianIG_loc; 60*9785fe93SJed Brown problem->setup_sur.qfunction = SetupBoundary; 61*9785fe93SJed Brown problem->setup_sur.qfunction_loc = SetupBoundary_loc; 62*9785fe93SJed Brown problem->apply_vol_rhs.qfunction = Newtonian; 63*9785fe93SJed Brown problem->apply_vol_rhs.qfunction_loc = Newtonian_loc; 64*9785fe93SJed Brown problem->apply_vol_ifunction.qfunction = IFunction_Newtonian; 65*9785fe93SJed Brown problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_loc; 663a8779fbSJames Wright problem->setup_ctx = SetupContext_DENSITY_CURRENT; 673a8779fbSJames Wright problem->non_zero_time = PETSC_FALSE; 683a8779fbSJames Wright problem->print_info = PRINT_DENSITY_CURRENT; 693a8779fbSJames Wright 703a8779fbSJames Wright // ------------------------------------------------------ 713a8779fbSJames Wright // Create the libCEED context 723a8779fbSJames Wright // ------------------------------------------------------ 733a8779fbSJames Wright CeedScalar cv = 717.; // J/(kg K) 743a8779fbSJames Wright CeedScalar cp = 1004.; // J/(kg K) 75bb8a0c61SJames Wright CeedScalar g[3] = {0, 0, -9.81}; // m/s^2 763a8779fbSJames Wright CeedScalar lambda = -2./3.; // - 77bb8a0c61SJames Wright CeedScalar mu = 1.8e-5; // Pa s, dynamic viscosity 783a8779fbSJames Wright CeedScalar k = 0.02638; // W/(m K) 793a8779fbSJames Wright CeedScalar c_tau = 0.5; // - 80bb8a0c61SJames Wright CeedScalar Ctau_t = 1.0; // - 81bb8a0c61SJames Wright CeedScalar Ctau_v = 36.0; // TODO make function of degree 82bb8a0c61SJames Wright CeedScalar Ctau_C = 1.0; // TODO make function of degree 83bb8a0c61SJames Wright CeedScalar Ctau_M = 1.0; // TODO make function of degree 84bb8a0c61SJames Wright CeedScalar Ctau_E = 1.0; // TODO make function of degree 853a8779fbSJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 863a8779fbSJames Wright ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 873a8779fbSJames Wright for (int i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 883a8779fbSJames Wright 893a8779fbSJames Wright // ------------------------------------------------------ 903a8779fbSJames Wright // Create the PETSc context 913a8779fbSJames Wright // ------------------------------------------------------ 92bb8a0c61SJames Wright PetscScalar meter = 1; // 1 meter in scaled length units 93bb8a0c61SJames Wright PetscScalar kilogram = 1; // 1 kilogram in scaled mass units 94bb8a0c61SJames Wright PetscScalar second = 1; // 1 second in scaled time units 953a8779fbSJames Wright PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 963a8779fbSJames Wright PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s; 973a8779fbSJames Wright 983a8779fbSJames Wright // ------------------------------------------------------ 993a8779fbSJames Wright // Command line Options 1003a8779fbSJames Wright // ------------------------------------------------------ 1011485969bSJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", 1021485969bSJeremy L Thompson NULL); 1031485969bSJeremy L Thompson 1043a8779fbSJames Wright // -- Physics 1053a8779fbSJames Wright ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume", 1063a8779fbSJames Wright NULL, cv, &cv, NULL); CHKERRQ(ierr); 1073a8779fbSJames Wright ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure", 1083a8779fbSJames Wright NULL, cp, &cp, NULL); CHKERRQ(ierr); 1093a8779fbSJames Wright ierr = PetscOptionsScalar("-lambda", 1103a8779fbSJames Wright "Stokes hypothesis second viscosity coefficient", 1113a8779fbSJames Wright NULL, lambda, &lambda, NULL); CHKERRQ(ierr); 1123a8779fbSJames Wright ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", 1133a8779fbSJames Wright NULL, mu, &mu, NULL); CHKERRQ(ierr); 1143a8779fbSJames Wright ierr = PetscOptionsScalar("-k", "Thermal conductivity", 1153a8779fbSJames Wright NULL, k, &k, NULL); CHKERRQ(ierr); 1163a8779fbSJames Wright 117bb8a0c61SJames Wright PetscInt dim = problem->dim; 118bb8a0c61SJames Wright ierr = PetscOptionsRealArray("-g", "Gravitational acceleration", 119bb8a0c61SJames Wright NULL, g, &dim, NULL); CHKERRQ(ierr); 1203a8779fbSJames Wright ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL, 1213a8779fbSJames Wright StabilizationTypes, (PetscEnum)(stab = STAB_NONE), 1223a8779fbSJames Wright (PetscEnum *)&stab, NULL); CHKERRQ(ierr); 1233a8779fbSJames Wright ierr = PetscOptionsScalar("-c_tau", "Stabilization constant", 1243a8779fbSJames Wright NULL, c_tau, &c_tau, NULL); CHKERRQ(ierr); 125bb8a0c61SJames Wright ierr = PetscOptionsScalar("-Ctau_t", "Stabilization time constant", 126bb8a0c61SJames Wright NULL, Ctau_t, &Ctau_t, NULL); CHKERRQ(ierr); 127bb8a0c61SJames Wright ierr = PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", 128bb8a0c61SJames Wright NULL, Ctau_v, &Ctau_v, NULL); CHKERRQ(ierr); 129bb8a0c61SJames Wright ierr = PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", 130bb8a0c61SJames Wright NULL, Ctau_C, &Ctau_C, NULL); CHKERRQ(ierr); 131bb8a0c61SJames Wright ierr = PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", 132bb8a0c61SJames Wright NULL, Ctau_M, &Ctau_M, NULL); CHKERRQ(ierr); 133bb8a0c61SJames Wright ierr = PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", 134bb8a0c61SJames Wright NULL, Ctau_E, &Ctau_E, NULL); CHKERRQ(ierr); 1353a8779fbSJames Wright ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", 1363a8779fbSJames Wright NULL, implicit=PETSC_FALSE, &implicit, NULL); 1373a8779fbSJames Wright CHKERRQ(ierr); 1383a8779fbSJames Wright 1393a8779fbSJames Wright // -- Units 1403a8779fbSJames Wright ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units", 1413a8779fbSJames Wright NULL, meter, &meter, NULL); CHKERRQ(ierr); 1423a8779fbSJames Wright meter = fabs(meter); 1433a8779fbSJames Wright ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units", 1443a8779fbSJames Wright NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr); 1453a8779fbSJames Wright kilogram = fabs(kilogram); 1463a8779fbSJames Wright ierr = PetscOptionsScalar("-units_second","1 second in scaled time units", 1473a8779fbSJames Wright NULL, second, &second, NULL); CHKERRQ(ierr); 1483a8779fbSJames Wright second = fabs(second); 1493a8779fbSJames Wright ierr = PetscOptionsScalar("-units_Kelvin", 1503a8779fbSJames Wright "1 Kelvin in scaled temperature units", 1513a8779fbSJames Wright NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr); 1523a8779fbSJames Wright Kelvin = fabs(Kelvin); 1533a8779fbSJames Wright 1543a8779fbSJames Wright // -- Warnings 1553a8779fbSJames Wright if (stab == STAB_SUPG && !implicit) { 1563a8779fbSJames Wright ierr = PetscPrintf(comm, 1573a8779fbSJames Wright "Warning! Use -stab supg only with -implicit\n"); 1583a8779fbSJames Wright CHKERRQ(ierr); 1593a8779fbSJames Wright } 1601485969bSJeremy L Thompson PetscOptionsEnd(); 1613a8779fbSJames Wright 1623a8779fbSJames Wright // ------------------------------------------------------ 1633a8779fbSJames Wright // Set up the PETSc context 1643a8779fbSJames Wright // ------------------------------------------------------ 1653a8779fbSJames Wright // -- Define derived units 1663a8779fbSJames Wright Pascal = kilogram / (meter * PetscSqr(second)); 1673a8779fbSJames Wright J_per_kg_K = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 1683a8779fbSJames Wright m_per_squared_s = meter / PetscSqr(second); 1693a8779fbSJames Wright W_per_m_K = kilogram * meter / (pow(second,3) * Kelvin); 1703a8779fbSJames Wright 1713a8779fbSJames Wright user->units->meter = meter; 1723a8779fbSJames Wright user->units->kilogram = kilogram; 1733a8779fbSJames Wright user->units->second = second; 1743a8779fbSJames Wright user->units->Kelvin = Kelvin; 1753a8779fbSJames Wright user->units->Pascal = Pascal; 1763a8779fbSJames Wright user->units->J_per_kg_K = J_per_kg_K; 1773a8779fbSJames Wright user->units->m_per_squared_s = m_per_squared_s; 1783a8779fbSJames Wright user->units->W_per_m_K = W_per_m_K; 1793a8779fbSJames Wright 1803a8779fbSJames Wright // ------------------------------------------------------ 1813a8779fbSJames Wright // Set up the libCEED context 1823a8779fbSJames Wright // ------------------------------------------------------ 1833a8779fbSJames Wright // -- Scale variables to desired units 1843a8779fbSJames Wright cv *= J_per_kg_K; 1853a8779fbSJames Wright cp *= J_per_kg_K; 1863a8779fbSJames Wright mu *= Pascal * second; 1873a8779fbSJames Wright k *= W_per_m_K; 1883a8779fbSJames Wright for (int i=0; i<3; i++) domain_size[i] *= meter; 189bb8a0c61SJames Wright for (int i=0; i<3; i++) g[i] *= m_per_squared_s; 1903a8779fbSJames Wright problem->dm_scale = meter; 1913a8779fbSJames Wright 1923a8779fbSJames Wright // -- Setup Context 1933a8779fbSJames Wright setup_context->cv = cv; 1943a8779fbSJames Wright setup_context->cp = cp; 1953a8779fbSJames Wright setup_context->lx = domain_size[0]; 1963a8779fbSJames Wright setup_context->ly = domain_size[1]; 1973a8779fbSJames Wright setup_context->lz = domain_size[2]; 1983a8779fbSJames Wright setup_context->time = 0; 199bb8a0c61SJames Wright ierr = PetscArraycpy(setup_context->g, g, 3); CHKERRQ(ierr); 2003a8779fbSJames Wright 2013a8779fbSJames Wright // -- Solver Settings 2023a8779fbSJames Wright user->phys->stab = stab; 2033a8779fbSJames Wright user->phys->implicit = implicit; 2043a8779fbSJames Wright user->phys->has_curr_time = has_curr_time; 2053a8779fbSJames Wright 2063a8779fbSJames Wright // -- QFunction Context 2073a8779fbSJames Wright user->phys->newtonian_ig_ctx->lambda = lambda; 2083a8779fbSJames Wright user->phys->newtonian_ig_ctx->mu = mu; 2093a8779fbSJames Wright user->phys->newtonian_ig_ctx->k = k; 2103a8779fbSJames Wright user->phys->newtonian_ig_ctx->cv = cv; 2113a8779fbSJames Wright user->phys->newtonian_ig_ctx->cp = cp; 2123a8779fbSJames Wright user->phys->newtonian_ig_ctx->c_tau = c_tau; 213bb8a0c61SJames Wright user->phys->newtonian_ig_ctx->Ctau_t = Ctau_t; 214bb8a0c61SJames Wright user->phys->newtonian_ig_ctx->Ctau_v = Ctau_v; 215bb8a0c61SJames Wright user->phys->newtonian_ig_ctx->Ctau_C = Ctau_C; 216bb8a0c61SJames Wright user->phys->newtonian_ig_ctx->Ctau_M = Ctau_M; 217bb8a0c61SJames Wright user->phys->newtonian_ig_ctx->Ctau_E = Ctau_E; 2183a8779fbSJames Wright user->phys->newtonian_ig_ctx->stabilization = stab; 219bb8a0c61SJames Wright ierr = PetscArraycpy(user->phys->newtonian_ig_ctx->g, g, 3); CHKERRQ(ierr); 2203a8779fbSJames Wright 2213a8779fbSJames Wright PetscFunctionReturn(0); 2223a8779fbSJames Wright } 2233a8779fbSJames Wright 2243a8779fbSJames Wright PetscErrorCode SetupContext_NEWTONIAN_IG(Ceed ceed, CeedData ceed_data, 2253a8779fbSJames Wright AppCtx app_ctx, SetupContext setup_ctx, Physics phys) { 2263a8779fbSJames Wright PetscFunctionBeginUser; 2273a8779fbSJames Wright CeedQFunctionContextCreate(ceed, &ceed_data->setup_context); 2283a8779fbSJames Wright CeedQFunctionContextSetData(ceed_data->setup_context, CEED_MEM_HOST, 2293a8779fbSJames Wright CEED_USE_POINTER, sizeof(*setup_ctx), setup_ctx); 2303a8779fbSJames Wright CeedQFunctionSetContext(ceed_data->qf_ics, ceed_data->setup_context); 231bb8a0c61SJames Wright 2323a8779fbSJames Wright CeedQFunctionContextCreate(ceed, &ceed_data->newt_ig_context); 2333a8779fbSJames Wright CeedQFunctionContextSetData(ceed_data->newt_ig_context, CEED_MEM_HOST, 2343a8779fbSJames Wright CEED_USE_POINTER, 2353a8779fbSJames Wright sizeof(*phys->newtonian_ig_ctx), phys->newtonian_ig_ctx); 236bb8a0c61SJames Wright CeedQFunctionContextRegisterDouble(ceed_data->newt_ig_context, "timestep size", 237bb8a0c61SJames Wright offsetof(struct NewtonianIdealGasContext_, dt), 1, "Size of timestep, delta t"); 238bb8a0c61SJames Wright 2393a8779fbSJames Wright if (ceed_data->qf_rhs_vol) 2403a8779fbSJames Wright CeedQFunctionSetContext(ceed_data->qf_rhs_vol, ceed_data->newt_ig_context); 2413a8779fbSJames Wright if (ceed_data->qf_ifunction_vol) 2423a8779fbSJames Wright CeedQFunctionSetContext(ceed_data->qf_ifunction_vol, 2433a8779fbSJames Wright ceed_data->newt_ig_context); 2443a8779fbSJames Wright PetscFunctionReturn(0); 2453a8779fbSJames Wright } 246