xref: /honee/problems/eulervortex.c (revision d6cac2202708b595f7f0c5ad1c28f21a66db1a8c)
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 static PetscErrorCode EulerVortexOutflowBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) {
15   HoneeBCStruct honee_bc;
16 
17   PetscFunctionBeginUser;
18   PetscCall(BCDefinitionGetContext(bc_def, &honee_bc));
19   PetscCall(HoneeBCCreateIFunctionQF(bc_def, Euler_Outflow, Euler_Outflow_loc, honee_bc->qfctx, qf));
20   PetscFunctionReturn(PETSC_SUCCESS);
21 }
22 
23 static PetscErrorCode EulerVortexInflowBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) {
24   HoneeBCStruct honee_bc;
25 
26   PetscFunctionBeginUser;
27   PetscCall(BCDefinitionGetContext(bc_def, &honee_bc));
28   PetscCall(HoneeBCCreateIFunctionQF(bc_def, TravelingVortex_Inflow, TravelingVortex_Inflow_loc, honee_bc->qfctx, qf));
29   PetscFunctionReturn(PETSC_SUCCESS);
30 }
31 
32 PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx, SimpleBC bc) {
33   EulerTestType        euler_test;
34   Honee                honee = *(Honee *)ctx;
35   StabilizationType    stab;
36   MPI_Comm             comm = honee->comm;
37   Ceed                 ceed = honee->ceed;
38   PetscBool            implicit;
39   EulerContext         euler_ctx;
40   CeedQFunctionContext euler_qfctx;
41   PetscInt             dim;
42 
43   PetscFunctionBeginUser;
44   PetscCall(PetscCalloc1(1, &euler_ctx));
45 
46   // ------------------------------------------------------
47   //               SET UP EULER VORTEX
48   // ------------------------------------------------------
49   problem->ics.qf_func_ptr                 = ICsEuler;
50   problem->ics.qf_loc                      = ICsEuler_loc;
51   problem->apply_vol_rhs.qf_func_ptr       = Euler;
52   problem->apply_vol_rhs.qf_loc            = Euler_loc;
53   problem->apply_vol_ifunction.qf_func_ptr = IFunction_Euler;
54   problem->apply_vol_ifunction.qf_loc      = IFunction_Euler_loc;
55   problem->jac_data_size_vol               = 0;
56   problem->jac_data_size_sur               = 0;
57   problem->compute_exact_solution_error    = PETSC_TRUE;
58   problem->print_info                      = PRINT_EULER_VORTEX;
59 
60   // ------------------------------------------------------
61   //             Create the libCEED context
62   // ------------------------------------------------------
63   CeedScalar vortex_strength = 5.;   // -
64   CeedScalar c_tau           = 0.5;  // -
65   // c_tau = 0.5 is reported as "optimal" in Hughes et al 2010
66   PetscReal center[3]  = {0.},         // m
67       mean_velocity[3] = {1., 1., 0};  // m/s
68   PetscReal domain_min[3], domain_max[3], domain_size[3] = {0.};
69   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
70   PetscCall(DMGetDimension(dm, &dim));
71   for (PetscInt i = 0; i < dim; i++) domain_size[i] = domain_max[i] - domain_min[i];
72 
73   // ------------------------------------------------------
74   //             Create the PETSc context
75   // ------------------------------------------------------
76   PetscScalar meter  = 1e-2;  // 1 meter in scaled length units
77   PetscScalar second = 1e-2;  // 1 second in scaled time units
78 
79   // ------------------------------------------------------
80   //              Command line Options
81   // ------------------------------------------------------
82   PetscOptionsBegin(comm, NULL, "Options for EULER_VORTEX problem", NULL);
83   // -- Physics
84   PetscCall(PetscOptionsScalar("-vortex_strength", "Strength of Vortex", NULL, vortex_strength, &vortex_strength, NULL));
85   PetscInt  n = dim;
86   PetscBool user_velocity;
87   PetscCall(PetscOptionsRealArray("-mean_velocity", "Background velocity vector", NULL, mean_velocity, &n, &user_velocity));
88   for (PetscInt i = 0; i < dim; i++) center[i] = .5 * domain_size[i];
89   n = dim;
90   PetscCall(PetscOptionsRealArray("-center", "Location of vortex center", NULL, center, &n, NULL));
91   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
92   PetscCall(PetscOptionsEnum("-euler_test", "Euler test option", NULL, EulerTestTypes, (PetscEnum)(euler_test = EULER_TEST_ISENTROPIC_VORTEX),
93                              (PetscEnum *)&euler_test, NULL));
94   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
95   PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL));
96   // -- Units
97   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
98   meter = fabs(meter);
99   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
100   second = fabs(second);
101 
102   // -- Warnings
103   if (stab == STAB_SUPG && !implicit) {
104     PetscCall(PetscPrintf(comm, "Warning! Use -stab supg only with -implicit\n"));
105   }
106   if (user_velocity && (euler_test == EULER_TEST_1 || euler_test == EULER_TEST_3)) {
107     PetscCall(PetscPrintf(comm, "Warning! Background velocity vector for -euler_test t1 and -euler_test t3 is (0,0,0)\n"));
108   }
109 
110   PetscOptionsEnd();
111 
112   // ------------------------------------------------------
113   //           Set up the PETSc context
114   // ------------------------------------------------------
115   honee->units->meter  = meter;
116   honee->units->second = second;
117 
118   // ------------------------------------------------------
119   //           Set up the libCEED context
120   // ------------------------------------------------------
121   // -- Scale variables to desired units
122   for (PetscInt i = 0; i < 3; i++) {
123     center[i] *= meter;
124     domain_size[i] *= meter;
125     mean_velocity[i] *= (meter / second);
126   }
127 
128   // -- QFunction Context
129   honee->phys->implicit       = implicit;
130   euler_ctx->curr_time        = 0.;
131   euler_ctx->implicit         = implicit;
132   euler_ctx->euler_test       = euler_test;
133   euler_ctx->center[0]        = center[0];
134   euler_ctx->center[1]        = center[1];
135   euler_ctx->center[2]        = center[2];
136   euler_ctx->vortex_strength  = vortex_strength;
137   euler_ctx->c_tau            = c_tau;
138   euler_ctx->mean_velocity[0] = mean_velocity[0];
139   euler_ctx->mean_velocity[1] = mean_velocity[1];
140   euler_ctx->mean_velocity[2] = mean_velocity[2];
141   euler_ctx->stabilization    = stab;
142 
143   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &euler_qfctx));
144   PetscCallCeed(ceed, CeedQFunctionContextSetData(euler_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*euler_ctx), euler_ctx));
145   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(euler_qfctx, CEED_MEM_HOST, FreeContextPetsc));
146   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(euler_qfctx, "solution time", offsetof(struct EulerContext_, curr_time), 1,
147                                                          "Physical time of the solution"));
148   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->ics.qfctx));
149   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->apply_vol_rhs.qfctx));
150   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->apply_vol_ifunction.qfctx));
151 
152   for (PetscCount b = 0; b < problem->num_bc_defs; b++) {
153     BCDefinition bc_def = problem->bc_defs[b];
154     const char  *name;
155 
156     PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL));
157     if (!strcmp(name, "outflow")) {
158       HoneeBCStruct honee_bc;
159 
160       PetscCall(PetscNew(&honee_bc));
161       PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &honee_bc->qfctx));
162       honee_bc->honee             = honee;
163       honee_bc->jac_data_size_sur = honee->phys->implicit ? problem->jac_data_size_sur : 0;
164       PetscCall(BCDefinitionSetContext(bc_def, HoneeBCDestroy, honee_bc));
165 
166       PetscCall(BCDefinitionSetIFunction(bc_def, EulerVortexOutflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp));
167       PetscCall(BCDefinitionSetIJacobian(bc_def, NULL, NULL));
168     } else if (!strcmp(name, "inflow")) {
169       HoneeBCStruct honee_bc;
170 
171       PetscCall(PetscNew(&honee_bc));
172       PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &honee_bc->qfctx));
173       honee_bc->honee             = honee;
174       honee_bc->jac_data_size_sur = honee->phys->implicit ? problem->jac_data_size_sur : 0;
175       PetscCall(BCDefinitionSetContext(bc_def, HoneeBCDestroy, honee_bc));
176 
177       PetscCall(BCDefinitionSetIFunction(bc_def, EulerVortexInflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp));
178       PetscCall(BCDefinitionSetIJacobian(bc_def, NULL, NULL));
179     }
180   }
181   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&euler_qfctx));
182   PetscFunctionReturn(PETSC_SUCCESS);
183 }
184 
185 PetscErrorCode PRINT_EULER_VORTEX(Honee honee, ProblemData problem, AppCtx app_ctx) {
186   MPI_Comm     comm = honee->comm;
187   Ceed         ceed = honee->ceed;
188   EulerContext euler_ctx;
189 
190   PetscFunctionBeginUser;
191   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfctx, CEED_MEM_HOST, &euler_ctx));
192   PetscCall(PetscPrintf(comm,
193                         "  Problem:\n"
194                         "    Problem Name                       : %s\n"
195                         "    Test Case                          : %s\n"
196                         "    Background Velocity                : %f,%f,%f\n"
197                         "    Stabilization                      : %s\n",
198                         app_ctx->problem_name, EulerTestTypes[euler_ctx->euler_test], euler_ctx->mean_velocity[0], euler_ctx->mean_velocity[1],
199                         euler_ctx->mean_velocity[2], StabilizationTypes[euler_ctx->stabilization]));
200 
201   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfctx, &euler_ctx));
202   PetscFunctionReturn(PETSC_SUCCESS);
203 }
204