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