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