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