xref: /libCEED/examples/fluids/problems/eulervortex.c (revision ba6664ae303f5b2ef46b3df96973d9bdc665107c)
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   EulerTestType     euler_test;
17   User              user = *(User *)ctx;
18   StabilizationType stab;
19   MPI_Comm          comm = PETSC_COMM_WORLD;
20   PetscBool         implicit;
21   PetscBool         has_curr_time = PETSC_TRUE;
22   PetscBool         has_neumann = PETSC_TRUE;
23   PetscInt          ierr;
24   EulerContext      euler_ctx;
25   CeedQFunctionContext euler_context;
26 
27   PetscFunctionBeginUser;
28   ierr = PetscCalloc1(1, &euler_ctx); CHKERRQ(ierr);
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   ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
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   ierr = PetscOptionsScalar("-vortex_strength", "Strength of Vortex",
79                             NULL, vortex_strength, &vortex_strength, NULL);
80   CHKERRQ(ierr);
81   PetscInt n = problem->dim;
82   PetscBool user_velocity;
83   ierr = PetscOptionsRealArray("-mean_velocity", "Background velocity vector",
84                                NULL, mean_velocity, &n, &user_velocity);
85   CHKERRQ(ierr);
86   for (PetscInt i=0; i<3; i++) center[i] = .5*domain_size[i];
87   n = problem->dim;
88   ierr = PetscOptionsRealArray("-center", "Location of vortex center",
89                                NULL, center, &n, NULL); CHKERRQ(ierr);
90   ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation",
91                           NULL, implicit=PETSC_FALSE, &implicit, NULL);
92   CHKERRQ(ierr);
93   ierr = PetscOptionsEnum("-euler_test", "Euler test option", NULL,
94                           EulerTestTypes, (PetscEnum)(euler_test = EULER_TEST_ISENTROPIC_VORTEX),
95                           (PetscEnum *)&euler_test, NULL); CHKERRQ(ierr);
96   ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL,
97                           StabilizationTypes, (PetscEnum)(stab = STAB_NONE),
98                           (PetscEnum *)&stab, NULL); CHKERRQ(ierr);
99   ierr = PetscOptionsScalar("-c_tau", "Stabilization constant",
100                             NULL, c_tau, &c_tau, NULL); CHKERRQ(ierr);
101   // -- Units
102   ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units",
103                             NULL, meter, &meter, NULL); CHKERRQ(ierr);
104   meter = fabs(meter);
105   ierr = PetscOptionsScalar("-units_second","1 second in scaled time units",
106                             NULL, second, &second, NULL); CHKERRQ(ierr);
107   second = fabs(second);
108 
109   // -- Warnings
110   if (stab == STAB_SUPG && !implicit) {
111     ierr = PetscPrintf(comm,
112                        "Warning! Use -stab supg only with -implicit\n");
113     CHKERRQ(ierr);
114   }
115   if (user_velocity && (euler_test == EULER_TEST_1
116                         || euler_test == EULER_TEST_3)) {
117     ierr = PetscPrintf(comm,
118                        "Warning! Background velocity vector for -euler_test t1 and -euler_test t3 is (0,0,0)\n");
119     CHKERRQ(ierr);
120   }
121 
122   PetscOptionsEnd();
123 
124   // ------------------------------------------------------
125   //           Set up the PETSc context
126   // ------------------------------------------------------
127   user->units->meter  = meter;
128   user->units->second = second;
129 
130   // ------------------------------------------------------
131   //           Set up the libCEED context
132   // ------------------------------------------------------
133   // -- Scale variables to desired units
134   for (PetscInt i=0; i<3; i++) {
135     center[i] *= meter;
136     domain_size[i] *= meter;
137     mean_velocity[i] *= (meter/second);
138   }
139   problem->dm_scale = meter;
140 
141   // -- QFunction Context
142   user->phys->stab                        = stab;
143   user->phys->euler_test                  = euler_test;
144   user->phys->implicit                    = implicit;
145   user->phys->has_curr_time               = has_curr_time;
146   user->phys->has_neumann                 = has_neumann;
147   euler_ctx->curr_time        = 0.;
148   euler_ctx->implicit         = implicit;
149   euler_ctx->euler_test       = euler_test;
150   euler_ctx->center[0]        = center[0];
151   euler_ctx->center[1]        = center[1];
152   euler_ctx->center[2]        = center[2];
153   euler_ctx->vortex_strength  = vortex_strength;
154   euler_ctx->c_tau            = c_tau;
155   euler_ctx->mean_velocity[0] = mean_velocity[0];
156   euler_ctx->mean_velocity[1] = mean_velocity[1];
157   euler_ctx->mean_velocity[2] = mean_velocity[2];
158   euler_ctx->stabilization    = stab;
159 
160   CeedQFunctionContextCreate(user->ceed, &euler_context);
161   CeedQFunctionContextSetData(euler_context, CEED_MEM_HOST, CEED_USE_POINTER,
162                               sizeof(*euler_ctx), euler_ctx);
163   CeedQFunctionContextSetDataDestroy(euler_context, CEED_MEM_HOST,
164                                      FreeContextPetsc);
165   CeedQFunctionContextRegisterDouble(euler_context, "solution time",
166                                      offsetof(struct EulerContext_, curr_time), 1, "Phyiscal time of the solution");
167   CeedQFunctionContextReferenceCopy(euler_context,
168                                     &problem->ics.qfunction_context);
169   CeedQFunctionContextReferenceCopy(euler_context,
170                                     &problem->apply_vol_rhs.qfunction_context);
171   CeedQFunctionContextReferenceCopy(euler_context,
172                                     &problem->apply_vol_ifunction.qfunction_context);
173   CeedQFunctionContextReferenceCopy(euler_context,
174                                     &problem->apply_inflow.qfunction_context);
175   CeedQFunctionContextReferenceCopy(euler_context,
176                                     &problem->apply_outflow.qfunction_context);
177   PetscFunctionReturn(0);
178 }
179 
180 PetscErrorCode PRINT_EULER_VORTEX(ProblemData *problem,
181                                   AppCtx app_ctx) {
182   MPI_Comm       comm = PETSC_COMM_WORLD;
183   PetscErrorCode ierr;
184   EulerContext   euler_ctx;
185 
186   PetscFunctionBeginUser;
187   CeedQFunctionContextGetData(problem->ics.qfunction_context, CEED_MEM_HOST,
188                               &euler_ctx);
189   ierr = PetscPrintf(comm,
190                      "  Problem:\n"
191                      "    Problem Name                       : %s\n"
192                      "    Test Case                          : %s\n"
193                      "    Background Velocity                : %f,%f,%f\n"
194                      "    Stabilization                      : %s\n",
195                      app_ctx->problem_name, EulerTestTypes[euler_ctx->euler_test],
196                      euler_ctx->mean_velocity[0],
197                      euler_ctx->mean_velocity[1],
198                      euler_ctx->mean_velocity[2],
199                      StabilizationTypes[euler_ctx->stabilization]); CHKERRQ(ierr);
200 
201   CeedQFunctionContextRestoreData(problem->ics.qfunction_context, &euler_ctx);
202   PetscFunctionReturn(0);
203 }
204