xref: /honee/problems/eulervortex.c (revision 3182272fe6b846b1236f33efc33a9283274e9953)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 /// @file
9 /// Utility functions for setting up EULER_VORTEX
10 
11 #include "../qfunctions/eulervortex.h"
12 
13 #include "../navierstokes.h"
14 #include "../qfunctions/setupgeo.h"
15 
16 PetscErrorCode NS_EULER_VORTEX(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) {
17   EulerTestType        euler_test;
18   User                 user = *(User *)ctx;
19   StabilizationType    stab;
20   MPI_Comm             comm = PETSC_COMM_WORLD;
21   PetscBool            implicit;
22   PetscBool            has_curr_time = PETSC_TRUE;
23   PetscBool            has_neumann   = PETSC_TRUE;
24   EulerContext         euler_ctx;
25   CeedQFunctionContext euler_context;
26 
27   PetscFunctionBeginUser;
28   PetscCall(PetscCalloc1(1, &euler_ctx));
29 
30   // ------------------------------------------------------
31   //               SET UP DENSITY_CURRENT
32   // ------------------------------------------------------
33   problem->dim                               = 3;
34   problem->q_data_size_vol                   = 10;
35   problem->q_data_size_sur                   = 10;
36   problem->setup_vol.qfunction               = Setup;
37   problem->setup_vol.qfunction_loc           = Setup_loc;
38   problem->setup_sur.qfunction               = SetupBoundary;
39   problem->setup_sur.qfunction_loc           = SetupBoundary_loc;
40   problem->ics.qfunction                     = ICsEuler;
41   problem->ics.qfunction_loc                 = ICsEuler_loc;
42   problem->apply_vol_rhs.qfunction           = Euler;
43   problem->apply_vol_rhs.qfunction_loc       = Euler_loc;
44   problem->apply_vol_ifunction.qfunction     = IFunction_Euler;
45   problem->apply_vol_ifunction.qfunction_loc = IFunction_Euler_loc;
46   problem->apply_inflow.qfunction            = TravelingVortex_Inflow;
47   problem->apply_inflow.qfunction_loc        = TravelingVortex_Inflow_loc;
48   problem->apply_outflow.qfunction           = Euler_Outflow;
49   problem->apply_outflow.qfunction_loc       = Euler_Outflow_loc;
50   problem->bc                                = Exact_Euler;
51   problem->bc_ctx                            = euler_ctx;
52   problem->non_zero_time                     = PETSC_TRUE;
53   problem->print_info                        = PRINT_EULER_VORTEX;
54 
55   // ------------------------------------------------------
56   //             Create the libCEED context
57   // ------------------------------------------------------
58   CeedScalar vortex_strength = 5.;   // -
59   CeedScalar c_tau           = 0.5;  // -
60   // c_tau = 0.5 is reported as "optimal" in Hughes et al 2010
61   PetscReal center[3],                 // m
62       mean_velocity[3] = {1., 1., 0};  // m/s
63   PetscReal domain_min[3], domain_max[3], domain_size[3];
64   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
65   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
66 
67   // ------------------------------------------------------
68   //             Create the PETSc context
69   // ------------------------------------------------------
70   PetscScalar meter  = 1e-2;  // 1 meter in scaled length units
71   PetscScalar second = 1e-2;  // 1 second in scaled time units
72 
73   // ------------------------------------------------------
74   //              Command line Options
75   // ------------------------------------------------------
76   PetscOptionsBegin(comm, NULL, "Options for EULER_VORTEX problem", NULL);
77   // -- Physics
78   PetscCall(PetscOptionsScalar("-vortex_strength", "Strength of Vortex", NULL, vortex_strength, &vortex_strength, NULL));
79   PetscInt  n = problem->dim;
80   PetscBool user_velocity;
81   PetscCall(PetscOptionsRealArray("-mean_velocity", "Background velocity vector", NULL, mean_velocity, &n, &user_velocity));
82   for (PetscInt i = 0; i < 3; i++) center[i] = .5 * domain_size[i];
83   n = problem->dim;
84   PetscCall(PetscOptionsRealArray("-center", "Location of vortex center", NULL, center, &n, NULL));
85   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
86   PetscCall(PetscOptionsEnum("-euler_test", "Euler test option", NULL, EulerTestTypes, (PetscEnum)(euler_test = EULER_TEST_ISENTROPIC_VORTEX),
87                              (PetscEnum *)&euler_test, NULL));
88   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
89   PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL));
90   // -- Units
91   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
92   meter = fabs(meter);
93   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
94   second = fabs(second);
95 
96   // -- Warnings
97   if (stab == STAB_SUPG && !implicit) {
98     PetscCall(PetscPrintf(comm, "Warning! Use -stab supg only with -implicit\n"));
99   }
100   if (user_velocity && (euler_test == EULER_TEST_1 || euler_test == EULER_TEST_3)) {
101     PetscCall(PetscPrintf(comm, "Warning! Background velocity vector for -euler_test t1 and -euler_test t3 is (0,0,0)\n"));
102   }
103 
104   PetscOptionsEnd();
105 
106   // ------------------------------------------------------
107   //           Set up the PETSc context
108   // ------------------------------------------------------
109   user->units->meter  = meter;
110   user->units->second = second;
111 
112   // ------------------------------------------------------
113   //           Set up the libCEED context
114   // ------------------------------------------------------
115   // -- Scale variables to desired units
116   for (PetscInt i = 0; i < 3; i++) {
117     center[i] *= meter;
118     domain_size[i] *= meter;
119     mean_velocity[i] *= (meter / second);
120   }
121   problem->dm_scale = meter;
122 
123   // -- QFunction Context
124   user->phys->stab            = stab;
125   user->phys->euler_test      = euler_test;
126   user->phys->implicit        = implicit;
127   user->phys->has_curr_time   = has_curr_time;
128   user->phys->has_neumann     = has_neumann;
129   euler_ctx->curr_time        = 0.;
130   euler_ctx->implicit         = implicit;
131   euler_ctx->euler_test       = euler_test;
132   euler_ctx->center[0]        = center[0];
133   euler_ctx->center[1]        = center[1];
134   euler_ctx->center[2]        = center[2];
135   euler_ctx->vortex_strength  = vortex_strength;
136   euler_ctx->c_tau            = c_tau;
137   euler_ctx->mean_velocity[0] = mean_velocity[0];
138   euler_ctx->mean_velocity[1] = mean_velocity[1];
139   euler_ctx->mean_velocity[2] = mean_velocity[2];
140   euler_ctx->stabilization    = stab;
141 
142   CeedQFunctionContextCreate(user->ceed, &euler_context);
143   CeedQFunctionContextSetData(euler_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*euler_ctx), euler_ctx);
144   CeedQFunctionContextSetDataDestroy(euler_context, CEED_MEM_HOST, FreeContextPetsc);
145   CeedQFunctionContextRegisterDouble(euler_context, "solution time", offsetof(struct EulerContext_, curr_time), 1, "Physical time of the solution");
146   CeedQFunctionContextReferenceCopy(euler_context, &problem->ics.qfunction_context);
147   CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_vol_rhs.qfunction_context);
148   CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_vol_ifunction.qfunction_context);
149   CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_inflow.qfunction_context);
150   CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_outflow.qfunction_context);
151   PetscFunctionReturn(0);
152 }
153 
154 PetscErrorCode PRINT_EULER_VORTEX(ProblemData *problem, AppCtx app_ctx) {
155   MPI_Comm     comm = PETSC_COMM_WORLD;
156   EulerContext euler_ctx;
157 
158   PetscFunctionBeginUser;
159   CeedQFunctionContextGetData(problem->ics.qfunction_context, CEED_MEM_HOST, &euler_ctx);
160   PetscCall(PetscPrintf(comm,
161                         "  Problem:\n"
162                         "    Problem Name                       : %s\n"
163                         "    Test Case                          : %s\n"
164                         "    Background Velocity                : %f,%f,%f\n"
165                         "    Stabilization                      : %s\n",
166                         app_ctx->problem_name, EulerTestTypes[euler_ctx->euler_test], euler_ctx->mean_velocity[0], euler_ctx->mean_velocity[1],
167                         euler_ctx->mean_velocity[2], StabilizationTypes[euler_ctx->stabilization]));
168 
169   CeedQFunctionContextRestoreData(problem->ics.qfunction_context, &euler_ctx);
170   PetscFunctionReturn(0);
171 }
172