xref: /libCEED/examples/fluids/problems/newtonian.c (revision 91e5af17cc38f7d4d16af29119518d78e5298d3f)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
388b783a1SJames Wright //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
588b783a1SJames Wright //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
788b783a1SJames Wright 
888b783a1SJames Wright /// @file
988b783a1SJames Wright /// Utility functions for setting up problems using the Newtonian Qfunction
1088b783a1SJames Wright 
1188b783a1SJames Wright #include "../navierstokes.h"
1288b783a1SJames Wright #include "../qfunctions/setupgeo.h"
1388b783a1SJames Wright #include "../qfunctions/newtonian.h"
1488b783a1SJames Wright 
1588626eedSJames Wright 
1688626eedSJames Wright #ifndef newtonian_context_struct
1788626eedSJames Wright #define newtonian_context_struct
1888626eedSJames Wright typedef struct NewtonianIdealGasContext_ *NewtonianIdealGasContext;
1988626eedSJames Wright struct NewtonianIdealGasContext_ {
2088626eedSJames Wright   CeedScalar lambda;
2188626eedSJames Wright   CeedScalar mu;
2288626eedSJames Wright   CeedScalar k;
2388626eedSJames Wright   CeedScalar cv;
2488626eedSJames Wright   CeedScalar cp;
2588626eedSJames Wright   CeedScalar g[3];
2688626eedSJames Wright   CeedScalar c_tau;
2788626eedSJames Wright   CeedScalar Ctau_t;
2888626eedSJames Wright   CeedScalar Ctau_v;
2988626eedSJames Wright   CeedScalar Ctau_C;
3088626eedSJames Wright   CeedScalar Ctau_M;
3188626eedSJames Wright   CeedScalar Ctau_E;
3288626eedSJames Wright   CeedScalar dt;
3388626eedSJames Wright   StabilizationType stabilization;
3488626eedSJames Wright };
3588626eedSJames Wright #endif
3688626eedSJames Wright 
3788b783a1SJames Wright PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *setup_ctx,
3888b783a1SJames Wright                                void *ctx) {
3988b783a1SJames Wright   SetupContext      setup_context = *(SetupContext *)setup_ctx;
4088b783a1SJames Wright   User              user = *(User *)ctx;
4188b783a1SJames Wright   StabilizationType stab;
4288b783a1SJames Wright   MPI_Comm          comm = PETSC_COMM_WORLD;
4388b783a1SJames Wright   PetscBool         implicit;
4488b783a1SJames Wright   PetscBool         has_curr_time = PETSC_FALSE;
4588b783a1SJames Wright   PetscInt          ierr;
4688b783a1SJames Wright   PetscFunctionBeginUser;
4788b783a1SJames Wright 
4888b783a1SJames Wright   ierr = PetscCalloc1(1, &user->phys->newtonian_ig_ctx); CHKERRQ(ierr);
4988b783a1SJames Wright 
5088b783a1SJames Wright   // ------------------------------------------------------
5188b783a1SJames Wright   //           Setup Generic Newtonian IG Problem
5288b783a1SJames Wright   // ------------------------------------------------------
5388b783a1SJames Wright   problem->dim                               = 3;
5488b783a1SJames Wright   problem->q_data_size_vol                   = 10;
5588b783a1SJames Wright   problem->q_data_size_sur                   = 4;
56*91e5af17SJed Brown   problem->setup_vol.qfunction               = Setup;
57*91e5af17SJed Brown   problem->setup_vol.qfunction_loc           = Setup_loc;
58*91e5af17SJed Brown   problem->ics.qfunction                     = ICsNewtonianIG;
59*91e5af17SJed Brown   problem->ics.qfunction_loc                 = ICsNewtonianIG_loc;
60*91e5af17SJed Brown   problem->setup_sur.qfunction               = SetupBoundary;
61*91e5af17SJed Brown   problem->setup_sur.qfunction_loc           = SetupBoundary_loc;
62*91e5af17SJed Brown   problem->apply_vol_rhs.qfunction           = Newtonian;
63*91e5af17SJed Brown   problem->apply_vol_rhs.qfunction_loc       = Newtonian_loc;
64*91e5af17SJed Brown   problem->apply_vol_ifunction.qfunction     = IFunction_Newtonian;
65*91e5af17SJed Brown   problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_loc;
6688b783a1SJames Wright   problem->setup_ctx                         = SetupContext_DENSITY_CURRENT;
6788b783a1SJames Wright   problem->non_zero_time                     = PETSC_FALSE;
6888b783a1SJames Wright   problem->print_info                        = PRINT_DENSITY_CURRENT;
6988b783a1SJames Wright 
7088b783a1SJames Wright   // ------------------------------------------------------
7188b783a1SJames Wright   //             Create the libCEED context
7288b783a1SJames Wright   // ------------------------------------------------------
7388b783a1SJames Wright   CeedScalar cv     = 717.;          // J/(kg K)
7488b783a1SJames Wright   CeedScalar cp     = 1004.;         // J/(kg K)
7588626eedSJames Wright   CeedScalar g[3]   = {0, 0, -9.81}; // m/s^2
7688b783a1SJames Wright   CeedScalar lambda = -2./3.;        // -
7788626eedSJames Wright   CeedScalar mu     = 1.8e-5;        // Pa s, dynamic viscosity
7888b783a1SJames Wright   CeedScalar k      = 0.02638;       // W/(m K)
7988b783a1SJames Wright   CeedScalar c_tau  = 0.5;           // -
8088626eedSJames Wright   CeedScalar Ctau_t  = 1.0;          // -
8188626eedSJames Wright   CeedScalar Ctau_v  = 36.0;         // TODO make function of degree
8288626eedSJames Wright   CeedScalar Ctau_C  = 1.0;          // TODO make function of degree
8388626eedSJames Wright   CeedScalar Ctau_M  = 1.0;          // TODO make function of degree
8488626eedSJames Wright   CeedScalar Ctau_E  = 1.0;          // TODO make function of degree
8588b783a1SJames Wright   PetscReal domain_min[3], domain_max[3], domain_size[3];
8688b783a1SJames Wright   ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
8788b783a1SJames Wright   for (int i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i];
8888b783a1SJames Wright 
8988b783a1SJames Wright   // ------------------------------------------------------
9088b783a1SJames Wright   //             Create the PETSc context
9188b783a1SJames Wright   // ------------------------------------------------------
9288626eedSJames Wright   PetscScalar meter    = 1;  // 1 meter in scaled length units
9388626eedSJames Wright   PetscScalar kilogram = 1;  // 1 kilogram in scaled mass units
9488626eedSJames Wright   PetscScalar second   = 1;  // 1 second in scaled time units
9588b783a1SJames Wright   PetscScalar Kelvin   = 1;     // 1 Kelvin in scaled temperature units
9688b783a1SJames Wright   PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s;
9788b783a1SJames Wright 
9888b783a1SJames Wright   // ------------------------------------------------------
9988b783a1SJames Wright   //              Command line Options
10088b783a1SJames Wright   // ------------------------------------------------------
10167490bc6SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem",
10267490bc6SJeremy L Thompson                     NULL);
10367490bc6SJeremy L Thompson 
10488b783a1SJames Wright   // -- Physics
10588b783a1SJames Wright   ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume",
10688b783a1SJames Wright                             NULL, cv, &cv, NULL); CHKERRQ(ierr);
10788b783a1SJames Wright   ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure",
10888b783a1SJames Wright                             NULL, cp, &cp, NULL); CHKERRQ(ierr);
10988b783a1SJames Wright   ierr = PetscOptionsScalar("-lambda",
11088b783a1SJames Wright                             "Stokes hypothesis second viscosity coefficient",
11188b783a1SJames Wright                             NULL, lambda, &lambda, NULL); CHKERRQ(ierr);
11288b783a1SJames Wright   ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient",
11388b783a1SJames Wright                             NULL, mu, &mu, NULL); CHKERRQ(ierr);
11488b783a1SJames Wright   ierr = PetscOptionsScalar("-k", "Thermal conductivity",
11588b783a1SJames Wright                             NULL, k, &k, NULL); CHKERRQ(ierr);
11688b783a1SJames Wright 
11788626eedSJames Wright   PetscInt dim = problem->dim;
11888626eedSJames Wright   ierr = PetscOptionsRealArray("-g", "Gravitational acceleration",
11988626eedSJames Wright                                NULL, g, &dim, NULL); CHKERRQ(ierr);
12088b783a1SJames Wright   ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL,
12188b783a1SJames Wright                           StabilizationTypes, (PetscEnum)(stab = STAB_NONE),
12288b783a1SJames Wright                           (PetscEnum *)&stab, NULL); CHKERRQ(ierr);
12388b783a1SJames Wright   ierr = PetscOptionsScalar("-c_tau", "Stabilization constant",
12488b783a1SJames Wright                             NULL, c_tau, &c_tau, NULL); CHKERRQ(ierr);
12588626eedSJames Wright   ierr = PetscOptionsScalar("-Ctau_t", "Stabilization time constant",
12688626eedSJames Wright                             NULL, Ctau_t, &Ctau_t, NULL); CHKERRQ(ierr);
12788626eedSJames Wright   ierr = PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant",
12888626eedSJames Wright                             NULL, Ctau_v, &Ctau_v, NULL); CHKERRQ(ierr);
12988626eedSJames Wright   ierr = PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant",
13088626eedSJames Wright                             NULL, Ctau_C, &Ctau_C, NULL); CHKERRQ(ierr);
13188626eedSJames Wright   ierr = PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant",
13288626eedSJames Wright                             NULL, Ctau_M, &Ctau_M, NULL); CHKERRQ(ierr);
13388626eedSJames Wright   ierr = PetscOptionsScalar("-Ctau_E", "Stabilization energy constant",
13488626eedSJames Wright                             NULL, Ctau_E, &Ctau_E, NULL); CHKERRQ(ierr);
13588b783a1SJames Wright   ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation",
13688b783a1SJames Wright                           NULL, implicit=PETSC_FALSE, &implicit, NULL);
13788b783a1SJames Wright   CHKERRQ(ierr);
13888b783a1SJames Wright 
13988b783a1SJames Wright   // -- Units
14088b783a1SJames Wright   ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units",
14188b783a1SJames Wright                             NULL, meter, &meter, NULL); CHKERRQ(ierr);
14288b783a1SJames Wright   meter = fabs(meter);
14388b783a1SJames Wright   ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units",
14488b783a1SJames Wright                             NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr);
14588b783a1SJames Wright   kilogram = fabs(kilogram);
14688b783a1SJames Wright   ierr = PetscOptionsScalar("-units_second","1 second in scaled time units",
14788b783a1SJames Wright                             NULL, second, &second, NULL); CHKERRQ(ierr);
14888b783a1SJames Wright   second = fabs(second);
14988b783a1SJames Wright   ierr = PetscOptionsScalar("-units_Kelvin",
15088b783a1SJames Wright                             "1 Kelvin in scaled temperature units",
15188b783a1SJames Wright                             NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr);
15288b783a1SJames Wright   Kelvin = fabs(Kelvin);
15388b783a1SJames Wright 
15488b783a1SJames Wright   // -- Warnings
15588b783a1SJames Wright   if (stab == STAB_SUPG && !implicit) {
15688b783a1SJames Wright     ierr = PetscPrintf(comm,
15788b783a1SJames Wright                        "Warning! Use -stab supg only with -implicit\n");
15888b783a1SJames Wright     CHKERRQ(ierr);
15988b783a1SJames Wright   }
16067490bc6SJeremy L Thompson   PetscOptionsEnd();
16188b783a1SJames Wright 
16288b783a1SJames Wright   // ------------------------------------------------------
16388b783a1SJames Wright   //           Set up the PETSc context
16488b783a1SJames Wright   // ------------------------------------------------------
16588b783a1SJames Wright   // -- Define derived units
16688b783a1SJames Wright   Pascal          = kilogram / (meter * PetscSqr(second));
16788b783a1SJames Wright   J_per_kg_K      =  PetscSqr(meter) / (PetscSqr(second) * Kelvin);
16888b783a1SJames Wright   m_per_squared_s = meter / PetscSqr(second);
16988b783a1SJames Wright   W_per_m_K       = kilogram * meter / (pow(second,3) * Kelvin);
17088b783a1SJames Wright 
17188b783a1SJames Wright   user->units->meter           = meter;
17288b783a1SJames Wright   user->units->kilogram        = kilogram;
17388b783a1SJames Wright   user->units->second          = second;
17488b783a1SJames Wright   user->units->Kelvin          = Kelvin;
17588b783a1SJames Wright   user->units->Pascal          = Pascal;
17688b783a1SJames Wright   user->units->J_per_kg_K      = J_per_kg_K;
17788b783a1SJames Wright   user->units->m_per_squared_s = m_per_squared_s;
17888b783a1SJames Wright   user->units->W_per_m_K       = W_per_m_K;
17988b783a1SJames Wright 
18088b783a1SJames Wright   // ------------------------------------------------------
18188b783a1SJames Wright   //           Set up the libCEED context
18288b783a1SJames Wright   // ------------------------------------------------------
18388b783a1SJames Wright   // -- Scale variables to desired units
18488b783a1SJames Wright   cv     *= J_per_kg_K;
18588b783a1SJames Wright   cp     *= J_per_kg_K;
18688b783a1SJames Wright   mu     *= Pascal * second;
18788b783a1SJames Wright   k      *= W_per_m_K;
18888b783a1SJames Wright   for (int i=0; i<3; i++) domain_size[i] *= meter;
18988626eedSJames Wright   for (int i=0; i<3; i++) g[i]           *= m_per_squared_s;
19088b783a1SJames Wright   problem->dm_scale = meter;
19188b783a1SJames Wright 
19288b783a1SJames Wright   // -- Setup Context
19388b783a1SJames Wright   setup_context->cv         = cv;
19488b783a1SJames Wright   setup_context->cp         = cp;
19588b783a1SJames Wright   setup_context->lx         = domain_size[0];
19688b783a1SJames Wright   setup_context->ly         = domain_size[1];
19788b783a1SJames Wright   setup_context->lz         = domain_size[2];
19888b783a1SJames Wright   setup_context->time       = 0;
19988626eedSJames Wright   ierr = PetscArraycpy(setup_context->g, g, 3); CHKERRQ(ierr);
20088b783a1SJames Wright 
20188b783a1SJames Wright   // -- Solver Settings
20288b783a1SJames Wright   user->phys->stab          = stab;
20388b783a1SJames Wright   user->phys->implicit      = implicit;
20488b783a1SJames Wright   user->phys->has_curr_time = has_curr_time;
20588b783a1SJames Wright 
20688b783a1SJames Wright   // -- QFunction Context
20788b783a1SJames Wright   user->phys->newtonian_ig_ctx->lambda        = lambda;
20888b783a1SJames Wright   user->phys->newtonian_ig_ctx->mu            = mu;
20988b783a1SJames Wright   user->phys->newtonian_ig_ctx->k             = k;
21088b783a1SJames Wright   user->phys->newtonian_ig_ctx->cv            = cv;
21188b783a1SJames Wright   user->phys->newtonian_ig_ctx->cp            = cp;
21288b783a1SJames Wright   user->phys->newtonian_ig_ctx->c_tau         = c_tau;
21388626eedSJames Wright   user->phys->newtonian_ig_ctx->Ctau_t        = Ctau_t;
21488626eedSJames Wright   user->phys->newtonian_ig_ctx->Ctau_v        = Ctau_v;
21588626eedSJames Wright   user->phys->newtonian_ig_ctx->Ctau_C        = Ctau_C;
21688626eedSJames Wright   user->phys->newtonian_ig_ctx->Ctau_M        = Ctau_M;
21788626eedSJames Wright   user->phys->newtonian_ig_ctx->Ctau_E        = Ctau_E;
21888b783a1SJames Wright   user->phys->newtonian_ig_ctx->stabilization = stab;
21988626eedSJames Wright   ierr = PetscArraycpy(user->phys->newtonian_ig_ctx->g, g, 3); CHKERRQ(ierr);
22088b783a1SJames Wright 
22188b783a1SJames Wright   PetscFunctionReturn(0);
22288b783a1SJames Wright }
22388b783a1SJames Wright 
22488b783a1SJames Wright PetscErrorCode SetupContext_NEWTONIAN_IG(Ceed ceed, CeedData ceed_data,
22588b783a1SJames Wright     AppCtx app_ctx, SetupContext setup_ctx, Physics phys) {
22688b783a1SJames Wright   PetscFunctionBeginUser;
22788b783a1SJames Wright   CeedQFunctionContextCreate(ceed, &ceed_data->setup_context);
22888b783a1SJames Wright   CeedQFunctionContextSetData(ceed_data->setup_context, CEED_MEM_HOST,
22988b783a1SJames Wright                               CEED_USE_POINTER, sizeof(*setup_ctx), setup_ctx);
23088b783a1SJames Wright   CeedQFunctionSetContext(ceed_data->qf_ics, ceed_data->setup_context);
23188626eedSJames Wright 
23288b783a1SJames Wright   CeedQFunctionContextCreate(ceed, &ceed_data->newt_ig_context);
23388b783a1SJames Wright   CeedQFunctionContextSetData(ceed_data->newt_ig_context, CEED_MEM_HOST,
23488b783a1SJames Wright                               CEED_USE_POINTER,
23588b783a1SJames Wright                               sizeof(*phys->newtonian_ig_ctx), phys->newtonian_ig_ctx);
23688626eedSJames Wright   CeedQFunctionContextRegisterDouble(ceed_data->newt_ig_context, "timestep size",
23788626eedSJames Wright                                      offsetof(struct NewtonianIdealGasContext_, dt), 1, "Size of timestep, delta t");
23888626eedSJames Wright 
23988b783a1SJames Wright   if (ceed_data->qf_rhs_vol)
24088b783a1SJames Wright     CeedQFunctionSetContext(ceed_data->qf_rhs_vol, ceed_data->newt_ig_context);
24188b783a1SJames Wright   if (ceed_data->qf_ifunction_vol)
24288b783a1SJames Wright     CeedQFunctionSetContext(ceed_data->qf_ifunction_vol,
24388b783a1SJames Wright                             ceed_data->newt_ig_context);
24488b783a1SJames Wright   PetscFunctionReturn(0);
24588b783a1SJames Wright }
246