xref: /honee/problems/eulervortex.c (revision 2b916ea7fa53b5c2584160b9274b1b14ca18ff4f)
1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3a515125bSLeila Ghaffari //
4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5a515125bSLeila Ghaffari //
6727da7e7SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7a515125bSLeila Ghaffari 
8a515125bSLeila Ghaffari /// @file
9a515125bSLeila Ghaffari /// Utility functions for setting up EULER_VORTEX
10a515125bSLeila Ghaffari 
11a515125bSLeila Ghaffari #include "../qfunctions/eulervortex.h"
12a515125bSLeila Ghaffari 
13*2b916ea7SJeremy L Thompson #include "../navierstokes.h"
14*2b916ea7SJeremy L Thompson #include "../qfunctions/setupgeo.h"
156dd99bedSLeila Ghaffari 
16*2b916ea7SJeremy L Thompson PetscErrorCode NS_EULER_VORTEX(ProblemData *problem, DM dm, void *ctx) {
17a515125bSLeila Ghaffari   EulerTestType        euler_test;
18a515125bSLeila Ghaffari   User                 user = *(User *)ctx;
19139613f2SLeila Ghaffari   StabilizationType    stab;
20a515125bSLeila Ghaffari   MPI_Comm             comm = PETSC_COMM_WORLD;
21a515125bSLeila Ghaffari   PetscBool            implicit;
22a515125bSLeila Ghaffari   PetscBool            has_curr_time = PETSC_TRUE;
23a515125bSLeila Ghaffari   PetscBool            has_neumann   = PETSC_TRUE;
2415a3537eSJed Brown   EulerContext         euler_ctx;
2515a3537eSJed Brown   CeedQFunctionContext euler_context;
26a515125bSLeila Ghaffari 
2715a3537eSJed Brown   PetscFunctionBeginUser;
28*2b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &euler_ctx));
29a515125bSLeila Ghaffari 
30a515125bSLeila Ghaffari   // ------------------------------------------------------
31a515125bSLeila Ghaffari   //               SET UP DENSITY_CURRENT
32a515125bSLeila Ghaffari   // ------------------------------------------------------
33a515125bSLeila Ghaffari   problem->dim                               = 3;
34a515125bSLeila Ghaffari   problem->q_data_size_vol                   = 10;
35493642f1SJames Wright   problem->q_data_size_sur                   = 10;
369785fe93SJed Brown   problem->setup_vol.qfunction               = Setup;
379785fe93SJed Brown   problem->setup_vol.qfunction_loc           = Setup_loc;
389785fe93SJed Brown   problem->setup_sur.qfunction               = SetupBoundary;
399785fe93SJed Brown   problem->setup_sur.qfunction_loc           = SetupBoundary_loc;
409785fe93SJed Brown   problem->ics.qfunction                     = ICsEuler;
419785fe93SJed Brown   problem->ics.qfunction_loc                 = ICsEuler_loc;
429785fe93SJed Brown   problem->apply_vol_rhs.qfunction           = Euler;
439785fe93SJed Brown   problem->apply_vol_rhs.qfunction_loc       = Euler_loc;
449785fe93SJed Brown   problem->apply_vol_ifunction.qfunction     = IFunction_Euler;
459785fe93SJed Brown   problem->apply_vol_ifunction.qfunction_loc = IFunction_Euler_loc;
469785fe93SJed Brown   problem->apply_inflow.qfunction            = TravelingVortex_Inflow;
479785fe93SJed Brown   problem->apply_inflow.qfunction_loc        = TravelingVortex_Inflow_loc;
489785fe93SJed Brown   problem->apply_outflow.qfunction           = Euler_Outflow;
499785fe93SJed Brown   problem->apply_outflow.qfunction_loc       = Euler_Outflow_loc;
50a515125bSLeila Ghaffari   problem->bc                                = Exact_Euler;
51b7f03f12SJed Brown   problem->bc_ctx                            = euler_ctx;
52a515125bSLeila Ghaffari   problem->non_zero_time                     = PETSC_TRUE;
53a515125bSLeila Ghaffari   problem->print_info                        = PRINT_EULER_VORTEX;
54a515125bSLeila Ghaffari 
55a515125bSLeila Ghaffari   // ------------------------------------------------------
56a515125bSLeila Ghaffari   //             Create the libCEED context
57a515125bSLeila Ghaffari   // ------------------------------------------------------
58a515125bSLeila Ghaffari   CeedScalar vortex_strength = 5.;   // -
59f821ee77SLeila Ghaffari   CeedScalar c_tau           = 0.5;  // -
60d8a22b9eSJed Brown   // c_tau = 0.5 is reported as "optimal" in Hughes et al 2010
613b451b28SLeila Ghaffari   PetscReal center[3],                 // m
623b451b28SLeila Ghaffari       mean_velocity[3] = {1., 1., 0};  // m/s
6305a512bdSLeila Ghaffari   PetscReal domain_min[3], domain_max[3], domain_size[3];
64*2b916ea7SJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
65493642f1SJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
66a515125bSLeila Ghaffari 
67a515125bSLeila Ghaffari   // ------------------------------------------------------
68a515125bSLeila Ghaffari   //             Create the PETSc context
69a515125bSLeila Ghaffari   // ------------------------------------------------------
70a515125bSLeila Ghaffari   PetscScalar meter  = 1e-2;  // 1 meter in scaled length units
71a515125bSLeila Ghaffari   PetscScalar second = 1e-2;  // 1 second in scaled time units
72a515125bSLeila Ghaffari 
73a515125bSLeila Ghaffari   // ------------------------------------------------------
74a515125bSLeila Ghaffari   //              Command line Options
75a515125bSLeila Ghaffari   // ------------------------------------------------------
761485969bSJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for EULER_VORTEX problem", NULL);
77a515125bSLeila Ghaffari   // -- Physics
78*2b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-vortex_strength", "Strength of Vortex", NULL, vortex_strength, &vortex_strength, NULL));
79a515125bSLeila Ghaffari   PetscInt  n = problem->dim;
80a515125bSLeila Ghaffari   PetscBool user_velocity;
81*2b916ea7SJeremy L Thompson   PetscCall(PetscOptionsRealArray("-mean_velocity", "Background velocity vector", NULL, mean_velocity, &n, &user_velocity));
82493642f1SJames Wright   for (PetscInt i = 0; i < 3; i++) center[i] = .5 * domain_size[i];
83a515125bSLeila Ghaffari   n = problem->dim;
84*2b916ea7SJeremy L Thompson   PetscCall(PetscOptionsRealArray("-center", "Location of vortex center", NULL, center, &n, NULL));
85*2b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
86*2b916ea7SJeremy L Thompson   PetscCall(PetscOptionsEnum("-euler_test", "Euler test option", NULL, EulerTestTypes, (PetscEnum)(euler_test = EULER_TEST_ISENTROPIC_VORTEX),
87*2b916ea7SJeremy L Thompson                              (PetscEnum *)&euler_test, NULL));
88*2b916ea7SJeremy L Thompson   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
89*2b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL));
90a515125bSLeila Ghaffari   // -- Units
91*2b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
92a515125bSLeila Ghaffari   meter = fabs(meter);
93*2b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
94a515125bSLeila Ghaffari   second = fabs(second);
95a515125bSLeila Ghaffari 
96a515125bSLeila Ghaffari   // -- Warnings
97139613f2SLeila Ghaffari   if (stab == STAB_SUPG && !implicit) {
98*2b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! Use -stab supg only with -implicit\n"));
99139613f2SLeila Ghaffari   }
100*2b916ea7SJeremy L Thompson   if (user_velocity && (euler_test == EULER_TEST_1 || euler_test == EULER_TEST_3)) {
101*2b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! Background velocity vector for -euler_test t1 and -euler_test t3 is (0,0,0)\n"));
102a515125bSLeila Ghaffari   }
103a515125bSLeila Ghaffari 
1041485969bSJeremy L Thompson   PetscOptionsEnd();
105a515125bSLeila Ghaffari 
106a515125bSLeila Ghaffari   // ------------------------------------------------------
107a515125bSLeila Ghaffari   //           Set up the PETSc context
108a515125bSLeila Ghaffari   // ------------------------------------------------------
109a515125bSLeila Ghaffari   user->units->meter  = meter;
110a515125bSLeila Ghaffari   user->units->second = second;
111a515125bSLeila Ghaffari 
112a515125bSLeila Ghaffari   // ------------------------------------------------------
113a515125bSLeila Ghaffari   //           Set up the libCEED context
114a515125bSLeila Ghaffari   // ------------------------------------------------------
115a515125bSLeila Ghaffari   // -- Scale variables to desired units
116493642f1SJames Wright   for (PetscInt i = 0; i < 3; i++) {
1173b451b28SLeila Ghaffari     center[i] *= meter;
11805a512bdSLeila Ghaffari     domain_size[i] *= meter;
11905a512bdSLeila Ghaffari     mean_velocity[i] *= (meter / second);
1203b451b28SLeila Ghaffari   }
12105a512bdSLeila Ghaffari   problem->dm_scale = meter;
122a515125bSLeila Ghaffari 
123a515125bSLeila Ghaffari   // -- QFunction Context
124139613f2SLeila Ghaffari   user->phys->stab            = stab;
125a515125bSLeila Ghaffari   user->phys->euler_test      = euler_test;
126a515125bSLeila Ghaffari   user->phys->implicit        = implicit;
127a515125bSLeila Ghaffari   user->phys->has_curr_time   = has_curr_time;
128a515125bSLeila Ghaffari   user->phys->has_neumann     = has_neumann;
12915a3537eSJed Brown   euler_ctx->curr_time        = 0.;
13015a3537eSJed Brown   euler_ctx->implicit         = implicit;
13115a3537eSJed Brown   euler_ctx->euler_test       = euler_test;
13215a3537eSJed Brown   euler_ctx->center[0]        = center[0];
13315a3537eSJed Brown   euler_ctx->center[1]        = center[1];
13415a3537eSJed Brown   euler_ctx->center[2]        = center[2];
13515a3537eSJed Brown   euler_ctx->vortex_strength  = vortex_strength;
13615a3537eSJed Brown   euler_ctx->c_tau            = c_tau;
13715a3537eSJed Brown   euler_ctx->mean_velocity[0] = mean_velocity[0];
13815a3537eSJed Brown   euler_ctx->mean_velocity[1] = mean_velocity[1];
13915a3537eSJed Brown   euler_ctx->mean_velocity[2] = mean_velocity[2];
14015a3537eSJed Brown   euler_ctx->stabilization    = stab;
141a515125bSLeila Ghaffari 
14215a3537eSJed Brown   CeedQFunctionContextCreate(user->ceed, &euler_context);
143*2b916ea7SJeremy L Thompson   CeedQFunctionContextSetData(euler_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*euler_ctx), euler_ctx);
144*2b916ea7SJeremy L Thompson   CeedQFunctionContextSetDataDestroy(euler_context, CEED_MEM_HOST, FreeContextPetsc);
145*2b916ea7SJeremy L Thompson   CeedQFunctionContextRegisterDouble(euler_context, "solution time", offsetof(struct EulerContext_, curr_time), 1, "Physical time of the solution");
146*2b916ea7SJeremy L Thompson   CeedQFunctionContextReferenceCopy(euler_context, &problem->ics.qfunction_context);
147*2b916ea7SJeremy L Thompson   CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_vol_rhs.qfunction_context);
148*2b916ea7SJeremy L Thompson   CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_vol_ifunction.qfunction_context);
149*2b916ea7SJeremy L Thompson   CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_inflow.qfunction_context);
150*2b916ea7SJeremy L Thompson   CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_outflow.qfunction_context);
151a515125bSLeila Ghaffari   PetscFunctionReturn(0);
152a515125bSLeila Ghaffari }
153a515125bSLeila Ghaffari 
154*2b916ea7SJeremy L Thompson PetscErrorCode PRINT_EULER_VORTEX(ProblemData *problem, AppCtx app_ctx) {
155a515125bSLeila Ghaffari   MPI_Comm     comm = PETSC_COMM_WORLD;
15615a3537eSJed Brown   EulerContext euler_ctx;
157a515125bSLeila Ghaffari 
15815a3537eSJed Brown   PetscFunctionBeginUser;
159*2b916ea7SJeremy L Thompson   CeedQFunctionContextGetData(problem->ics.qfunction_context, CEED_MEM_HOST, &euler_ctx);
160*2b916ea7SJeremy L Thompson   PetscCall(PetscPrintf(comm,
161a515125bSLeila Ghaffari                         "  Problem:\n"
162a515125bSLeila Ghaffari                         "    Problem Name                       : %s\n"
163a515125bSLeila Ghaffari                         "    Test Case                          : %s\n"
164139613f2SLeila Ghaffari                         "    Background Velocity                : %f,%f,%f\n"
165139613f2SLeila Ghaffari                         "    Stabilization                      : %s\n",
166*2b916ea7SJeremy L Thompson                         app_ctx->problem_name, EulerTestTypes[euler_ctx->euler_test], euler_ctx->mean_velocity[0], euler_ctx->mean_velocity[1],
167*2b916ea7SJeremy L Thompson                         euler_ctx->mean_velocity[2], StabilizationTypes[euler_ctx->stabilization]));
168a515125bSLeila Ghaffari 
16915a3537eSJed Brown   CeedQFunctionContextRestoreData(problem->ics.qfunction_context, &euler_ctx);
170a515125bSLeila Ghaffari   PetscFunctionReturn(0);
171a515125bSLeila Ghaffari }
172