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