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