xref: /honee/problems/newtonian.c (revision 727da7e70824329c32c70307643ee0d0901369c8)
1*727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2*727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
33a8779fbSJames Wright //
4*727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
53a8779fbSJames Wright //
6*727da7e7SJeremy 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 
153a8779fbSJames Wright PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *setup_ctx,
163a8779fbSJames Wright                                void *ctx) {
173a8779fbSJames Wright   SetupContext      setup_context = *(SetupContext *)setup_ctx;
183a8779fbSJames Wright   User              user = *(User *)ctx;
193a8779fbSJames Wright   StabilizationType stab;
203a8779fbSJames Wright   MPI_Comm          comm = PETSC_COMM_WORLD;
213a8779fbSJames Wright   PetscBool         implicit;
223a8779fbSJames Wright   PetscBool         has_curr_time = PETSC_FALSE;
233a8779fbSJames Wright   PetscInt          ierr;
243a8779fbSJames Wright   PetscFunctionBeginUser;
253a8779fbSJames Wright 
263a8779fbSJames Wright   ierr = PetscCalloc1(1, &user->phys->newtonian_ig_ctx); CHKERRQ(ierr);
273a8779fbSJames Wright 
283a8779fbSJames Wright   // ------------------------------------------------------
293a8779fbSJames Wright   //           Setup Generic Newtonian IG Problem
303a8779fbSJames Wright   // ------------------------------------------------------
313a8779fbSJames Wright   problem->dim                     = 3;
323a8779fbSJames Wright   problem->q_data_size_vol         = 10;
333a8779fbSJames Wright   problem->q_data_size_sur         = 4;
343a8779fbSJames Wright   problem->setup_vol               = Setup;
353a8779fbSJames Wright   problem->setup_vol_loc           = Setup_loc;
363a8779fbSJames Wright   problem->ics                     = ICsNewtonianIG;
373a8779fbSJames Wright   problem->ics_loc                 = ICsNewtonianIG_loc;
383a8779fbSJames Wright   problem->setup_sur               = SetupBoundary;
393a8779fbSJames Wright   problem->setup_sur_loc           = SetupBoundary_loc;
403a8779fbSJames Wright   problem->apply_vol_rhs           = Newtonian;
413a8779fbSJames Wright   problem->apply_vol_rhs_loc       = Newtonian_loc;
423a8779fbSJames Wright   problem->apply_vol_ifunction     = IFunction_Newtonian;
433a8779fbSJames Wright   problem->apply_vol_ifunction_loc = IFunction_Newtonian_loc;
443a8779fbSJames Wright   problem->setup_ctx               = SetupContext_DENSITY_CURRENT;
453a8779fbSJames Wright   problem->non_zero_time           = PETSC_FALSE;
463a8779fbSJames Wright   problem->print_info              = PRINT_DENSITY_CURRENT;
473a8779fbSJames Wright 
483a8779fbSJames Wright   // ------------------------------------------------------
493a8779fbSJames Wright   //             Create the libCEED context
503a8779fbSJames Wright   // ------------------------------------------------------
513a8779fbSJames Wright   CeedScalar theta0 = 300.;    // K
523a8779fbSJames Wright   CeedScalar thetaC = -15.;    // K
533a8779fbSJames Wright   CeedScalar P0     = 1.e5;    // Pa
543a8779fbSJames Wright   CeedScalar N      = 0.01;    // 1/s
553a8779fbSJames Wright   CeedScalar cv     = 717.;    // J/(kg K)
563a8779fbSJames Wright   CeedScalar cp     = 1004.;   // J/(kg K)
573a8779fbSJames Wright   CeedScalar g      = 9.81;    // m/s^2
583a8779fbSJames Wright   CeedScalar lambda = -2./3.;  // -
593a8779fbSJames Wright   CeedScalar mu     = 75.;     // Pa s, dynamic viscosity
603a8779fbSJames Wright   // mu = 75 is not physical for air, but is good for numerical stability
613a8779fbSJames Wright   CeedScalar k      = 0.02638; // W/(m K)
623a8779fbSJames Wright   CeedScalar c_tau  = 0.5;     // -
633a8779fbSJames Wright   // c_tau = 0.5 is reported as "optimal" in Hughes et al 2010
643a8779fbSJames Wright   PetscReal domain_min[3], domain_max[3], domain_size[3];
653a8779fbSJames Wright   ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
663a8779fbSJames Wright   for (int i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i];
673a8779fbSJames Wright 
683a8779fbSJames Wright   // ------------------------------------------------------
693a8779fbSJames Wright   //             Create the PETSc context
703a8779fbSJames Wright   // ------------------------------------------------------
713a8779fbSJames Wright   PetscScalar meter    = 1e-2;  // 1 meter in scaled length units
723a8779fbSJames Wright   PetscScalar kilogram = 1e-6;  // 1 kilogram in scaled mass units
733a8779fbSJames Wright   PetscScalar second   = 1e-2;  // 1 second in scaled time units
743a8779fbSJames Wright   PetscScalar Kelvin   = 1;     // 1 Kelvin in scaled temperature units
753a8779fbSJames Wright   PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s;
763a8779fbSJames Wright 
773a8779fbSJames Wright   // ------------------------------------------------------
783a8779fbSJames Wright   //              Command line Options
793a8779fbSJames Wright   // ------------------------------------------------------
803a8779fbSJames Wright   ierr = PetscOptionsBegin(comm, NULL,
813a8779fbSJames Wright                            "Options for Newtonian Ideal Gas based problem",
823a8779fbSJames Wright                            NULL); CHKERRQ(ierr);
833a8779fbSJames Wright   // -- Physics
843a8779fbSJames Wright   ierr = PetscOptionsScalar("-theta0", "Reference potential temperature",
853a8779fbSJames Wright                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
863a8779fbSJames Wright   ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature",
873a8779fbSJames Wright                             NULL, thetaC, &thetaC, NULL); CHKERRQ(ierr);
883a8779fbSJames Wright   ierr = PetscOptionsScalar("-P0", "Atmospheric pressure",
893a8779fbSJames Wright                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
903a8779fbSJames Wright   ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency",
913a8779fbSJames Wright                             NULL, N, &N, NULL); CHKERRQ(ierr);
923a8779fbSJames Wright   ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume",
933a8779fbSJames Wright                             NULL, cv, &cv, NULL); CHKERRQ(ierr);
943a8779fbSJames Wright   ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure",
953a8779fbSJames Wright                             NULL, cp, &cp, NULL); CHKERRQ(ierr);
963a8779fbSJames Wright   ierr = PetscOptionsScalar("-g", "Gravitational acceleration",
973a8779fbSJames Wright                             NULL, g, &g, NULL); CHKERRQ(ierr);
983a8779fbSJames Wright   ierr = PetscOptionsScalar("-lambda",
993a8779fbSJames Wright                             "Stokes hypothesis second viscosity coefficient",
1003a8779fbSJames Wright                             NULL, lambda, &lambda, NULL); CHKERRQ(ierr);
1013a8779fbSJames Wright   ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient",
1023a8779fbSJames Wright                             NULL, mu, &mu, NULL); CHKERRQ(ierr);
1033a8779fbSJames Wright   ierr = PetscOptionsScalar("-k", "Thermal conductivity",
1043a8779fbSJames Wright                             NULL, k, &k, NULL); CHKERRQ(ierr);
1053a8779fbSJames Wright 
1063a8779fbSJames Wright   ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL,
1073a8779fbSJames Wright                           StabilizationTypes, (PetscEnum)(stab = STAB_NONE),
1083a8779fbSJames Wright                           (PetscEnum *)&stab, NULL); CHKERRQ(ierr);
1093a8779fbSJames Wright   ierr = PetscOptionsScalar("-c_tau", "Stabilization constant",
1103a8779fbSJames Wright                             NULL, c_tau, &c_tau, NULL); CHKERRQ(ierr);
1113a8779fbSJames Wright   ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation",
1123a8779fbSJames Wright                           NULL, implicit=PETSC_FALSE, &implicit, NULL);
1133a8779fbSJames Wright   CHKERRQ(ierr);
1143a8779fbSJames Wright 
1153a8779fbSJames Wright   // -- Units
1163a8779fbSJames Wright   ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units",
1173a8779fbSJames Wright                             NULL, meter, &meter, NULL); CHKERRQ(ierr);
1183a8779fbSJames Wright   meter = fabs(meter);
1193a8779fbSJames Wright   ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units",
1203a8779fbSJames Wright                             NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr);
1213a8779fbSJames Wright   kilogram = fabs(kilogram);
1223a8779fbSJames Wright   ierr = PetscOptionsScalar("-units_second","1 second in scaled time units",
1233a8779fbSJames Wright                             NULL, second, &second, NULL); CHKERRQ(ierr);
1243a8779fbSJames Wright   second = fabs(second);
1253a8779fbSJames Wright   ierr = PetscOptionsScalar("-units_Kelvin",
1263a8779fbSJames Wright                             "1 Kelvin in scaled temperature units",
1273a8779fbSJames Wright                             NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr);
1283a8779fbSJames Wright   Kelvin = fabs(Kelvin);
1293a8779fbSJames Wright 
1303a8779fbSJames Wright   // -- Warnings
1313a8779fbSJames Wright   if (stab == STAB_SUPG && !implicit) {
1323a8779fbSJames Wright     ierr = PetscPrintf(comm,
1333a8779fbSJames Wright                        "Warning! Use -stab supg only with -implicit\n");
1343a8779fbSJames Wright     CHKERRQ(ierr);
1353a8779fbSJames Wright   }
1363a8779fbSJames Wright   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
1373a8779fbSJames Wright 
1383a8779fbSJames Wright   // ------------------------------------------------------
1393a8779fbSJames Wright   //           Set up the PETSc context
1403a8779fbSJames Wright   // ------------------------------------------------------
1413a8779fbSJames Wright   // -- Define derived units
1423a8779fbSJames Wright   Pascal          = kilogram / (meter * PetscSqr(second));
1433a8779fbSJames Wright   J_per_kg_K      =  PetscSqr(meter) / (PetscSqr(second) * Kelvin);
1443a8779fbSJames Wright   m_per_squared_s = meter / PetscSqr(second);
1453a8779fbSJames Wright   W_per_m_K       = kilogram * meter / (pow(second,3) * Kelvin);
1463a8779fbSJames Wright 
1473a8779fbSJames Wright   user->units->meter           = meter;
1483a8779fbSJames Wright   user->units->kilogram        = kilogram;
1493a8779fbSJames Wright   user->units->second          = second;
1503a8779fbSJames Wright   user->units->Kelvin          = Kelvin;
1513a8779fbSJames Wright   user->units->Pascal          = Pascal;
1523a8779fbSJames Wright   user->units->J_per_kg_K      = J_per_kg_K;
1533a8779fbSJames Wright   user->units->m_per_squared_s = m_per_squared_s;
1543a8779fbSJames Wright   user->units->W_per_m_K       = W_per_m_K;
1553a8779fbSJames Wright 
1563a8779fbSJames Wright   // ------------------------------------------------------
1573a8779fbSJames Wright   //           Set up the libCEED context
1583a8779fbSJames Wright   // ------------------------------------------------------
1593a8779fbSJames Wright   // -- Scale variables to desired units
1603a8779fbSJames Wright   theta0 *= Kelvin;
1613a8779fbSJames Wright   thetaC *= Kelvin;
1623a8779fbSJames Wright   P0     *= Pascal;
1633a8779fbSJames Wright   N      *= (1./second);
1643a8779fbSJames Wright   cv     *= J_per_kg_K;
1653a8779fbSJames Wright   cp     *= J_per_kg_K;
1663a8779fbSJames Wright   g      *= m_per_squared_s;
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;
1703a8779fbSJames Wright   problem->dm_scale = meter;
1713a8779fbSJames Wright 
1723a8779fbSJames Wright   // -- Setup Context
1733a8779fbSJames Wright   setup_context->theta0     = theta0;
1743a8779fbSJames Wright   setup_context->thetaC     = thetaC;
1753a8779fbSJames Wright   setup_context->P0         = P0;
1763a8779fbSJames Wright   setup_context->N          = N;
1773a8779fbSJames Wright   setup_context->cv         = cv;
1783a8779fbSJames Wright   setup_context->cp         = cp;
1793a8779fbSJames Wright   setup_context->g          = g;
1803a8779fbSJames Wright   setup_context->lx         = domain_size[0];
1813a8779fbSJames Wright   setup_context->ly         = domain_size[1];
1823a8779fbSJames Wright   setup_context->lz         = domain_size[2];
1833a8779fbSJames Wright   setup_context->time       = 0;
1843a8779fbSJames Wright 
1853a8779fbSJames Wright   // -- Solver Settings
1863a8779fbSJames Wright   user->phys->stab          = stab;
1873a8779fbSJames Wright   user->phys->implicit      = implicit;
1883a8779fbSJames Wright   user->phys->has_curr_time = has_curr_time;
1893a8779fbSJames Wright 
1903a8779fbSJames Wright   // -- QFunction Context
1913a8779fbSJames Wright   user->phys->newtonian_ig_ctx->lambda        = lambda;
1923a8779fbSJames Wright   user->phys->newtonian_ig_ctx->mu            = mu;
1933a8779fbSJames Wright   user->phys->newtonian_ig_ctx->k             = k;
1943a8779fbSJames Wright   user->phys->newtonian_ig_ctx->cv            = cv;
1953a8779fbSJames Wright   user->phys->newtonian_ig_ctx->cp            = cp;
1963a8779fbSJames Wright   user->phys->newtonian_ig_ctx->g             = g;
1973a8779fbSJames Wright   user->phys->newtonian_ig_ctx->c_tau         = c_tau;
1983a8779fbSJames Wright   user->phys->newtonian_ig_ctx->stabilization = stab;
1993a8779fbSJames Wright 
2003a8779fbSJames Wright   PetscFunctionReturn(0);
2013a8779fbSJames Wright }
2023a8779fbSJames Wright 
2033a8779fbSJames Wright PetscErrorCode SetupContext_NEWTONIAN_IG(Ceed ceed, CeedData ceed_data,
2043a8779fbSJames Wright     AppCtx app_ctx, SetupContext setup_ctx, Physics phys) {
2053a8779fbSJames Wright   PetscFunctionBeginUser;
2063a8779fbSJames Wright   CeedQFunctionContextCreate(ceed, &ceed_data->setup_context);
2073a8779fbSJames Wright   CeedQFunctionContextSetData(ceed_data->setup_context, CEED_MEM_HOST,
2083a8779fbSJames Wright                               CEED_USE_POINTER, sizeof(*setup_ctx), setup_ctx);
2093a8779fbSJames Wright   CeedQFunctionSetContext(ceed_data->qf_ics, ceed_data->setup_context);
2103a8779fbSJames Wright   CeedQFunctionContextCreate(ceed, &ceed_data->newt_ig_context);
2113a8779fbSJames Wright   CeedQFunctionContextSetData(ceed_data->newt_ig_context, CEED_MEM_HOST,
2123a8779fbSJames Wright                               CEED_USE_POINTER,
2133a8779fbSJames Wright                               sizeof(*phys->newtonian_ig_ctx), phys->newtonian_ig_ctx);
2143a8779fbSJames Wright   if (ceed_data->qf_rhs_vol)
2153a8779fbSJames Wright     CeedQFunctionSetContext(ceed_data->qf_rhs_vol, ceed_data->newt_ig_context);
2163a8779fbSJames Wright   if (ceed_data->qf_ifunction_vol)
2173a8779fbSJames Wright     CeedQFunctionSetContext(ceed_data->qf_ifunction_vol,
2183a8779fbSJames Wright                             ceed_data->newt_ig_context);
2193a8779fbSJames Wright   PetscFunctionReturn(0);
2203a8779fbSJames Wright }
221