xref: /libCEED/examples/fluids/problems/eulervortex.c (revision 9c34f28efacf2a1982e0e33d5357c64f609fddc7)
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   EulerContext         euler_ctx;
27   CeedQFunctionContext euler_context;
28 
29   PetscFunctionBeginUser;
30   PetscCall(PetscCalloc1(1, &euler_ctx));
31 
32   // ------------------------------------------------------
33   //               SET UP DENSITY_CURRENT
34   // ------------------------------------------------------
35   problem->dim                               = 3;
36   problem->q_data_size_vol                   = 10;
37   problem->q_data_size_sur                   = 10;
38   problem->setup_vol.qfunction               = Setup;
39   problem->setup_vol.qfunction_loc           = Setup_loc;
40   problem->setup_sur.qfunction               = SetupBoundary;
41   problem->setup_sur.qfunction_loc           = SetupBoundary_loc;
42   problem->ics.qfunction                     = ICsEuler;
43   problem->ics.qfunction_loc                 = ICsEuler_loc;
44   problem->apply_vol_rhs.qfunction           = Euler;
45   problem->apply_vol_rhs.qfunction_loc       = Euler_loc;
46   problem->apply_vol_ifunction.qfunction     = IFunction_Euler;
47   problem->apply_vol_ifunction.qfunction_loc = IFunction_Euler_loc;
48   problem->apply_inflow.qfunction            = TravelingVortex_Inflow;
49   problem->apply_inflow.qfunction_loc        = TravelingVortex_Inflow_loc;
50   problem->apply_outflow.qfunction           = Euler_Outflow;
51   problem->apply_outflow.qfunction_loc       = Euler_Outflow_loc;
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->implicit        = implicit;
125   euler_ctx->curr_time        = 0.;
126   euler_ctx->implicit         = implicit;
127   euler_ctx->euler_test       = euler_test;
128   euler_ctx->center[0]        = center[0];
129   euler_ctx->center[1]        = center[1];
130   euler_ctx->center[2]        = center[2];
131   euler_ctx->vortex_strength  = vortex_strength;
132   euler_ctx->c_tau            = c_tau;
133   euler_ctx->mean_velocity[0] = mean_velocity[0];
134   euler_ctx->mean_velocity[1] = mean_velocity[1];
135   euler_ctx->mean_velocity[2] = mean_velocity[2];
136   euler_ctx->stabilization    = stab;
137 
138   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &euler_context));
139   PetscCallCeed(ceed, CeedQFunctionContextSetData(euler_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*euler_ctx), euler_ctx));
140   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(euler_context, CEED_MEM_HOST, FreeContextPetsc));
141   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(euler_context, "solution time", offsetof(struct EulerContext_, curr_time), 1,
142                                                          "Physical time of the solution"));
143   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_context, &problem->ics.qfunction_context));
144   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_vol_rhs.qfunction_context));
145   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_vol_ifunction.qfunction_context));
146   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_inflow.qfunction_context));
147   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_outflow.qfunction_context));
148   PetscFunctionReturn(PETSC_SUCCESS);
149 }
150 
151 PetscErrorCode PRINT_EULER_VORTEX(User user, ProblemData *problem, AppCtx app_ctx) {
152   MPI_Comm     comm = user->comm;
153   Ceed         ceed = user->ceed;
154   EulerContext euler_ctx;
155 
156   PetscFunctionBeginUser;
157   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfunction_context, CEED_MEM_HOST, &euler_ctx));
158   PetscCall(PetscPrintf(comm,
159                         "  Problem:\n"
160                         "    Problem Name                       : %s\n"
161                         "    Test Case                          : %s\n"
162                         "    Background Velocity                : %f,%f,%f\n"
163                         "    Stabilization                      : %s\n",
164                         app_ctx->problem_name, EulerTestTypes[euler_ctx->euler_test], euler_ctx->mean_velocity[0], euler_ctx->mean_velocity[1],
165                         euler_ctx->mean_velocity[2], StabilizationTypes[euler_ctx->stabilization]));
166 
167   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfunction_context, &euler_ctx));
168   PetscFunctionReturn(PETSC_SUCCESS);
169 }
170