xref: /honee/problems/densitycurrent.c (revision 05a512bda999ca8e4cad9319b5f88d1ff577e38b)
1a515125bSLeila Ghaffari // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2a515125bSLeila Ghaffari // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3a515125bSLeila Ghaffari // reserved. See files LICENSE and NOTICE for details.
4a515125bSLeila Ghaffari //
5a515125bSLeila Ghaffari // This file is part of CEED, a collection of benchmarks, miniapps, software
6a515125bSLeila Ghaffari // libraries and APIs for efficient high-order finite element and spectral
7a515125bSLeila Ghaffari // element discretizations for exascale applications. For more information and
8a515125bSLeila Ghaffari // source code availability see http://github.com/ceed.
9a515125bSLeila Ghaffari //
10a515125bSLeila Ghaffari // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11a515125bSLeila Ghaffari // a collaborative effort of two U.S. Department of Energy organizations (Office
12a515125bSLeila Ghaffari // of Science and the National Nuclear Security Administration) responsible for
13a515125bSLeila Ghaffari // the planning and preparation of a capable exascale ecosystem, including
14a515125bSLeila Ghaffari // software, applications, hardware, advanced system engineering and early
15a515125bSLeila Ghaffari // testbed platforms, in support of the nation's exascale computing imperative.
16a515125bSLeila Ghaffari 
17a515125bSLeila Ghaffari /// @file
18a515125bSLeila Ghaffari /// Utility functions for setting up DENSITY_CURRENT
19a515125bSLeila Ghaffari 
20a515125bSLeila Ghaffari #include "../navierstokes.h"
21a515125bSLeila Ghaffari #include "../qfunctions/setupgeo.h"
22a515125bSLeila Ghaffari #include "../qfunctions/densitycurrent.h"
23a515125bSLeila Ghaffari 
24*05a512bdSLeila Ghaffari PetscErrorCode NS_DENSITY_CURRENT(ProblemData *problem, DM dm, void *setup_ctx,
25a515125bSLeila Ghaffari                                   void *ctx) {
26a515125bSLeila Ghaffari   SetupContext      setup_context = *(SetupContext *)setup_ctx;
27a515125bSLeila Ghaffari   User              user = *(User *)ctx;
28a515125bSLeila Ghaffari   StabilizationType stab;
29a515125bSLeila Ghaffari   MPI_Comm          comm = PETSC_COMM_WORLD;
30a515125bSLeila Ghaffari   PetscBool         implicit;
31a515125bSLeila Ghaffari   PetscBool         has_curr_time = PETSC_FALSE;
32a515125bSLeila Ghaffari   PetscInt          ierr;
33a515125bSLeila Ghaffari   PetscFunctionBeginUser;
34a515125bSLeila Ghaffari 
35a515125bSLeila Ghaffari   ierr = PetscCalloc1(1, &user->phys->dc_ctx); CHKERRQ(ierr);
36a515125bSLeila Ghaffari 
37a515125bSLeila Ghaffari   // ------------------------------------------------------
38a515125bSLeila Ghaffari   //               SET UP DENSITY_CURRENT
39a515125bSLeila Ghaffari   // ------------------------------------------------------
40a515125bSLeila Ghaffari   problem->dim                     = 3;
41a515125bSLeila Ghaffari   problem->q_data_size_vol         = 10;
42a515125bSLeila Ghaffari   problem->q_data_size_sur         = 4;
43a515125bSLeila Ghaffari   problem->setup_vol               = Setup;
44a515125bSLeila Ghaffari   problem->setup_vol_loc           = Setup_loc;
45a515125bSLeila Ghaffari   problem->setup_sur               = SetupBoundary;
46a515125bSLeila Ghaffari   problem->setup_sur_loc           = SetupBoundary_loc;
47a515125bSLeila Ghaffari   problem->ics                     = ICsDC;
48a515125bSLeila Ghaffari   problem->ics_loc                 = ICsDC_loc;
49a515125bSLeila Ghaffari   problem->apply_vol_rhs           = DC;
50a515125bSLeila Ghaffari   problem->apply_vol_rhs_loc       = DC_loc;
51a515125bSLeila Ghaffari   problem->apply_vol_ifunction     = IFunction_DC;
52a515125bSLeila Ghaffari   problem->apply_vol_ifunction_loc = IFunction_DC_loc;
53a515125bSLeila Ghaffari   problem->bc                      = Exact_DC;
54ba5420e5SLeila Ghaffari   problem->setup_ctx               = SetupContext_DENSITY_CURRENT;
55a515125bSLeila Ghaffari   problem->bc_func                 = BC_DENSITY_CURRENT;
56a515125bSLeila Ghaffari   problem->non_zero_time           = PETSC_FALSE;
57a515125bSLeila Ghaffari   problem->print_info              = PRINT_DENSITY_CURRENT;
58a515125bSLeila Ghaffari 
59a515125bSLeila Ghaffari   // ------------------------------------------------------
60a515125bSLeila Ghaffari   //             Create the libCEED context
61a515125bSLeila Ghaffari   // ------------------------------------------------------
62a515125bSLeila Ghaffari   CeedScalar theta0 = 300.;    // K
63a515125bSLeila Ghaffari   CeedScalar thetaC = -15.;    // K
64a515125bSLeila Ghaffari   CeedScalar P0     = 1.e5;    // Pa
65a515125bSLeila Ghaffari   CeedScalar N      = 0.01;    // 1/s
66a515125bSLeila Ghaffari   CeedScalar cv     = 717.;    // J/(kg K)
67a515125bSLeila Ghaffari   CeedScalar cp     = 1004.;   // J/(kg K)
68a515125bSLeila Ghaffari   CeedScalar g      = 9.81;    // m/s^2
69a515125bSLeila Ghaffari   CeedScalar lambda = -2./3.;  // -
70a515125bSLeila Ghaffari   CeedScalar mu     = 75.;     // Pa s, dynamic viscosity
71a515125bSLeila Ghaffari   // mu = 75 is not physical for air, but is good for numerical stability
72a515125bSLeila Ghaffari   CeedScalar k      = 0.02638; // W/(m K)
73f821ee77SLeila Ghaffari   CeedScalar c_tau  = 0.5;     // -
74d8a22b9eSJed Brown   // c_tau = 0.5 is reported as "optimal" in Hughes et al 2010
75a515125bSLeila Ghaffari   CeedScalar rc     = 1000.;   // m (Radius of bubble)
76a515125bSLeila Ghaffari   PetscReal center[3], dc_axis[3] = {0, 0, 0};
77*05a512bdSLeila Ghaffari   PetscReal domain_min[3], domain_max[3], domain_size[3];
78*05a512bdSLeila Ghaffari   ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
79*05a512bdSLeila Ghaffari   for (int i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i];
80a515125bSLeila Ghaffari 
81a515125bSLeila Ghaffari   // ------------------------------------------------------
82a515125bSLeila Ghaffari   //             Create the PETSc context
83a515125bSLeila Ghaffari   // ------------------------------------------------------
84a515125bSLeila Ghaffari   PetscScalar meter    = 1e-2;  // 1 meter in scaled length units
85a515125bSLeila Ghaffari   PetscScalar kilogram = 1e-6;  // 1 kilogram in scaled mass units
86a515125bSLeila Ghaffari   PetscScalar second   = 1e-2;  // 1 second in scaled time units
87a515125bSLeila Ghaffari   PetscScalar Kelvin   = 1;     // 1 Kelvin in scaled temperature units
88a515125bSLeila Ghaffari   PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s;
89a515125bSLeila Ghaffari 
90a515125bSLeila Ghaffari   // ------------------------------------------------------
91a515125bSLeila Ghaffari   //              Command line Options
92a515125bSLeila Ghaffari   // ------------------------------------------------------
93a515125bSLeila Ghaffari   ierr = PetscOptionsBegin(comm, NULL, "Options for DENSITY_CURRENT problem",
94a515125bSLeila Ghaffari                            NULL); CHKERRQ(ierr);
95a515125bSLeila Ghaffari   // -- Physics
96a515125bSLeila Ghaffari   ierr = PetscOptionsScalar("-theta0", "Reference potential temperature",
97a515125bSLeila Ghaffari                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
98a515125bSLeila Ghaffari   ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature",
99a515125bSLeila Ghaffari                             NULL, thetaC, &thetaC, NULL); CHKERRQ(ierr);
100a515125bSLeila Ghaffari   ierr = PetscOptionsScalar("-P0", "Atmospheric pressure",
101a515125bSLeila Ghaffari                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
102a515125bSLeila Ghaffari   ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency",
103a515125bSLeila Ghaffari                             NULL, N, &N, NULL); CHKERRQ(ierr);
104a515125bSLeila Ghaffari   ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume",
105a515125bSLeila Ghaffari                             NULL, cv, &cv, NULL); CHKERRQ(ierr);
106a515125bSLeila Ghaffari   ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure",
107a515125bSLeila Ghaffari                             NULL, cp, &cp, NULL); CHKERRQ(ierr);
108a515125bSLeila Ghaffari   ierr = PetscOptionsScalar("-g", "Gravitational acceleration",
109a515125bSLeila Ghaffari                             NULL, g, &g, NULL); CHKERRQ(ierr);
110a515125bSLeila Ghaffari   ierr = PetscOptionsScalar("-lambda",
111a515125bSLeila Ghaffari                             "Stokes hypothesis second viscosity coefficient",
112a515125bSLeila Ghaffari                             NULL, lambda, &lambda, NULL); CHKERRQ(ierr);
113a515125bSLeila Ghaffari   ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient",
114a515125bSLeila Ghaffari                             NULL, mu, &mu, NULL); CHKERRQ(ierr);
115a515125bSLeila Ghaffari   ierr = PetscOptionsScalar("-k", "Thermal conductivity",
116a515125bSLeila Ghaffari                             NULL, k, &k, NULL); CHKERRQ(ierr);
117a515125bSLeila Ghaffari   ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble",
118a515125bSLeila Ghaffari                             NULL, rc, &rc, NULL); CHKERRQ(ierr);
119*05a512bdSLeila Ghaffari   for (int i=0; i<3; i++) center[i] = .5*domain_size[i];
120a515125bSLeila Ghaffari   PetscInt n = problem->dim;
121a515125bSLeila Ghaffari   ierr = PetscOptionsRealArray("-center", "Location of bubble center",
122a515125bSLeila Ghaffari                                NULL, center, &n, NULL); CHKERRQ(ierr);
123a515125bSLeila Ghaffari   n = problem->dim;
124a515125bSLeila Ghaffari   ierr = PetscOptionsRealArray("-dc_axis",
125a515125bSLeila Ghaffari                                "Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric",
126a515125bSLeila Ghaffari                                NULL, dc_axis, &n, NULL); CHKERRQ(ierr);
127a515125bSLeila Ghaffari   {
128a515125bSLeila Ghaffari     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) + PetscSqr(dc_axis[1]) +
129a515125bSLeila Ghaffari                                    PetscSqr(dc_axis[2]));
130a515125bSLeila Ghaffari     if (norm > 0) {
131a515125bSLeila Ghaffari       for (int i=0; i<3; i++)  dc_axis[i] /= norm;
132a515125bSLeila Ghaffari     }
133a515125bSLeila Ghaffari   }
134a515125bSLeila Ghaffari   ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL,
135a515125bSLeila Ghaffari                           StabilizationTypes, (PetscEnum)(stab = STAB_NONE),
136a515125bSLeila Ghaffari                           (PetscEnum *)&stab, NULL); CHKERRQ(ierr);
137d8a22b9eSJed Brown   ierr = PetscOptionsScalar("-c_tau", "Stabilization constant",
138d8a22b9eSJed Brown                             NULL, c_tau, &c_tau, NULL); CHKERRQ(ierr);
139a515125bSLeila Ghaffari   ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation",
140a515125bSLeila Ghaffari                           NULL, implicit=PETSC_FALSE, &implicit, NULL);
141a515125bSLeila Ghaffari   CHKERRQ(ierr);
142a515125bSLeila Ghaffari 
143a515125bSLeila Ghaffari   // -- Units
144a515125bSLeila Ghaffari   ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units",
145a515125bSLeila Ghaffari                             NULL, meter, &meter, NULL); CHKERRQ(ierr);
146a515125bSLeila Ghaffari   meter = fabs(meter);
147a515125bSLeila Ghaffari   ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units",
148a515125bSLeila Ghaffari                             NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr);
149a515125bSLeila Ghaffari   kilogram = fabs(kilogram);
150a515125bSLeila Ghaffari   ierr = PetscOptionsScalar("-units_second","1 second in scaled time units",
151a515125bSLeila Ghaffari                             NULL, second, &second, NULL); CHKERRQ(ierr);
152a515125bSLeila Ghaffari   second = fabs(second);
153a515125bSLeila Ghaffari   ierr = PetscOptionsScalar("-units_Kelvin",
154a515125bSLeila Ghaffari                             "1 Kelvin in scaled temperature units",
155a515125bSLeila Ghaffari                             NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr);
156a515125bSLeila Ghaffari   Kelvin = fabs(Kelvin);
157a515125bSLeila Ghaffari 
158a515125bSLeila Ghaffari   // -- Warnings
159a515125bSLeila Ghaffari   if (stab == STAB_SUPG && !implicit) {
160a515125bSLeila Ghaffari     ierr = PetscPrintf(comm,
161a515125bSLeila Ghaffari                        "Warning! Use -stab supg only with -implicit\n");
162a515125bSLeila Ghaffari     CHKERRQ(ierr);
163a515125bSLeila Ghaffari   }
164a515125bSLeila Ghaffari 
165a515125bSLeila Ghaffari   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
166a515125bSLeila Ghaffari 
167a515125bSLeila Ghaffari   // ------------------------------------------------------
168a515125bSLeila Ghaffari   //           Set up the PETSc context
169a515125bSLeila Ghaffari   // ------------------------------------------------------
170a515125bSLeila Ghaffari   // -- Define derived units
171a515125bSLeila Ghaffari   Pascal          = kilogram / (meter * PetscSqr(second));
172a515125bSLeila Ghaffari   J_per_kg_K      =  PetscSqr(meter) / (PetscSqr(second) * Kelvin);
173a515125bSLeila Ghaffari   m_per_squared_s = meter / PetscSqr(second);
174a515125bSLeila Ghaffari   W_per_m_K       = kilogram * meter / (pow(second,3) * Kelvin);
175a515125bSLeila Ghaffari 
176a515125bSLeila Ghaffari   user->units->meter           = meter;
177a515125bSLeila Ghaffari   user->units->kilogram        = kilogram;
178a515125bSLeila Ghaffari   user->units->second          = second;
179a515125bSLeila Ghaffari   user->units->Kelvin          = Kelvin;
180a515125bSLeila Ghaffari   user->units->Pascal          = Pascal;
181a515125bSLeila Ghaffari   user->units->J_per_kg_K      = J_per_kg_K;
182a515125bSLeila Ghaffari   user->units->m_per_squared_s = m_per_squared_s;
183a515125bSLeila Ghaffari   user->units->W_per_m_K       = W_per_m_K;
184a515125bSLeila Ghaffari 
185a515125bSLeila Ghaffari   // ------------------------------------------------------
186a515125bSLeila Ghaffari   //           Set up the libCEED context
187a515125bSLeila Ghaffari   // ------------------------------------------------------
188a515125bSLeila Ghaffari   // -- Scale variables to desired units
189a515125bSLeila Ghaffari   theta0 *= Kelvin;
190a515125bSLeila Ghaffari   thetaC *= Kelvin;
191a515125bSLeila Ghaffari   P0     *= Pascal;
192a515125bSLeila Ghaffari   N      *= (1./second);
193a515125bSLeila Ghaffari   cv     *= J_per_kg_K;
194a515125bSLeila Ghaffari   cp     *= J_per_kg_K;
195a515125bSLeila Ghaffari   g      *= m_per_squared_s;
196a515125bSLeila Ghaffari   mu     *= Pascal * second;
197a515125bSLeila Ghaffari   k      *= W_per_m_K;
198a515125bSLeila Ghaffari   rc     = fabs(rc) * meter;
199*05a512bdSLeila Ghaffari   for (int i=0; i<3; i++) domain_size[i] *= meter;
200a515125bSLeila Ghaffari   for (int i=0; i<3; i++) center[i] *= meter;
201*05a512bdSLeila Ghaffari   problem->dm_scale = meter;
202a515125bSLeila Ghaffari 
203a515125bSLeila Ghaffari   // -- Setup Context
204a515125bSLeila Ghaffari   setup_context->theta0     = theta0;
205a515125bSLeila Ghaffari   setup_context->thetaC     = thetaC;
206a515125bSLeila Ghaffari   setup_context->P0         = P0;
207a515125bSLeila Ghaffari   setup_context->N          = N;
208a515125bSLeila Ghaffari   setup_context->cv         = cv;
209a515125bSLeila Ghaffari   setup_context->cp         = cp;
210a515125bSLeila Ghaffari   setup_context->g          = g;
211a515125bSLeila Ghaffari   setup_context->rc         = rc;
212*05a512bdSLeila Ghaffari   setup_context->lx         = domain_size[0];
213*05a512bdSLeila Ghaffari   setup_context->ly         = domain_size[1];
214*05a512bdSLeila Ghaffari   setup_context->lz         = domain_size[2];
215a515125bSLeila Ghaffari   setup_context->center[0]  = center[0];
216a515125bSLeila Ghaffari   setup_context->center[1]  = center[1];
217a515125bSLeila Ghaffari   setup_context->center[2]  = center[2];
218a515125bSLeila Ghaffari   setup_context->dc_axis[0] = dc_axis[0];
219a515125bSLeila Ghaffari   setup_context->dc_axis[1] = dc_axis[1];
220a515125bSLeila Ghaffari   setup_context->dc_axis[2] = dc_axis[2];
221a515125bSLeila Ghaffari   setup_context->time       = 0;
222a515125bSLeila Ghaffari 
223a515125bSLeila Ghaffari   // -- QFunction Context
224a515125bSLeila Ghaffari   user->phys->stab             = stab;
225a515125bSLeila Ghaffari   user->phys->implicit         = implicit;
226a515125bSLeila Ghaffari   user->phys->has_curr_time    = has_curr_time;
227a515125bSLeila Ghaffari   user->phys->dc_ctx->lambda   = lambda;
228a515125bSLeila Ghaffari   user->phys->dc_ctx->mu       = mu;
229a515125bSLeila Ghaffari   user->phys->dc_ctx->k        = k;
230a515125bSLeila Ghaffari   user->phys->dc_ctx->cv       = cv;
231a515125bSLeila Ghaffari   user->phys->dc_ctx->cp       = cp;
232a515125bSLeila Ghaffari   user->phys->dc_ctx->g        = g;
233d8a22b9eSJed Brown   user->phys->dc_ctx->c_tau    = c_tau;
234a515125bSLeila Ghaffari   user->phys->dc_ctx->stabilization = stab;
235a515125bSLeila Ghaffari 
236a515125bSLeila Ghaffari   PetscFunctionReturn(0);
237a515125bSLeila Ghaffari }
238a515125bSLeila Ghaffari 
239ba5420e5SLeila Ghaffari PetscErrorCode SetupContext_DENSITY_CURRENT(Ceed ceed, CeedData ceed_data,
240ba5420e5SLeila Ghaffari     AppCtx app_ctx, SetupContext setup_ctx,
241ba5420e5SLeila Ghaffari     Physics phys) {
242ba5420e5SLeila Ghaffari   PetscFunctionBeginUser;
243ba5420e5SLeila Ghaffari 
244ba5420e5SLeila Ghaffari   CeedQFunctionContextCreate(ceed, &ceed_data->setup_context);
245ba5420e5SLeila Ghaffari   CeedQFunctionContextSetData(ceed_data->setup_context, CEED_MEM_HOST,
246ba5420e5SLeila Ghaffari                               CEED_USE_POINTER, sizeof(*setup_ctx), setup_ctx);
247ba5420e5SLeila Ghaffari   CeedQFunctionSetContext(ceed_data->qf_ics, ceed_data->setup_context);
248ba5420e5SLeila Ghaffari   CeedQFunctionContextCreate(ceed, &ceed_data->dc_context);
249ba5420e5SLeila Ghaffari   CeedQFunctionContextSetData(ceed_data->dc_context, CEED_MEM_HOST,
250ba5420e5SLeila Ghaffari                               CEED_USE_POINTER,
251ba5420e5SLeila Ghaffari                               sizeof(*phys->dc_ctx), phys->dc_ctx);
252ba5420e5SLeila Ghaffari   if (ceed_data->qf_rhs_vol)
253ba5420e5SLeila Ghaffari     CeedQFunctionSetContext(ceed_data->qf_rhs_vol, ceed_data->dc_context);
254ba5420e5SLeila Ghaffari   if (ceed_data->qf_ifunction_vol)
255ba5420e5SLeila Ghaffari     CeedQFunctionSetContext(ceed_data->qf_ifunction_vol, ceed_data->dc_context);
256ba5420e5SLeila Ghaffari 
257ba5420e5SLeila Ghaffari   PetscFunctionReturn(0);
258ba5420e5SLeila Ghaffari }
259ba5420e5SLeila Ghaffari 
260a515125bSLeila Ghaffari PetscErrorCode BC_DENSITY_CURRENT(DM dm, SimpleBC bc, Physics phys,
261a515125bSLeila Ghaffari                                   void *setup_ctx) {
262a515125bSLeila Ghaffari 
263a515125bSLeila Ghaffari   PetscInt       len;
264a515125bSLeila Ghaffari   PetscBool      flg;
265a515125bSLeila Ghaffari   MPI_Comm       comm = PETSC_COMM_WORLD;
266a515125bSLeila Ghaffari   PetscErrorCode ierr;
267a515125bSLeila Ghaffari   PetscFunctionBeginUser;
268a515125bSLeila Ghaffari 
269a515125bSLeila Ghaffari   // Default boundary conditions
270a515125bSLeila Ghaffari   //   slip bc on all faces and no wall bc
271a515125bSLeila Ghaffari   bc->num_slip[0] = bc->num_slip[1] = bc->num_slip[2] = 2;
272a515125bSLeila Ghaffari   bc->slips[0][0] = 5;
273a515125bSLeila Ghaffari   bc->slips[0][1] = 6;
274a515125bSLeila Ghaffari   bc->slips[1][0] = 3;
275a515125bSLeila Ghaffari   bc->slips[1][1] = 4;
276a515125bSLeila Ghaffari   bc->slips[2][0] = 1;
277a515125bSLeila Ghaffari   bc->slips[2][1] = 2;
278a515125bSLeila Ghaffari 
279a515125bSLeila Ghaffari   // Parse command line options
280a515125bSLeila Ghaffari   ierr = PetscOptionsBegin(comm, NULL, "Options for DENSITY_CURRENT BCs ",
281a515125bSLeila Ghaffari                            NULL); CHKERRQ(ierr);
282a515125bSLeila Ghaffari   ierr = PetscOptionsIntArray("-bc_wall",
283a515125bSLeila Ghaffari                               "Use wall boundary conditions on this list of faces",
284a515125bSLeila Ghaffari                               NULL, bc->walls,
285a515125bSLeila Ghaffari                               (len = sizeof(bc->walls) / sizeof(bc->walls[0]),
286a515125bSLeila Ghaffari                                &len), &flg); CHKERRQ(ierr);
287a515125bSLeila Ghaffari   if (flg) {
288a515125bSLeila Ghaffari     bc->num_wall = len;
289a515125bSLeila Ghaffari     // Using a no-slip wall disables automatic slip walls (they must be set explicitly)
290a515125bSLeila Ghaffari     bc->num_slip[0] = bc->num_slip[1] = bc->num_slip[2] = 0;
291a515125bSLeila Ghaffari   }
292a515125bSLeila Ghaffari   for (PetscInt j=0; j<3; j++) {
293a515125bSLeila Ghaffari     const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"};
294a515125bSLeila Ghaffari     ierr = PetscOptionsIntArray(flags[j],
295a515125bSLeila Ghaffari                                 "Use slip boundary conditions on this list of faces",
296a515125bSLeila Ghaffari                                 NULL, bc->slips[j],
297a515125bSLeila Ghaffari                                 (len = sizeof(bc->slips[j]) / sizeof(bc->slips[j][0]),
298a515125bSLeila Ghaffari                                  &len), &flg); CHKERRQ(ierr);
299a515125bSLeila Ghaffari     if (flg) {
300a515125bSLeila Ghaffari       bc->num_slip[j] = len;
301a515125bSLeila Ghaffari       bc->user_bc = PETSC_TRUE;
302a515125bSLeila Ghaffari     }
303a515125bSLeila Ghaffari   }
304a515125bSLeila Ghaffari   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
305a515125bSLeila Ghaffari 
306a515125bSLeila Ghaffari   {
307a515125bSLeila Ghaffari     // Set slip boundary conditions
308a515125bSLeila Ghaffari     DMLabel label;
309a515125bSLeila Ghaffari     ierr = DMGetLabel(dm, "Face Sets", &label); CHKERRQ(ierr);
310a515125bSLeila Ghaffari     PetscInt comps[1] = {1};
311b9842f74SJeremy L Thompson     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", label,
312a515125bSLeila Ghaffari                          bc->num_slip[0], bc->slips[0], 0, 1, comps,
313a515125bSLeila Ghaffari                          (void(*)(void))NULL, NULL, setup_ctx, NULL);
314a515125bSLeila Ghaffari     CHKERRQ(ierr);
315a515125bSLeila Ghaffari     comps[0] = 2;
316b9842f74SJeremy L Thompson     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", label,
317a515125bSLeila Ghaffari                          bc->num_slip[1], bc->slips[1], 0, 1, comps,
318a515125bSLeila Ghaffari                          (void(*)(void))NULL, NULL, setup_ctx, NULL);
319a515125bSLeila Ghaffari     CHKERRQ(ierr);
320a515125bSLeila Ghaffari     comps[0] = 3;
321b9842f74SJeremy L Thompson     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", label,
322a515125bSLeila Ghaffari                          bc->num_slip[2], bc->slips[2], 0, 1, comps,
323a515125bSLeila Ghaffari                          (void(*)(void))NULL, NULL, setup_ctx, NULL);
324a515125bSLeila Ghaffari     CHKERRQ(ierr);
325a515125bSLeila Ghaffari   }
326a515125bSLeila Ghaffari 
327139613f2SLeila Ghaffari   if (bc->user_bc) {
328a515125bSLeila Ghaffari     for (PetscInt c = 0; c < 3; c++) {
329a515125bSLeila Ghaffari       for (PetscInt s = 0; s < bc->num_slip[c]; s++) {
330a515125bSLeila Ghaffari         for (PetscInt w = 0; w < bc->num_wall; w++) {
331a515125bSLeila Ghaffari           if (bc->slips[c][s] == bc->walls[w])
332a515125bSLeila Ghaffari             SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG,
333a515125bSLeila Ghaffari                      "Boundary condition already set on face %D!\n",
334a515125bSLeila Ghaffari                      bc->walls[w]);
335a515125bSLeila Ghaffari         }
336a515125bSLeila Ghaffari       }
337a515125bSLeila Ghaffari     }
338a515125bSLeila Ghaffari   }
339a515125bSLeila Ghaffari 
340a515125bSLeila Ghaffari   // Set wall boundary conditions
341a515125bSLeila Ghaffari   //   zero velocity and zero flux for mass density and energy density
342a515125bSLeila Ghaffari   {
343a515125bSLeila Ghaffari     DMLabel  label;
344a515125bSLeila Ghaffari     PetscInt comps[3] = {1, 2, 3};
345a515125bSLeila Ghaffari     ierr = DMGetLabel(dm, "Face Sets", &label); CHKERRQ(ierr);
346b9842f74SJeremy L Thompson     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label,
347b9842f74SJeremy L Thompson                          bc->num_wall, bc->walls, 0, 3, comps,
348b9842f74SJeremy L Thompson                          (void(*)(void))Exact_DC, NULL,
349a515125bSLeila Ghaffari                          setup_ctx, NULL); CHKERRQ(ierr);
350a515125bSLeila Ghaffari   }
351a515125bSLeila Ghaffari 
352a515125bSLeila Ghaffari   PetscFunctionReturn(0);
353a515125bSLeila Ghaffari }
354a515125bSLeila Ghaffari 
355a515125bSLeila Ghaffari PetscErrorCode PRINT_DENSITY_CURRENT(Physics phys, SetupContext setup_ctx,
356a515125bSLeila Ghaffari                                      AppCtx app_ctx) {
357a515125bSLeila Ghaffari   MPI_Comm       comm = PETSC_COMM_WORLD;
358a515125bSLeila Ghaffari   PetscErrorCode ierr;
359a515125bSLeila Ghaffari   PetscFunctionBeginUser;
360a515125bSLeila Ghaffari 
361a515125bSLeila Ghaffari   ierr = PetscPrintf(comm,
362a515125bSLeila Ghaffari                      "  Problem:\n"
363a515125bSLeila Ghaffari                      "    Problem Name                       : %s\n"
364a515125bSLeila Ghaffari                      "    Stabilization                      : %s\n",
365a515125bSLeila Ghaffari                      app_ctx->problem_name, StabilizationTypes[phys->stab]);
366a515125bSLeila Ghaffari   CHKERRQ(ierr);
367a515125bSLeila Ghaffari 
368a515125bSLeila Ghaffari   PetscFunctionReturn(0);
369a515125bSLeila Ghaffari }
370