xref: /honee/problems/eulervortex.c (revision 9b05e62e4a81b29e40cfcd7da7a7fab1171c8eb3)
1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3a515125bSLeila Ghaffari 
4a515125bSLeila Ghaffari /// @file
5a515125bSLeila Ghaffari /// Utility functions for setting up EULER_VORTEX
6a515125bSLeila Ghaffari 
7a515125bSLeila Ghaffari #include "../qfunctions/eulervortex.h"
8a515125bSLeila Ghaffari 
9e419654dSJeremy L Thompson #include <ceed.h>
10e419654dSJeremy L Thompson #include <petscdm.h>
11e419654dSJeremy L Thompson 
12149fb536SJames Wright #include <navierstokes.h>
136dd99bedSLeila Ghaffari 
14*9b05e62eSJames Wright static PetscErrorCode PRINT_EULER_VORTEX(Honee honee, ProblemData problem, AppCtx app_ctx) {
15*9b05e62eSJames Wright   MPI_Comm     comm = honee->comm;
16*9b05e62eSJames Wright   Ceed         ceed = honee->ceed;
17*9b05e62eSJames Wright   EulerContext euler_ctx;
18*9b05e62eSJames Wright 
19*9b05e62eSJames Wright   PetscFunctionBeginUser;
20*9b05e62eSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfctx, CEED_MEM_HOST, &euler_ctx));
21*9b05e62eSJames Wright   PetscCall(PetscPrintf(comm,
22*9b05e62eSJames Wright                         "  Problem:\n"
23*9b05e62eSJames Wright                         "    Problem Name                       : %s\n"
24*9b05e62eSJames Wright                         "    Test Case                          : %s\n"
25*9b05e62eSJames Wright                         "    Background Velocity                : %f,%f,%f\n"
26*9b05e62eSJames Wright                         "    Stabilization                      : %s\n",
27*9b05e62eSJames Wright                         app_ctx->problem_name, EulerTestTypes[euler_ctx->euler_test], euler_ctx->mean_velocity[0], euler_ctx->mean_velocity[1],
28*9b05e62eSJames Wright                         euler_ctx->mean_velocity[2], StabilizationTypes[euler_ctx->stabilization]));
29*9b05e62eSJames Wright 
30*9b05e62eSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfctx, &euler_ctx));
31*9b05e62eSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
32*9b05e62eSJames Wright }
33*9b05e62eSJames Wright 
34f978755dSJames Wright static PetscErrorCode EulerVortexOutflowBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) {
35f978755dSJames Wright   HoneeBCStruct honee_bc;
36f978755dSJames Wright 
37f978755dSJames Wright   PetscFunctionBeginUser;
38f978755dSJames Wright   PetscCall(BCDefinitionGetContext(bc_def, &honee_bc));
39f978755dSJames Wright   PetscCall(HoneeBCCreateIFunctionQF(bc_def, Euler_Outflow, Euler_Outflow_loc, honee_bc->qfctx, qf));
40f978755dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
41f978755dSJames Wright }
42f978755dSJames Wright 
43d6cac220SJames Wright static PetscErrorCode EulerVortexInflowBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) {
44d6cac220SJames Wright   HoneeBCStruct honee_bc;
45d6cac220SJames Wright 
46d6cac220SJames Wright   PetscFunctionBeginUser;
47d6cac220SJames Wright   PetscCall(BCDefinitionGetContext(bc_def, &honee_bc));
48d6cac220SJames Wright   PetscCall(HoneeBCCreateIFunctionQF(bc_def, TravelingVortex_Inflow, TravelingVortex_Inflow_loc, honee_bc->qfctx, qf));
49d6cac220SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
50d6cac220SJames Wright }
51d6cac220SJames Wright 
52d3c60affSJames Wright PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx) {
53a515125bSLeila Ghaffari   EulerTestType        euler_test;
540c373b74SJames Wright   Honee                honee = *(Honee *)ctx;
55139613f2SLeila Ghaffari   StabilizationType    stab;
560c373b74SJames Wright   MPI_Comm             comm = honee->comm;
570c373b74SJames Wright   Ceed                 ceed = honee->ceed;
58a515125bSLeila Ghaffari   PetscBool            implicit;
5915a3537eSJed Brown   EulerContext         euler_ctx;
60e07531f7SJames Wright   CeedQFunctionContext euler_qfctx;
6166d54740SJames Wright   PetscInt             dim;
62a515125bSLeila Ghaffari 
6315a3537eSJed Brown   PetscFunctionBeginUser;
642b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &euler_ctx));
65a515125bSLeila Ghaffari 
66a515125bSLeila Ghaffari   // ------------------------------------------------------
67f978755dSJames Wright   //               SET UP EULER VORTEX
68a515125bSLeila Ghaffari   // ------------------------------------------------------
69e07531f7SJames Wright   problem->ics.qf_func_ptr                 = ICsEuler;
70e07531f7SJames Wright   problem->ics.qf_loc                      = ICsEuler_loc;
71e07531f7SJames Wright   problem->apply_vol_rhs.qf_func_ptr       = Euler;
72e07531f7SJames Wright   problem->apply_vol_rhs.qf_loc            = Euler_loc;
73e07531f7SJames Wright   problem->apply_vol_ifunction.qf_func_ptr = IFunction_Euler;
74e07531f7SJames Wright   problem->apply_vol_ifunction.qf_loc      = IFunction_Euler_loc;
751abc2837SJames Wright   problem->num_comps_jac_data              = 0;
7658ce1233SJames Wright   problem->compute_exact_solution_error    = PETSC_TRUE;
77a515125bSLeila Ghaffari   problem->print_info                      = PRINT_EULER_VORTEX;
78a515125bSLeila Ghaffari 
79a515125bSLeila Ghaffari   // ------------------------------------------------------
80a515125bSLeila Ghaffari   //             Create the libCEED context
81a515125bSLeila Ghaffari   // ------------------------------------------------------
82a515125bSLeila Ghaffari   CeedScalar vortex_strength = 5.;   // -
83f821ee77SLeila Ghaffari   CeedScalar c_tau           = 0.5;  // -
84d8a22b9eSJed Brown   // c_tau = 0.5 is reported as "optimal" in Hughes et al 2010
8566d54740SJames Wright   PetscReal center[3]  = {0.},         // m
863b451b28SLeila Ghaffari       mean_velocity[3] = {1., 1., 0};  // m/s
8766d54740SJames Wright   PetscReal domain_min[3], domain_max[3], domain_size[3] = {0.};
882b916ea7SJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
8966d54740SJames Wright   PetscCall(DMGetDimension(dm, &dim));
9066d54740SJames Wright   for (PetscInt i = 0; i < dim; i++) domain_size[i] = domain_max[i] - domain_min[i];
91a515125bSLeila Ghaffari 
92a515125bSLeila Ghaffari   // ------------------------------------------------------
93a515125bSLeila Ghaffari   //             Create the PETSc context
94a515125bSLeila Ghaffari   // ------------------------------------------------------
95a515125bSLeila Ghaffari   PetscScalar meter  = 1e-2;  // 1 meter in scaled length units
96a515125bSLeila Ghaffari   PetscScalar second = 1e-2;  // 1 second in scaled time units
97a515125bSLeila Ghaffari 
98a515125bSLeila Ghaffari   // ------------------------------------------------------
99a515125bSLeila Ghaffari   //              Command line Options
100a515125bSLeila Ghaffari   // ------------------------------------------------------
1011485969bSJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for EULER_VORTEX problem", NULL);
102a515125bSLeila Ghaffari   // -- Physics
1032b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-vortex_strength", "Strength of Vortex", NULL, vortex_strength, &vortex_strength, NULL));
10466d54740SJames Wright   PetscInt  n = dim;
105a515125bSLeila Ghaffari   PetscBool user_velocity;
1062b916ea7SJeremy L Thompson   PetscCall(PetscOptionsRealArray("-mean_velocity", "Background velocity vector", NULL, mean_velocity, &n, &user_velocity));
10766d54740SJames Wright   for (PetscInt i = 0; i < dim; i++) center[i] = .5 * domain_size[i];
10866d54740SJames Wright   n = dim;
1092b916ea7SJeremy L Thompson   PetscCall(PetscOptionsRealArray("-center", "Location of vortex center", NULL, center, &n, NULL));
1102b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
1112b916ea7SJeremy L Thompson   PetscCall(PetscOptionsEnum("-euler_test", "Euler test option", NULL, EulerTestTypes, (PetscEnum)(euler_test = EULER_TEST_ISENTROPIC_VORTEX),
1122b916ea7SJeremy L Thompson                              (PetscEnum *)&euler_test, NULL));
1132b916ea7SJeremy L Thompson   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
1142b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL));
115a515125bSLeila Ghaffari   // -- Units
1162b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
117a515125bSLeila Ghaffari   meter = fabs(meter);
1182b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
119a515125bSLeila Ghaffari   second = fabs(second);
120a515125bSLeila Ghaffari 
121a515125bSLeila Ghaffari   // -- Warnings
122139613f2SLeila Ghaffari   if (stab == STAB_SUPG && !implicit) {
1232b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! Use -stab supg only with -implicit\n"));
124139613f2SLeila Ghaffari   }
1252b916ea7SJeremy L Thompson   if (user_velocity && (euler_test == EULER_TEST_1 || euler_test == EULER_TEST_3)) {
1262b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! Background velocity vector for -euler_test t1 and -euler_test t3 is (0,0,0)\n"));
127a515125bSLeila Ghaffari   }
128a515125bSLeila Ghaffari 
1291485969bSJeremy L Thompson   PetscOptionsEnd();
130a515125bSLeila Ghaffari 
131a515125bSLeila Ghaffari   // ------------------------------------------------------
132a515125bSLeila Ghaffari   //           Set up the PETSc context
133a515125bSLeila Ghaffari   // ------------------------------------------------------
1340c373b74SJames Wright   honee->units->meter  = meter;
1350c373b74SJames Wright   honee->units->second = second;
136a515125bSLeila Ghaffari 
137a515125bSLeila Ghaffari   // ------------------------------------------------------
138a515125bSLeila Ghaffari   //           Set up the libCEED context
139a515125bSLeila Ghaffari   // ------------------------------------------------------
140a515125bSLeila Ghaffari   // -- Scale variables to desired units
141493642f1SJames Wright   for (PetscInt i = 0; i < 3; i++) {
1423b451b28SLeila Ghaffari     center[i] *= meter;
14305a512bdSLeila Ghaffari     domain_size[i] *= meter;
14405a512bdSLeila Ghaffari     mean_velocity[i] *= (meter / second);
1453b451b28SLeila Ghaffari   }
146a515125bSLeila Ghaffari 
147a515125bSLeila Ghaffari   // -- QFunction Context
1480c373b74SJames Wright   honee->phys->implicit       = implicit;
14915a3537eSJed Brown   euler_ctx->curr_time        = 0.;
15015a3537eSJed Brown   euler_ctx->implicit         = implicit;
15115a3537eSJed Brown   euler_ctx->euler_test       = euler_test;
15215a3537eSJed Brown   euler_ctx->center[0]        = center[0];
15315a3537eSJed Brown   euler_ctx->center[1]        = center[1];
15415a3537eSJed Brown   euler_ctx->center[2]        = center[2];
15515a3537eSJed Brown   euler_ctx->vortex_strength  = vortex_strength;
15615a3537eSJed Brown   euler_ctx->c_tau            = c_tau;
15715a3537eSJed Brown   euler_ctx->mean_velocity[0] = mean_velocity[0];
15815a3537eSJed Brown   euler_ctx->mean_velocity[1] = mean_velocity[1];
15915a3537eSJed Brown   euler_ctx->mean_velocity[2] = mean_velocity[2];
16015a3537eSJed Brown   euler_ctx->stabilization    = stab;
161a515125bSLeila Ghaffari 
1620c373b74SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &euler_qfctx));
163e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(euler_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*euler_ctx), euler_ctx));
164e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(euler_qfctx, CEED_MEM_HOST, FreeContextPetsc));
165e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(euler_qfctx, "solution time", offsetof(struct EulerContext_, curr_time), 1,
166b4c37c5cSJames Wright                                                          "Physical time of the solution"));
167e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->ics.qfctx));
168e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->apply_vol_rhs.qfctx));
169e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->apply_vol_ifunction.qfctx));
170f978755dSJames Wright 
171f978755dSJames Wright   for (PetscCount b = 0; b < problem->num_bc_defs; b++) {
172f978755dSJames Wright     BCDefinition bc_def = problem->bc_defs[b];
173f978755dSJames Wright     const char  *name;
174f978755dSJames Wright 
175f978755dSJames Wright     PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL));
176f978755dSJames Wright     if (!strcmp(name, "outflow")) {
177f978755dSJames Wright       HoneeBCStruct honee_bc;
178f978755dSJames Wright 
179f978755dSJames Wright       PetscCall(PetscNew(&honee_bc));
180f978755dSJames Wright       PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &honee_bc->qfctx));
181f978755dSJames Wright       honee_bc->honee              = honee;
1821abc2837SJames Wright       honee_bc->num_comps_jac_data = 0;
183f978755dSJames Wright       PetscCall(BCDefinitionSetContext(bc_def, HoneeBCDestroy, honee_bc));
184f978755dSJames Wright 
185f978755dSJames Wright       PetscCall(BCDefinitionSetIFunction(bc_def, EulerVortexOutflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp));
186f978755dSJames Wright       PetscCall(BCDefinitionSetIJacobian(bc_def, NULL, NULL));
187d6cac220SJames Wright     } else if (!strcmp(name, "inflow")) {
188d6cac220SJames Wright       HoneeBCStruct honee_bc;
189d6cac220SJames Wright 
190d6cac220SJames Wright       PetscCall(PetscNew(&honee_bc));
191d6cac220SJames Wright       PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &honee_bc->qfctx));
192d6cac220SJames Wright       honee_bc->honee              = honee;
1931abc2837SJames Wright       honee_bc->num_comps_jac_data = 0;
194d6cac220SJames Wright       PetscCall(BCDefinitionSetContext(bc_def, HoneeBCDestroy, honee_bc));
195d6cac220SJames Wright 
196d6cac220SJames Wright       PetscCall(BCDefinitionSetIFunction(bc_def, EulerVortexInflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp));
197d6cac220SJames Wright       PetscCall(BCDefinitionSetIJacobian(bc_def, NULL, NULL));
198f978755dSJames Wright     }
199f978755dSJames Wright   }
200e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&euler_qfctx));
201d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
202a515125bSLeila Ghaffari }
203