xref: /libCEED/examples/fluids/problems/densitycurrent.c (revision 1864f1c2b4e770a2a9adc26a02ef77fc3a284256)
177841947SLeila Ghaffari // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
277841947SLeila Ghaffari // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
377841947SLeila Ghaffari // reserved. See files LICENSE and NOTICE for details.
477841947SLeila Ghaffari //
577841947SLeila Ghaffari // This file is part of CEED, a collection of benchmarks, miniapps, software
677841947SLeila Ghaffari // libraries and APIs for efficient high-order finite element and spectral
777841947SLeila Ghaffari // element discretizations for exascale applications. For more information and
877841947SLeila Ghaffari // source code availability see http://github.com/ceed.
977841947SLeila Ghaffari //
1077841947SLeila Ghaffari // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
1177841947SLeila Ghaffari // a collaborative effort of two U.S. Department of Energy organizations (Office
1277841947SLeila Ghaffari // of Science and the National Nuclear Security Administration) responsible for
1377841947SLeila Ghaffari // the planning and preparation of a capable exascale ecosystem, including
1477841947SLeila Ghaffari // software, applications, hardware, advanced system engineering and early
1577841947SLeila Ghaffari // testbed platforms, in support of the nation's exascale computing imperative.
1677841947SLeila Ghaffari 
1777841947SLeila Ghaffari /// @file
1877841947SLeila Ghaffari /// Utility functions for setting up DENSITY_CURRENT
1977841947SLeila Ghaffari 
2077841947SLeila Ghaffari #include "../navierstokes.h"
2177841947SLeila Ghaffari #include "../qfunctions/setupgeo.h"
2277841947SLeila Ghaffari #include "../qfunctions/densitycurrent.h"
2377841947SLeila Ghaffari 
24*1864f1c2SLeila Ghaffari PetscErrorCode NS_DENSITY_CURRENT(ProblemData *problem, DM dm, void *setup_ctx,
2577841947SLeila Ghaffari                                   void *ctx) {
2677841947SLeila Ghaffari   SetupContext      setup_context = *(SetupContext *)setup_ctx;
2777841947SLeila Ghaffari   User              user = *(User *)ctx;
2877841947SLeila Ghaffari   StabilizationType stab;
2977841947SLeila Ghaffari   MPI_Comm          comm = PETSC_COMM_WORLD;
3077841947SLeila Ghaffari   PetscBool         implicit;
3177841947SLeila Ghaffari   PetscBool         has_curr_time = PETSC_FALSE;
3277841947SLeila Ghaffari   PetscInt          ierr;
3377841947SLeila Ghaffari   PetscFunctionBeginUser;
3477841947SLeila Ghaffari 
3577841947SLeila Ghaffari   ierr = PetscCalloc1(1, &user->phys->dc_ctx); CHKERRQ(ierr);
3677841947SLeila Ghaffari 
3777841947SLeila Ghaffari   // ------------------------------------------------------
3877841947SLeila Ghaffari   //               SET UP DENSITY_CURRENT
3977841947SLeila Ghaffari   // ------------------------------------------------------
4077841947SLeila Ghaffari   problem->dim                     = 3;
4177841947SLeila Ghaffari   problem->q_data_size_vol         = 10;
4277841947SLeila Ghaffari   problem->q_data_size_sur         = 4;
4377841947SLeila Ghaffari   problem->setup_vol               = Setup;
4477841947SLeila Ghaffari   problem->setup_vol_loc           = Setup_loc;
4577841947SLeila Ghaffari   problem->setup_sur               = SetupBoundary;
4677841947SLeila Ghaffari   problem->setup_sur_loc           = SetupBoundary_loc;
4777841947SLeila Ghaffari   problem->ics                     = ICsDC;
4877841947SLeila Ghaffari   problem->ics_loc                 = ICsDC_loc;
4977841947SLeila Ghaffari   problem->apply_vol_rhs           = DC;
5077841947SLeila Ghaffari   problem->apply_vol_rhs_loc       = DC_loc;
5177841947SLeila Ghaffari   problem->apply_vol_ifunction     = IFunction_DC;
5277841947SLeila Ghaffari   problem->apply_vol_ifunction_loc = IFunction_DC_loc;
5377841947SLeila Ghaffari   problem->bc                      = Exact_DC;
54d0b732dbSLeila Ghaffari   problem->setup_ctx               = SetupContext_DENSITY_CURRENT;
5577841947SLeila Ghaffari   problem->bc_func                 = BC_DENSITY_CURRENT;
5677841947SLeila Ghaffari   problem->non_zero_time           = PETSC_FALSE;
5777841947SLeila Ghaffari   problem->print_info              = PRINT_DENSITY_CURRENT;
5877841947SLeila Ghaffari 
5977841947SLeila Ghaffari   // ------------------------------------------------------
6077841947SLeila Ghaffari   //             Create the libCEED context
6177841947SLeila Ghaffari   // ------------------------------------------------------
6277841947SLeila Ghaffari   CeedScalar theta0 = 300.;    // K
6377841947SLeila Ghaffari   CeedScalar thetaC = -15.;    // K
6477841947SLeila Ghaffari   CeedScalar P0     = 1.e5;    // Pa
6577841947SLeila Ghaffari   CeedScalar N      = 0.01;    // 1/s
6677841947SLeila Ghaffari   CeedScalar cv     = 717.;    // J/(kg K)
6777841947SLeila Ghaffari   CeedScalar cp     = 1004.;   // J/(kg K)
6877841947SLeila Ghaffari   CeedScalar g      = 9.81;    // m/s^2
6977841947SLeila Ghaffari   CeedScalar lambda = -2./3.;  // -
7077841947SLeila Ghaffari   CeedScalar mu     = 75.;     // Pa s, dynamic viscosity
7177841947SLeila Ghaffari   // mu = 75 is not physical for air, but is good for numerical stability
7277841947SLeila Ghaffari   CeedScalar k      = 0.02638; // W/(m K)
73504dc8e0SLeila Ghaffari   CeedScalar c_tau  = 0.5;     // -
74932417b3SJed Brown   // c_tau = 0.5 is reported as "optimal" in Hughes et al 2010
7577841947SLeila Ghaffari   CeedScalar rc     = 1000.;   // m (Radius of bubble)
7677841947SLeila Ghaffari   PetscReal center[3], dc_axis[3] = {0, 0, 0};
77*1864f1c2SLeila Ghaffari   PetscReal domain_min[3], domain_max[3], domain_size[3];
78*1864f1c2SLeila Ghaffari   ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
79*1864f1c2SLeila Ghaffari   for (int i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i];
8077841947SLeila Ghaffari 
8177841947SLeila Ghaffari   // ------------------------------------------------------
8277841947SLeila Ghaffari   //             Create the PETSc context
8377841947SLeila Ghaffari   // ------------------------------------------------------
8477841947SLeila Ghaffari   PetscScalar meter    = 1e-2;  // 1 meter in scaled length units
8577841947SLeila Ghaffari   PetscScalar kilogram = 1e-6;  // 1 kilogram in scaled mass units
8677841947SLeila Ghaffari   PetscScalar second   = 1e-2;  // 1 second in scaled time units
8777841947SLeila Ghaffari   PetscScalar Kelvin   = 1;     // 1 Kelvin in scaled temperature units
8877841947SLeila Ghaffari   PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s;
8977841947SLeila Ghaffari 
9077841947SLeila Ghaffari   // ------------------------------------------------------
9177841947SLeila Ghaffari   //              Command line Options
9277841947SLeila Ghaffari   // ------------------------------------------------------
9377841947SLeila Ghaffari   ierr = PetscOptionsBegin(comm, NULL, "Options for DENSITY_CURRENT problem",
9477841947SLeila Ghaffari                            NULL); CHKERRQ(ierr);
9577841947SLeila Ghaffari   // -- Physics
9677841947SLeila Ghaffari   ierr = PetscOptionsScalar("-theta0", "Reference potential temperature",
9777841947SLeila Ghaffari                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
9877841947SLeila Ghaffari   ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature",
9977841947SLeila Ghaffari                             NULL, thetaC, &thetaC, NULL); CHKERRQ(ierr);
10077841947SLeila Ghaffari   ierr = PetscOptionsScalar("-P0", "Atmospheric pressure",
10177841947SLeila Ghaffari                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
10277841947SLeila Ghaffari   ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency",
10377841947SLeila Ghaffari                             NULL, N, &N, NULL); CHKERRQ(ierr);
10477841947SLeila Ghaffari   ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume",
10577841947SLeila Ghaffari                             NULL, cv, &cv, NULL); CHKERRQ(ierr);
10677841947SLeila Ghaffari   ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure",
10777841947SLeila Ghaffari                             NULL, cp, &cp, NULL); CHKERRQ(ierr);
10877841947SLeila Ghaffari   ierr = PetscOptionsScalar("-g", "Gravitational acceleration",
10977841947SLeila Ghaffari                             NULL, g, &g, NULL); CHKERRQ(ierr);
11077841947SLeila Ghaffari   ierr = PetscOptionsScalar("-lambda",
11177841947SLeila Ghaffari                             "Stokes hypothesis second viscosity coefficient",
11277841947SLeila Ghaffari                             NULL, lambda, &lambda, NULL); CHKERRQ(ierr);
11377841947SLeila Ghaffari   ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient",
11477841947SLeila Ghaffari                             NULL, mu, &mu, NULL); CHKERRQ(ierr);
11577841947SLeila Ghaffari   ierr = PetscOptionsScalar("-k", "Thermal conductivity",
11677841947SLeila Ghaffari                             NULL, k, &k, NULL); CHKERRQ(ierr);
11777841947SLeila Ghaffari   ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble",
11877841947SLeila Ghaffari                             NULL, rc, &rc, NULL); CHKERRQ(ierr);
119*1864f1c2SLeila Ghaffari   for (int i=0; i<3; i++) center[i] = .5*domain_size[i];
12077841947SLeila Ghaffari   PetscInt n = problem->dim;
12177841947SLeila Ghaffari   ierr = PetscOptionsRealArray("-center", "Location of bubble center",
12277841947SLeila Ghaffari                                NULL, center, &n, NULL); CHKERRQ(ierr);
12377841947SLeila Ghaffari   n = problem->dim;
12477841947SLeila Ghaffari   ierr = PetscOptionsRealArray("-dc_axis",
12577841947SLeila Ghaffari                                "Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric",
12677841947SLeila Ghaffari                                NULL, dc_axis, &n, NULL); CHKERRQ(ierr);
12777841947SLeila Ghaffari   {
12877841947SLeila Ghaffari     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) + PetscSqr(dc_axis[1]) +
12977841947SLeila Ghaffari                                    PetscSqr(dc_axis[2]));
13077841947SLeila Ghaffari     if (norm > 0) {
13177841947SLeila Ghaffari       for (int i=0; i<3; i++)  dc_axis[i] /= norm;
13277841947SLeila Ghaffari     }
13377841947SLeila Ghaffari   }
13477841947SLeila Ghaffari   ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL,
13577841947SLeila Ghaffari                           StabilizationTypes, (PetscEnum)(stab = STAB_NONE),
13677841947SLeila Ghaffari                           (PetscEnum *)&stab, NULL); CHKERRQ(ierr);
137932417b3SJed Brown   ierr = PetscOptionsScalar("-c_tau", "Stabilization constant",
138932417b3SJed Brown                             NULL, c_tau, &c_tau, NULL); CHKERRQ(ierr);
13977841947SLeila Ghaffari   ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation",
14077841947SLeila Ghaffari                           NULL, implicit=PETSC_FALSE, &implicit, NULL);
14177841947SLeila Ghaffari   CHKERRQ(ierr);
14277841947SLeila Ghaffari 
14377841947SLeila Ghaffari   // -- Units
14477841947SLeila Ghaffari   ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units",
14577841947SLeila Ghaffari                             NULL, meter, &meter, NULL); CHKERRQ(ierr);
14677841947SLeila Ghaffari   meter = fabs(meter);
14777841947SLeila Ghaffari   ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units",
14877841947SLeila Ghaffari                             NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr);
14977841947SLeila Ghaffari   kilogram = fabs(kilogram);
15077841947SLeila Ghaffari   ierr = PetscOptionsScalar("-units_second","1 second in scaled time units",
15177841947SLeila Ghaffari                             NULL, second, &second, NULL); CHKERRQ(ierr);
15277841947SLeila Ghaffari   second = fabs(second);
15377841947SLeila Ghaffari   ierr = PetscOptionsScalar("-units_Kelvin",
15477841947SLeila Ghaffari                             "1 Kelvin in scaled temperature units",
15577841947SLeila Ghaffari                             NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr);
15677841947SLeila Ghaffari   Kelvin = fabs(Kelvin);
15777841947SLeila Ghaffari 
15877841947SLeila Ghaffari   // -- Warnings
15977841947SLeila Ghaffari   if (stab == STAB_SUPG && !implicit) {
16077841947SLeila Ghaffari     ierr = PetscPrintf(comm,
16177841947SLeila Ghaffari                        "Warning! Use -stab supg only with -implicit\n");
16277841947SLeila Ghaffari     CHKERRQ(ierr);
16377841947SLeila Ghaffari   }
16477841947SLeila Ghaffari 
16577841947SLeila Ghaffari   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
16677841947SLeila Ghaffari 
16777841947SLeila Ghaffari   // ------------------------------------------------------
16877841947SLeila Ghaffari   //           Set up the PETSc context
16977841947SLeila Ghaffari   // ------------------------------------------------------
17077841947SLeila Ghaffari   // -- Define derived units
17177841947SLeila Ghaffari   Pascal          = kilogram / (meter * PetscSqr(second));
17277841947SLeila Ghaffari   J_per_kg_K      =  PetscSqr(meter) / (PetscSqr(second) * Kelvin);
17377841947SLeila Ghaffari   m_per_squared_s = meter / PetscSqr(second);
17477841947SLeila Ghaffari   W_per_m_K       = kilogram * meter / (pow(second,3) * Kelvin);
17577841947SLeila Ghaffari 
17677841947SLeila Ghaffari   user->units->meter           = meter;
17777841947SLeila Ghaffari   user->units->kilogram        = kilogram;
17877841947SLeila Ghaffari   user->units->second          = second;
17977841947SLeila Ghaffari   user->units->Kelvin          = Kelvin;
18077841947SLeila Ghaffari   user->units->Pascal          = Pascal;
18177841947SLeila Ghaffari   user->units->J_per_kg_K      = J_per_kg_K;
18277841947SLeila Ghaffari   user->units->m_per_squared_s = m_per_squared_s;
18377841947SLeila Ghaffari   user->units->W_per_m_K       = W_per_m_K;
18477841947SLeila Ghaffari 
18577841947SLeila Ghaffari   // ------------------------------------------------------
18677841947SLeila Ghaffari   //           Set up the libCEED context
18777841947SLeila Ghaffari   // ------------------------------------------------------
18877841947SLeila Ghaffari   // -- Scale variables to desired units
18977841947SLeila Ghaffari   theta0 *= Kelvin;
19077841947SLeila Ghaffari   thetaC *= Kelvin;
19177841947SLeila Ghaffari   P0     *= Pascal;
19277841947SLeila Ghaffari   N      *= (1./second);
19377841947SLeila Ghaffari   cv     *= J_per_kg_K;
19477841947SLeila Ghaffari   cp     *= J_per_kg_K;
19577841947SLeila Ghaffari   g      *= m_per_squared_s;
19677841947SLeila Ghaffari   mu     *= Pascal * second;
19777841947SLeila Ghaffari   k      *= W_per_m_K;
19877841947SLeila Ghaffari   rc     = fabs(rc) * meter;
199*1864f1c2SLeila Ghaffari   for (int i=0; i<3; i++) domain_size[i] *= meter;
20077841947SLeila Ghaffari   for (int i=0; i<3; i++) center[i] *= meter;
201*1864f1c2SLeila Ghaffari   problem->dm_scale = meter;
20277841947SLeila Ghaffari 
20377841947SLeila Ghaffari   // -- Setup Context
20477841947SLeila Ghaffari   setup_context->theta0     = theta0;
20577841947SLeila Ghaffari   setup_context->thetaC     = thetaC;
20677841947SLeila Ghaffari   setup_context->P0         = P0;
20777841947SLeila Ghaffari   setup_context->N          = N;
20877841947SLeila Ghaffari   setup_context->cv         = cv;
20977841947SLeila Ghaffari   setup_context->cp         = cp;
21077841947SLeila Ghaffari   setup_context->g          = g;
21177841947SLeila Ghaffari   setup_context->rc         = rc;
212*1864f1c2SLeila Ghaffari   setup_context->lx         = domain_size[0];
213*1864f1c2SLeila Ghaffari   setup_context->ly         = domain_size[1];
214*1864f1c2SLeila Ghaffari   setup_context->lz         = domain_size[2];
21577841947SLeila Ghaffari   setup_context->center[0]  = center[0];
21677841947SLeila Ghaffari   setup_context->center[1]  = center[1];
21777841947SLeila Ghaffari   setup_context->center[2]  = center[2];
21877841947SLeila Ghaffari   setup_context->dc_axis[0] = dc_axis[0];
21977841947SLeila Ghaffari   setup_context->dc_axis[1] = dc_axis[1];
22077841947SLeila Ghaffari   setup_context->dc_axis[2] = dc_axis[2];
22177841947SLeila Ghaffari   setup_context->time       = 0;
22277841947SLeila Ghaffari 
22377841947SLeila Ghaffari   // -- QFunction Context
22477841947SLeila Ghaffari   user->phys->stab             = stab;
22577841947SLeila Ghaffari   user->phys->implicit         = implicit;
22677841947SLeila Ghaffari   user->phys->has_curr_time    = has_curr_time;
22777841947SLeila Ghaffari   user->phys->dc_ctx->lambda   = lambda;
22877841947SLeila Ghaffari   user->phys->dc_ctx->mu       = mu;
22977841947SLeila Ghaffari   user->phys->dc_ctx->k        = k;
23077841947SLeila Ghaffari   user->phys->dc_ctx->cv       = cv;
23177841947SLeila Ghaffari   user->phys->dc_ctx->cp       = cp;
23277841947SLeila Ghaffari   user->phys->dc_ctx->g        = g;
233932417b3SJed Brown   user->phys->dc_ctx->c_tau    = c_tau;
23477841947SLeila Ghaffari   user->phys->dc_ctx->stabilization = stab;
23577841947SLeila Ghaffari 
23677841947SLeila Ghaffari   PetscFunctionReturn(0);
23777841947SLeila Ghaffari }
23877841947SLeila Ghaffari 
239d0b732dbSLeila Ghaffari PetscErrorCode SetupContext_DENSITY_CURRENT(Ceed ceed, CeedData ceed_data,
240d0b732dbSLeila Ghaffari     AppCtx app_ctx, SetupContext setup_ctx,
241d0b732dbSLeila Ghaffari     Physics phys) {
242d0b732dbSLeila Ghaffari   PetscFunctionBeginUser;
243d0b732dbSLeila Ghaffari 
244d0b732dbSLeila Ghaffari   CeedQFunctionContextCreate(ceed, &ceed_data->setup_context);
245d0b732dbSLeila Ghaffari   CeedQFunctionContextSetData(ceed_data->setup_context, CEED_MEM_HOST,
246d0b732dbSLeila Ghaffari                               CEED_USE_POINTER, sizeof(*setup_ctx), setup_ctx);
247d0b732dbSLeila Ghaffari   CeedQFunctionSetContext(ceed_data->qf_ics, ceed_data->setup_context);
248d0b732dbSLeila Ghaffari   CeedQFunctionContextCreate(ceed, &ceed_data->dc_context);
249d0b732dbSLeila Ghaffari   CeedQFunctionContextSetData(ceed_data->dc_context, CEED_MEM_HOST,
250d0b732dbSLeila Ghaffari                               CEED_USE_POINTER,
251d0b732dbSLeila Ghaffari                               sizeof(*phys->dc_ctx), phys->dc_ctx);
252d0b732dbSLeila Ghaffari   if (ceed_data->qf_rhs_vol)
253d0b732dbSLeila Ghaffari     CeedQFunctionSetContext(ceed_data->qf_rhs_vol, ceed_data->dc_context);
254d0b732dbSLeila Ghaffari   if (ceed_data->qf_ifunction_vol)
255d0b732dbSLeila Ghaffari     CeedQFunctionSetContext(ceed_data->qf_ifunction_vol, ceed_data->dc_context);
256d0b732dbSLeila Ghaffari 
257d0b732dbSLeila Ghaffari   PetscFunctionReturn(0);
258d0b732dbSLeila Ghaffari }
259d0b732dbSLeila Ghaffari 
26077841947SLeila Ghaffari PetscErrorCode BC_DENSITY_CURRENT(DM dm, SimpleBC bc, Physics phys,
26177841947SLeila Ghaffari                                   void *setup_ctx) {
26277841947SLeila Ghaffari 
26377841947SLeila Ghaffari   PetscInt       len;
26477841947SLeila Ghaffari   PetscBool      flg;
26577841947SLeila Ghaffari   MPI_Comm       comm = PETSC_COMM_WORLD;
26677841947SLeila Ghaffari   PetscErrorCode ierr;
26777841947SLeila Ghaffari   PetscFunctionBeginUser;
26877841947SLeila Ghaffari 
26977841947SLeila Ghaffari   // Default boundary conditions
27077841947SLeila Ghaffari   //   slip bc on all faces and no wall bc
27177841947SLeila Ghaffari   bc->num_slip[0] = bc->num_slip[1] = bc->num_slip[2] = 2;
27277841947SLeila Ghaffari   bc->slips[0][0] = 5;
27377841947SLeila Ghaffari   bc->slips[0][1] = 6;
27477841947SLeila Ghaffari   bc->slips[1][0] = 3;
27577841947SLeila Ghaffari   bc->slips[1][1] = 4;
27677841947SLeila Ghaffari   bc->slips[2][0] = 1;
27777841947SLeila Ghaffari   bc->slips[2][1] = 2;
27877841947SLeila Ghaffari 
27977841947SLeila Ghaffari   // Parse command line options
28077841947SLeila Ghaffari   ierr = PetscOptionsBegin(comm, NULL, "Options for DENSITY_CURRENT BCs ",
28177841947SLeila Ghaffari                            NULL); CHKERRQ(ierr);
28277841947SLeila Ghaffari   ierr = PetscOptionsIntArray("-bc_wall",
28377841947SLeila Ghaffari                               "Use wall boundary conditions on this list of faces",
28477841947SLeila Ghaffari                               NULL, bc->walls,
28577841947SLeila Ghaffari                               (len = sizeof(bc->walls) / sizeof(bc->walls[0]),
28677841947SLeila Ghaffari                                &len), &flg); CHKERRQ(ierr);
28777841947SLeila Ghaffari   if (flg) {
28877841947SLeila Ghaffari     bc->num_wall = len;
28977841947SLeila Ghaffari     // Using a no-slip wall disables automatic slip walls (they must be set explicitly)
29077841947SLeila Ghaffari     bc->num_slip[0] = bc->num_slip[1] = bc->num_slip[2] = 0;
29177841947SLeila Ghaffari   }
29277841947SLeila Ghaffari   for (PetscInt j=0; j<3; j++) {
29377841947SLeila Ghaffari     const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"};
29477841947SLeila Ghaffari     ierr = PetscOptionsIntArray(flags[j],
29577841947SLeila Ghaffari                                 "Use slip boundary conditions on this list of faces",
29677841947SLeila Ghaffari                                 NULL, bc->slips[j],
29777841947SLeila Ghaffari                                 (len = sizeof(bc->slips[j]) / sizeof(bc->slips[j][0]),
29877841947SLeila Ghaffari                                  &len), &flg); CHKERRQ(ierr);
29977841947SLeila Ghaffari     if (flg) {
30077841947SLeila Ghaffari       bc->num_slip[j] = len;
30177841947SLeila Ghaffari       bc->user_bc = PETSC_TRUE;
30277841947SLeila Ghaffari     }
30377841947SLeila Ghaffari   }
30477841947SLeila Ghaffari   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
30577841947SLeila Ghaffari 
30677841947SLeila Ghaffari   {
30777841947SLeila Ghaffari     // Set slip boundary conditions
30877841947SLeila Ghaffari     DMLabel label;
30977841947SLeila Ghaffari     ierr = DMGetLabel(dm, "Face Sets", &label); CHKERRQ(ierr);
31077841947SLeila Ghaffari     PetscInt comps[1] = {1};
311b8962995SJeremy L Thompson     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", label,
31277841947SLeila Ghaffari                          bc->num_slip[0], bc->slips[0], 0, 1, comps,
31377841947SLeila Ghaffari                          (void(*)(void))NULL, NULL, setup_ctx, NULL);
31477841947SLeila Ghaffari     CHKERRQ(ierr);
31577841947SLeila Ghaffari     comps[0] = 2;
316b8962995SJeremy L Thompson     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", label,
31777841947SLeila Ghaffari                          bc->num_slip[1], bc->slips[1], 0, 1, comps,
31877841947SLeila Ghaffari                          (void(*)(void))NULL, NULL, setup_ctx, NULL);
31977841947SLeila Ghaffari     CHKERRQ(ierr);
32077841947SLeila Ghaffari     comps[0] = 3;
321b8962995SJeremy L Thompson     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", label,
32277841947SLeila Ghaffari                          bc->num_slip[2], bc->slips[2], 0, 1, comps,
32377841947SLeila Ghaffari                          (void(*)(void))NULL, NULL, setup_ctx, NULL);
32477841947SLeila Ghaffari     CHKERRQ(ierr);
32577841947SLeila Ghaffari   }
32677841947SLeila Ghaffari 
327e6225c47SLeila Ghaffari   if (bc->user_bc) {
32877841947SLeila Ghaffari     for (PetscInt c = 0; c < 3; c++) {
32977841947SLeila Ghaffari       for (PetscInt s = 0; s < bc->num_slip[c]; s++) {
33077841947SLeila Ghaffari         for (PetscInt w = 0; w < bc->num_wall; w++) {
33177841947SLeila Ghaffari           if (bc->slips[c][s] == bc->walls[w])
33277841947SLeila Ghaffari             SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG,
33377841947SLeila Ghaffari                      "Boundary condition already set on face %D!\n",
33477841947SLeila Ghaffari                      bc->walls[w]);
33577841947SLeila Ghaffari         }
33677841947SLeila Ghaffari       }
33777841947SLeila Ghaffari     }
33877841947SLeila Ghaffari   }
33977841947SLeila Ghaffari 
34077841947SLeila Ghaffari   // Set wall boundary conditions
34177841947SLeila Ghaffari   //   zero velocity and zero flux for mass density and energy density
34277841947SLeila Ghaffari   {
34377841947SLeila Ghaffari     DMLabel  label;
34477841947SLeila Ghaffari     PetscInt comps[3] = {1, 2, 3};
34577841947SLeila Ghaffari     ierr = DMGetLabel(dm, "Face Sets", &label); CHKERRQ(ierr);
346b8962995SJeremy L Thompson     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label,
347b8962995SJeremy L Thompson                          bc->num_wall, bc->walls, 0, 3, comps,
348b8962995SJeremy L Thompson                          (void(*)(void))Exact_DC, NULL,
34977841947SLeila Ghaffari                          setup_ctx, NULL); CHKERRQ(ierr);
35077841947SLeila Ghaffari   }
35177841947SLeila Ghaffari 
35277841947SLeila Ghaffari   PetscFunctionReturn(0);
35377841947SLeila Ghaffari }
35477841947SLeila Ghaffari 
35577841947SLeila Ghaffari PetscErrorCode PRINT_DENSITY_CURRENT(Physics phys, SetupContext setup_ctx,
35677841947SLeila Ghaffari                                      AppCtx app_ctx) {
35777841947SLeila Ghaffari   MPI_Comm       comm = PETSC_COMM_WORLD;
35877841947SLeila Ghaffari   PetscErrorCode ierr;
35977841947SLeila Ghaffari   PetscFunctionBeginUser;
36077841947SLeila Ghaffari 
36177841947SLeila Ghaffari   ierr = PetscPrintf(comm,
36277841947SLeila Ghaffari                      "  Problem:\n"
36377841947SLeila Ghaffari                      "    Problem Name                       : %s\n"
36477841947SLeila Ghaffari                      "    Stabilization                      : %s\n",
36577841947SLeila Ghaffari                      app_ctx->problem_name, StabilizationTypes[phys->stab]);
36677841947SLeila Ghaffari   CHKERRQ(ierr);
36777841947SLeila Ghaffari 
36877841947SLeila Ghaffari   PetscFunctionReturn(0);
36977841947SLeila Ghaffari }
370