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