xref: /honee/problems/eulervortex.c (revision 337840fc9f8b699b98947f09e30443f8c0c7f962)
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 
PRINT_EULER_VORTEX(Honee honee,ProblemData problem,AppCtx app_ctx)149b05e62eSJames Wright static PetscErrorCode PRINT_EULER_VORTEX(Honee honee, ProblemData problem, AppCtx app_ctx) {
159b05e62eSJames Wright   MPI_Comm     comm = honee->comm;
169b05e62eSJames Wright   Ceed         ceed = honee->ceed;
179b05e62eSJames Wright   EulerContext euler_ctx;
189b05e62eSJames Wright 
199b05e62eSJames Wright   PetscFunctionBeginUser;
209b05e62eSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfctx, CEED_MEM_HOST, &euler_ctx));
219b05e62eSJames Wright   PetscCall(PetscPrintf(comm,
229b05e62eSJames Wright                         "  Problem:\n"
239b05e62eSJames Wright                         "    Problem Name                       : %s\n"
249b05e62eSJames Wright                         "    Test Case                          : %s\n"
259b05e62eSJames Wright                         "    Background Velocity                : %f,%f,%f\n"
269b05e62eSJames Wright                         "    Stabilization                      : %s\n",
279b05e62eSJames Wright                         app_ctx->problem_name, EulerTestTypes[euler_ctx->euler_test], euler_ctx->mean_velocity[0], euler_ctx->mean_velocity[1],
289b05e62eSJames Wright                         euler_ctx->mean_velocity[2], StabilizationTypes[euler_ctx->stabilization]));
299b05e62eSJames Wright 
309b05e62eSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfctx, &euler_ctx));
319b05e62eSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
329b05e62eSJames Wright }
339b05e62eSJames Wright 
EulerVortexOutflowBCSetup_CreateIFunctionQF(BCDefinition bc_def,CeedQFunction * qf)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 
EulerVortexInflowBCSetup_CreateIFunctionQF(BCDefinition bc_def,CeedQFunction * qf)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 
52f27c2204SJames Wright static const char *const component_names[] = {"Density", "MomentumX", "MomentumY", "MomentumZ", "TotalEnergy"};
53f27c2204SJames Wright 
NS_EULER_VORTEX(ProblemData problem,DM dm,void * ctx)54d3c60affSJames Wright PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx) {
55a515125bSLeila Ghaffari   EulerTestType        euler_test;
560c373b74SJames Wright   Honee                honee = *(Honee *)ctx;
57139613f2SLeila Ghaffari   StabilizationType    stab;
580c373b74SJames Wright   MPI_Comm             comm = honee->comm;
590c373b74SJames Wright   Ceed                 ceed = honee->ceed;
60a515125bSLeila Ghaffari   PetscBool            implicit;
6115a3537eSJed Brown   EulerContext         euler_ctx;
62e07531f7SJames Wright   CeedQFunctionContext euler_qfctx;
6366d54740SJames Wright   PetscInt             dim;
64a515125bSLeila Ghaffari 
6515a3537eSJed Brown   PetscFunctionBeginUser;
662d898fa6SJames Wright   PetscCall(PetscNew(&euler_ctx));
67a515125bSLeila Ghaffari 
68a515125bSLeila Ghaffari   // ------------------------------------------------------
69f978755dSJames Wright   //               SET UP EULER VORTEX
70a515125bSLeila Ghaffari   // ------------------------------------------------------
71f5dc303cSJames Wright   problem->ics                          = (HoneeQFSpec){.qf_func_ptr = ICsEuler, .qf_loc = ICsEuler_loc};
72f5dc303cSJames Wright   problem->apply_vol_rhs                = (HoneeQFSpec){.qf_func_ptr = Euler, .qf_loc = Euler_loc};
73f5dc303cSJames Wright   problem->apply_vol_ifunction          = (HoneeQFSpec){.qf_func_ptr = IFunction_Euler, .qf_loc = IFunction_Euler_loc};
741abc2837SJames Wright   problem->num_comps_jac_data           = 0;
7558ce1233SJames Wright   problem->compute_exact_solution_error = PETSC_TRUE;
76a515125bSLeila Ghaffari   problem->print_info                   = PRINT_EULER_VORTEX;
77a515125bSLeila Ghaffari 
78f27c2204SJames Wright   problem->num_components = 5;
79f27c2204SJames Wright   PetscCall(PetscMalloc1(problem->num_components, &problem->component_names));
80f27c2204SJames Wright   for (PetscInt i = 0; i < 5; i++) PetscCall(PetscStrallocpy(component_names[i], &problem->component_names[i]));
81f27c2204SJames Wright 
82a515125bSLeila Ghaffari   // ------------------------------------------------------
83a515125bSLeila Ghaffari   //             Create the libCEED context
84a515125bSLeila Ghaffari   // ------------------------------------------------------
857e3656bdSJames Wright   Units      units           = honee->units;
86a515125bSLeila Ghaffari   CeedScalar vortex_strength = 5.;   // -
87f821ee77SLeila Ghaffari   CeedScalar c_tau           = 0.5;  // -
88d8a22b9eSJed Brown   // c_tau = 0.5 is reported as "optimal" in Hughes et al 2010
8966d54740SJames Wright   PetscReal center[3]  = {0.},         // m
903b451b28SLeila Ghaffari       mean_velocity[3] = {1., 1., 0};  // m/s
9166d54740SJames Wright   PetscReal domain_min[3], domain_max[3], domain_size[3] = {0.};
922b916ea7SJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
9366d54740SJames Wright   PetscCall(DMGetDimension(dm, &dim));
9466d54740SJames Wright   for (PetscInt i = 0; i < dim; i++) domain_size[i] = domain_max[i] - domain_min[i];
95a515125bSLeila Ghaffari 
96a515125bSLeila Ghaffari   // ------------------------------------------------------
97a515125bSLeila Ghaffari   //              Command line Options
98a515125bSLeila Ghaffari   // ------------------------------------------------------
991485969bSJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for EULER_VORTEX problem", NULL);
100a515125bSLeila Ghaffari   // -- Physics
1012b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-vortex_strength", "Strength of Vortex", NULL, vortex_strength, &vortex_strength, NULL));
10266d54740SJames Wright   PetscInt  n = dim;
103a515125bSLeila Ghaffari   PetscBool user_velocity;
1042b916ea7SJeremy L Thompson   PetscCall(PetscOptionsRealArray("-mean_velocity", "Background velocity vector", NULL, mean_velocity, &n, &user_velocity));
1057e3656bdSJames Wright   for (PetscInt i = 0; i < dim; i++) center[i] = .5 * domain_size[i] / units->meter;  // Redimensionalize domain
10666d54740SJames Wright   n = dim;
1072b916ea7SJeremy L Thompson   PetscCall(PetscOptionsRealArray("-center", "Location of vortex center", NULL, center, &n, NULL));
1082b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
1092b916ea7SJeremy L Thompson   PetscCall(PetscOptionsEnum("-euler_test", "Euler test option", NULL, EulerTestTypes, (PetscEnum)(euler_test = EULER_TEST_ISENTROPIC_VORTEX),
1102b916ea7SJeremy L Thompson                              (PetscEnum *)&euler_test, NULL));
1112b916ea7SJeremy L Thompson   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
1122b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL));
113a515125bSLeila Ghaffari 
114a515125bSLeila Ghaffari   // -- Warnings
115139613f2SLeila Ghaffari   if (stab == STAB_SUPG && !implicit) {
1162b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! Use -stab supg only with -implicit\n"));
117139613f2SLeila Ghaffari   }
1182b916ea7SJeremy L Thompson   if (user_velocity && (euler_test == EULER_TEST_1 || euler_test == EULER_TEST_3)) {
1192b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! Background velocity vector for -euler_test t1 and -euler_test t3 is (0,0,0)\n"));
120a515125bSLeila Ghaffari   }
121a515125bSLeila Ghaffari 
1221485969bSJeremy L Thompson   PetscOptionsEnd();
123a515125bSLeila Ghaffari 
124a515125bSLeila Ghaffari   // ------------------------------------------------------
125a515125bSLeila Ghaffari   //           Set up the libCEED context
126a515125bSLeila Ghaffari   // ------------------------------------------------------
127a515125bSLeila Ghaffari   // -- Scale variables to desired units
128493642f1SJames Wright   for (PetscInt i = 0; i < 3; i++) {
129c9f37605SMohammed Amin     center[i] *= units->meter;
130c9f37605SMohammed Amin     mean_velocity[i] *= (units->meter / units->second);
1313b451b28SLeila Ghaffari   }
132a515125bSLeila Ghaffari 
133a515125bSLeila Ghaffari   // -- QFunction Context
1340c373b74SJames Wright   honee->phys->implicit       = implicit;
13515a3537eSJed Brown   euler_ctx->curr_time        = 0.;
13615a3537eSJed Brown   euler_ctx->implicit         = implicit;
13715a3537eSJed Brown   euler_ctx->euler_test       = euler_test;
13815a3537eSJed Brown   euler_ctx->center[0]        = center[0];
13915a3537eSJed Brown   euler_ctx->center[1]        = center[1];
14015a3537eSJed Brown   euler_ctx->center[2]        = center[2];
14115a3537eSJed Brown   euler_ctx->vortex_strength  = vortex_strength;
14215a3537eSJed Brown   euler_ctx->c_tau            = c_tau;
14315a3537eSJed Brown   euler_ctx->mean_velocity[0] = mean_velocity[0];
14415a3537eSJed Brown   euler_ctx->mean_velocity[1] = mean_velocity[1];
14515a3537eSJed Brown   euler_ctx->mean_velocity[2] = mean_velocity[2];
14615a3537eSJed Brown   euler_ctx->stabilization    = stab;
147a515125bSLeila Ghaffari 
1480c373b74SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &euler_qfctx));
149e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(euler_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*euler_ctx), euler_ctx));
150e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(euler_qfctx, CEED_MEM_HOST, FreeContextPetsc));
151e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(euler_qfctx, "solution time", offsetof(struct EulerContext_, curr_time), 1,
152b4c37c5cSJames Wright                                                          "Physical time of the solution"));
153e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->ics.qfctx));
154e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->apply_vol_rhs.qfctx));
155e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->apply_vol_ifunction.qfctx));
156f978755dSJames Wright 
157f978755dSJames Wright   for (PetscCount b = 0; b < problem->num_bc_defs; b++) {
158f978755dSJames Wright     BCDefinition bc_def = problem->bc_defs[b];
159f978755dSJames Wright     const char  *name;
160f978755dSJames Wright 
161f978755dSJames Wright     PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL));
162f978755dSJames Wright     if (!strcmp(name, "outflow")) {
163f978755dSJames Wright       HoneeBCStruct honee_bc;
164f978755dSJames Wright 
165f978755dSJames Wright       PetscCall(PetscNew(&honee_bc));
166f978755dSJames Wright       PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &honee_bc->qfctx));
167f978755dSJames Wright       honee_bc->honee              = honee;
1681abc2837SJames Wright       honee_bc->num_comps_jac_data = 0;
169*26d401f3SJames Wright       PetscCall(BCDefinitionSetContext(bc_def, (PetscCtxDestroyFn *)HoneeBCDestroy, honee_bc));
170f978755dSJames Wright 
171f978755dSJames Wright       PetscCall(BCDefinitionSetIFunction(bc_def, EulerVortexOutflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp));
172f978755dSJames Wright       PetscCall(BCDefinitionSetIJacobian(bc_def, NULL, NULL));
173d6cac220SJames Wright     } else if (!strcmp(name, "inflow")) {
174d6cac220SJames Wright       HoneeBCStruct honee_bc;
175d6cac220SJames Wright 
176d6cac220SJames Wright       PetscCall(PetscNew(&honee_bc));
177d6cac220SJames Wright       PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &honee_bc->qfctx));
178d6cac220SJames Wright       honee_bc->honee              = honee;
1791abc2837SJames Wright       honee_bc->num_comps_jac_data = 0;
180*26d401f3SJames Wright       PetscCall(BCDefinitionSetContext(bc_def, (PetscCtxDestroyFn *)HoneeBCDestroy, honee_bc));
181d6cac220SJames Wright 
182d6cac220SJames Wright       PetscCall(BCDefinitionSetIFunction(bc_def, EulerVortexInflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp));
183d6cac220SJames Wright       PetscCall(BCDefinitionSetIJacobian(bc_def, NULL, NULL));
184f978755dSJames Wright     }
185f978755dSJames Wright   }
186e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&euler_qfctx));
187d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
188a515125bSLeila Ghaffari }
189