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