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 15*b7f03f12SJed Brown PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *ctx) { 16*b7f03f12SJed Brown SetupContext setup_context; 173a8779fbSJames Wright User user = *(User *)ctx; 183a8779fbSJames Wright StabilizationType stab; 193a8779fbSJames Wright MPI_Comm comm = PETSC_COMM_WORLD; 203a8779fbSJames Wright PetscBool implicit; 213a8779fbSJames Wright PetscBool has_curr_time = PETSC_FALSE; 223a8779fbSJames Wright PetscInt ierr; 2315a3537eSJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 2415a3537eSJed Brown CeedQFunctionContext newtonian_ig_context; 253a8779fbSJames Wright 2615a3537eSJed Brown PetscFunctionBeginUser; 27*b7f03f12SJed Brown ierr = PetscCalloc1(1, &setup_context); CHKERRQ(ierr); 2815a3537eSJed Brown ierr = PetscCalloc1(1, &newtonian_ig_ctx); CHKERRQ(ierr); 293a8779fbSJames Wright 303a8779fbSJames Wright // ------------------------------------------------------ 313a8779fbSJames Wright // Setup Generic Newtonian IG Problem 323a8779fbSJames Wright // ------------------------------------------------------ 333a8779fbSJames Wright problem->dim = 3; 343a8779fbSJames Wright problem->q_data_size_vol = 10; 353a8779fbSJames Wright problem->q_data_size_sur = 4; 369785fe93SJed Brown problem->setup_vol.qfunction = Setup; 379785fe93SJed Brown problem->setup_vol.qfunction_loc = Setup_loc; 389785fe93SJed Brown problem->ics.qfunction = ICsNewtonianIG; 399785fe93SJed Brown problem->ics.qfunction_loc = ICsNewtonianIG_loc; 409785fe93SJed Brown problem->setup_sur.qfunction = SetupBoundary; 419785fe93SJed Brown problem->setup_sur.qfunction_loc = SetupBoundary_loc; 429785fe93SJed Brown problem->apply_vol_rhs.qfunction = Newtonian; 439785fe93SJed Brown problem->apply_vol_rhs.qfunction_loc = Newtonian_loc; 449785fe93SJed Brown problem->apply_vol_ifunction.qfunction = IFunction_Newtonian; 459785fe93SJed Brown problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_loc; 46*b7f03f12SJed Brown problem->bc = NULL; 47*b7f03f12SJed Brown problem->bc_ctx = setup_context; 483a8779fbSJames Wright problem->non_zero_time = PETSC_FALSE; 493a8779fbSJames Wright problem->print_info = PRINT_DENSITY_CURRENT; 503a8779fbSJames Wright 513a8779fbSJames Wright // ------------------------------------------------------ 523a8779fbSJames Wright // Create the libCEED context 533a8779fbSJames Wright // ------------------------------------------------------ 543a8779fbSJames Wright CeedScalar cv = 717.; // J/(kg K) 553a8779fbSJames Wright CeedScalar cp = 1004.; // J/(kg K) 56bb8a0c61SJames Wright CeedScalar g[3] = {0, 0, -9.81}; // m/s^2 573a8779fbSJames Wright CeedScalar lambda = -2./3.; // - 58bb8a0c61SJames Wright CeedScalar mu = 1.8e-5; // Pa s, dynamic viscosity 593a8779fbSJames Wright CeedScalar k = 0.02638; // W/(m K) 603a8779fbSJames Wright CeedScalar c_tau = 0.5; // - 61bb8a0c61SJames Wright CeedScalar Ctau_t = 1.0; // - 62bb8a0c61SJames Wright CeedScalar Ctau_v = 36.0; // TODO make function of degree 63bb8a0c61SJames Wright CeedScalar Ctau_C = 1.0; // TODO make function of degree 64bb8a0c61SJames Wright CeedScalar Ctau_M = 1.0; // TODO make function of degree 65bb8a0c61SJames Wright CeedScalar Ctau_E = 1.0; // TODO make function of degree 663a8779fbSJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 673a8779fbSJames Wright ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 683a8779fbSJames Wright for (int i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 693a8779fbSJames Wright 703a8779fbSJames Wright // ------------------------------------------------------ 713a8779fbSJames Wright // Create the PETSc context 723a8779fbSJames Wright // ------------------------------------------------------ 73bb8a0c61SJames Wright PetscScalar meter = 1; // 1 meter in scaled length units 74bb8a0c61SJames Wright PetscScalar kilogram = 1; // 1 kilogram in scaled mass units 75bb8a0c61SJames Wright PetscScalar second = 1; // 1 second in scaled time units 763a8779fbSJames Wright PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 773a8779fbSJames Wright PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s; 783a8779fbSJames Wright 793a8779fbSJames Wright // ------------------------------------------------------ 803a8779fbSJames Wright // Command line Options 813a8779fbSJames Wright // ------------------------------------------------------ 821485969bSJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", 831485969bSJeremy L Thompson NULL); 841485969bSJeremy L Thompson 853a8779fbSJames Wright // -- Physics 863a8779fbSJames Wright ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume", 873a8779fbSJames Wright NULL, cv, &cv, NULL); CHKERRQ(ierr); 883a8779fbSJames Wright ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure", 893a8779fbSJames Wright NULL, cp, &cp, NULL); CHKERRQ(ierr); 903a8779fbSJames Wright ierr = PetscOptionsScalar("-lambda", 913a8779fbSJames Wright "Stokes hypothesis second viscosity coefficient", 923a8779fbSJames Wright NULL, lambda, &lambda, NULL); CHKERRQ(ierr); 933a8779fbSJames Wright ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", 943a8779fbSJames Wright NULL, mu, &mu, NULL); CHKERRQ(ierr); 953a8779fbSJames Wright ierr = PetscOptionsScalar("-k", "Thermal conductivity", 963a8779fbSJames Wright NULL, k, &k, NULL); CHKERRQ(ierr); 973a8779fbSJames Wright 98bb8a0c61SJames Wright PetscInt dim = problem->dim; 99bb8a0c61SJames Wright ierr = PetscOptionsRealArray("-g", "Gravitational acceleration", 100bb8a0c61SJames Wright NULL, g, &dim, NULL); CHKERRQ(ierr); 1013a8779fbSJames Wright ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL, 1023a8779fbSJames Wright StabilizationTypes, (PetscEnum)(stab = STAB_NONE), 1033a8779fbSJames Wright (PetscEnum *)&stab, NULL); CHKERRQ(ierr); 1043a8779fbSJames Wright ierr = PetscOptionsScalar("-c_tau", "Stabilization constant", 1053a8779fbSJames Wright NULL, c_tau, &c_tau, NULL); CHKERRQ(ierr); 106bb8a0c61SJames Wright ierr = PetscOptionsScalar("-Ctau_t", "Stabilization time constant", 107bb8a0c61SJames Wright NULL, Ctau_t, &Ctau_t, NULL); CHKERRQ(ierr); 108bb8a0c61SJames Wright ierr = PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", 109bb8a0c61SJames Wright NULL, Ctau_v, &Ctau_v, NULL); CHKERRQ(ierr); 110bb8a0c61SJames Wright ierr = PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", 111bb8a0c61SJames Wright NULL, Ctau_C, &Ctau_C, NULL); CHKERRQ(ierr); 112bb8a0c61SJames Wright ierr = PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", 113bb8a0c61SJames Wright NULL, Ctau_M, &Ctau_M, NULL); CHKERRQ(ierr); 114bb8a0c61SJames Wright ierr = PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", 115bb8a0c61SJames Wright NULL, Ctau_E, &Ctau_E, NULL); CHKERRQ(ierr); 1163a8779fbSJames Wright ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", 1173a8779fbSJames Wright NULL, implicit=PETSC_FALSE, &implicit, NULL); 1183a8779fbSJames Wright CHKERRQ(ierr); 1193a8779fbSJames Wright 1203a8779fbSJames Wright // -- Units 1213a8779fbSJames Wright ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units", 1223a8779fbSJames Wright NULL, meter, &meter, NULL); CHKERRQ(ierr); 1233a8779fbSJames Wright meter = fabs(meter); 1243a8779fbSJames Wright ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units", 1253a8779fbSJames Wright NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr); 1263a8779fbSJames Wright kilogram = fabs(kilogram); 1273a8779fbSJames Wright ierr = PetscOptionsScalar("-units_second","1 second in scaled time units", 1283a8779fbSJames Wright NULL, second, &second, NULL); CHKERRQ(ierr); 1293a8779fbSJames Wright second = fabs(second); 1303a8779fbSJames Wright ierr = PetscOptionsScalar("-units_Kelvin", 1313a8779fbSJames Wright "1 Kelvin in scaled temperature units", 1323a8779fbSJames Wright NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr); 1333a8779fbSJames Wright Kelvin = fabs(Kelvin); 1343a8779fbSJames Wright 1353a8779fbSJames Wright // -- Warnings 1363a8779fbSJames Wright if (stab == STAB_SUPG && !implicit) { 1373a8779fbSJames Wright ierr = PetscPrintf(comm, 1383a8779fbSJames Wright "Warning! Use -stab supg only with -implicit\n"); 1393a8779fbSJames Wright CHKERRQ(ierr); 1403a8779fbSJames Wright } 1411485969bSJeremy L Thompson PetscOptionsEnd(); 1423a8779fbSJames Wright 1433a8779fbSJames Wright // ------------------------------------------------------ 1443a8779fbSJames Wright // Set up the PETSc context 1453a8779fbSJames Wright // ------------------------------------------------------ 1463a8779fbSJames Wright // -- Define derived units 1473a8779fbSJames Wright Pascal = kilogram / (meter * PetscSqr(second)); 1483a8779fbSJames Wright J_per_kg_K = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 1493a8779fbSJames Wright m_per_squared_s = meter / PetscSqr(second); 1503a8779fbSJames Wright W_per_m_K = kilogram * meter / (pow(second,3) * Kelvin); 1513a8779fbSJames Wright 1523a8779fbSJames Wright user->units->meter = meter; 1533a8779fbSJames Wright user->units->kilogram = kilogram; 1543a8779fbSJames Wright user->units->second = second; 1553a8779fbSJames Wright user->units->Kelvin = Kelvin; 1563a8779fbSJames Wright user->units->Pascal = Pascal; 1573a8779fbSJames Wright user->units->J_per_kg_K = J_per_kg_K; 1583a8779fbSJames Wright user->units->m_per_squared_s = m_per_squared_s; 1593a8779fbSJames Wright user->units->W_per_m_K = W_per_m_K; 1603a8779fbSJames Wright 1613a8779fbSJames Wright // ------------------------------------------------------ 1623a8779fbSJames Wright // Set up the libCEED context 1633a8779fbSJames Wright // ------------------------------------------------------ 1643a8779fbSJames Wright // -- Scale variables to desired units 1653a8779fbSJames Wright cv *= J_per_kg_K; 1663a8779fbSJames Wright cp *= J_per_kg_K; 1673a8779fbSJames Wright mu *= Pascal * second; 1683a8779fbSJames Wright k *= W_per_m_K; 1693a8779fbSJames Wright for (int i=0; i<3; i++) domain_size[i] *= meter; 170bb8a0c61SJames Wright for (int i=0; i<3; i++) g[i] *= m_per_squared_s; 1713a8779fbSJames Wright problem->dm_scale = meter; 1723a8779fbSJames Wright 1733a8779fbSJames Wright // -- Setup Context 1743a8779fbSJames Wright setup_context->cv = cv; 1753a8779fbSJames Wright setup_context->cp = cp; 1763a8779fbSJames Wright setup_context->lx = domain_size[0]; 1773a8779fbSJames Wright setup_context->ly = domain_size[1]; 1783a8779fbSJames Wright setup_context->lz = domain_size[2]; 1793a8779fbSJames Wright setup_context->time = 0; 180bb8a0c61SJames Wright ierr = PetscArraycpy(setup_context->g, g, 3); CHKERRQ(ierr); 1813a8779fbSJames Wright 1823a8779fbSJames Wright // -- Solver Settings 1833a8779fbSJames Wright user->phys->stab = stab; 1843a8779fbSJames Wright user->phys->implicit = implicit; 1853a8779fbSJames Wright user->phys->has_curr_time = has_curr_time; 1863a8779fbSJames Wright 1873a8779fbSJames Wright // -- QFunction Context 18815a3537eSJed Brown newtonian_ig_ctx->lambda = lambda; 18915a3537eSJed Brown newtonian_ig_ctx->mu = mu; 19015a3537eSJed Brown newtonian_ig_ctx->k = k; 19115a3537eSJed Brown newtonian_ig_ctx->cv = cv; 19215a3537eSJed Brown newtonian_ig_ctx->cp = cp; 19315a3537eSJed Brown newtonian_ig_ctx->c_tau = c_tau; 19415a3537eSJed Brown newtonian_ig_ctx->Ctau_t = Ctau_t; 19515a3537eSJed Brown newtonian_ig_ctx->Ctau_v = Ctau_v; 19615a3537eSJed Brown newtonian_ig_ctx->Ctau_C = Ctau_C; 19715a3537eSJed Brown newtonian_ig_ctx->Ctau_M = Ctau_M; 19815a3537eSJed Brown newtonian_ig_ctx->Ctau_E = Ctau_E; 19915a3537eSJed Brown newtonian_ig_ctx->stabilization = stab; 20015a3537eSJed Brown ierr = PetscArraycpy(newtonian_ig_ctx->g, g, 3); CHKERRQ(ierr); 2013a8779fbSJames Wright 20215a3537eSJed Brown CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context); 20315a3537eSJed Brown CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, 20415a3537eSJed Brown CEED_USE_POINTER, sizeof(*setup_context), setup_context); 20515a3537eSJed Brown CeedQFunctionContextRegisterDouble(problem->ics.qfunction_context, 20615a3537eSJed Brown "evaluation time", 20715a3537eSJed Brown (char *)&setup_context->time - (char *)setup_context, 1, "Time of evaluation"); 20815a3537eSJed Brown 20915a3537eSJed Brown CeedQFunctionContextCreate(user->ceed, &newtonian_ig_context); 21015a3537eSJed Brown CeedQFunctionContextSetData(newtonian_ig_context, CEED_MEM_HOST, 21115a3537eSJed Brown CEED_USE_POINTER, 21215a3537eSJed Brown sizeof(*newtonian_ig_ctx), newtonian_ig_ctx); 21315a3537eSJed Brown CeedQFunctionContextSetDataDestroy(newtonian_ig_context, CEED_MEM_HOST, 21415a3537eSJed Brown FreeContextPetsc); 21515a3537eSJed Brown CeedQFunctionContextRegisterDouble(newtonian_ig_context, "timestep size", 21615a3537eSJed Brown offsetof(struct NewtonianIdealGasContext_, dt), 1, "Size of timestep, delta t"); 21715a3537eSJed Brown problem->apply_vol_rhs.qfunction_context = newtonian_ig_context; 21815a3537eSJed Brown CeedQFunctionContextReferenceCopy(newtonian_ig_context, 21915a3537eSJed Brown &problem->apply_vol_ifunction.qfunction_context); 2203a8779fbSJames Wright PetscFunctionReturn(0); 2213a8779fbSJames Wright } 2223a8779fbSJames Wright 22315a3537eSJed Brown PetscErrorCode PRINT_DENSITY_CURRENT(ProblemData *problem, 22415a3537eSJed Brown AppCtx app_ctx) { 22515a3537eSJed Brown MPI_Comm comm = PETSC_COMM_WORLD; 22615a3537eSJed Brown PetscErrorCode ierr; 22715a3537eSJed Brown NewtonianIdealGasContext newtonian_ctx; 22815a3537eSJed Brown 2293a8779fbSJames Wright PetscFunctionBeginUser; 23015a3537eSJed Brown CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, 23115a3537eSJed Brown CEED_MEM_HOST, &newtonian_ctx); 23215a3537eSJed Brown ierr = PetscPrintf(comm, 23315a3537eSJed Brown " Problem:\n" 23415a3537eSJed Brown " Problem Name : %s\n" 23515a3537eSJed Brown " Stabilization : %s\n", 23615a3537eSJed Brown app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization]); 23715a3537eSJed Brown CHKERRQ(ierr); 23815a3537eSJed Brown CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, 23915a3537eSJed Brown &newtonian_ctx); 2403a8779fbSJames Wright PetscFunctionReturn(0); 2413a8779fbSJames Wright } 242