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