xref: /libCEED/examples/fluids/problems/eulervortex.c (revision a424bcd0dd58ea8f8b80ddde5af211268f0b358a)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
377841947SLeila Ghaffari //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
577841947SLeila Ghaffari //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
777841947SLeila Ghaffari 
877841947SLeila Ghaffari /// @file
977841947SLeila Ghaffari /// Utility functions for setting up EULER_VORTEX
1077841947SLeila Ghaffari 
1177841947SLeila Ghaffari #include "../qfunctions/eulervortex.h"
1277841947SLeila Ghaffari 
1349aac155SJeremy L Thompson #include <ceed.h>
1449aac155SJeremy L Thompson #include <petscdm.h>
1549aac155SJeremy L Thompson 
162b730f8bSJeremy L Thompson #include "../navierstokes.h"
172b730f8bSJeremy L Thompson #include "../qfunctions/setupgeo.h"
186ef2784eSLeila Ghaffari 
1946de7363SJames Wright PetscErrorCode NS_EULER_VORTEX(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) {
2077841947SLeila Ghaffari   EulerTestType        euler_test;
2177841947SLeila Ghaffari   User                 user = *(User *)ctx;
22e6225c47SLeila Ghaffari   StabilizationType    stab;
23*a424bcd0SJames Wright   MPI_Comm             comm = user->comm;
24*a424bcd0SJames Wright   Ceed                 ceed = user->ceed;
2577841947SLeila Ghaffari   PetscBool            implicit;
2677841947SLeila Ghaffari   PetscBool            has_curr_time = PETSC_TRUE;
2777841947SLeila Ghaffari   PetscBool            has_neumann   = PETSC_TRUE;
28841e4c73SJed Brown   EulerContext         euler_ctx;
29841e4c73SJed Brown   CeedQFunctionContext euler_context;
3077841947SLeila Ghaffari 
31841e4c73SJed Brown   PetscFunctionBeginUser;
322b730f8bSJeremy L Thompson   PetscCall(PetscCalloc1(1, &euler_ctx));
3377841947SLeila Ghaffari 
3477841947SLeila Ghaffari   // ------------------------------------------------------
3577841947SLeila Ghaffari   //               SET UP DENSITY_CURRENT
3677841947SLeila Ghaffari   // ------------------------------------------------------
3777841947SLeila Ghaffari   problem->dim                               = 3;
3877841947SLeila Ghaffari   problem->q_data_size_vol                   = 10;
39ba6664aeSJames Wright   problem->q_data_size_sur                   = 10;
4091e5af17SJed Brown   problem->setup_vol.qfunction               = Setup;
4191e5af17SJed Brown   problem->setup_vol.qfunction_loc           = Setup_loc;
4291e5af17SJed Brown   problem->setup_sur.qfunction               = SetupBoundary;
4391e5af17SJed Brown   problem->setup_sur.qfunction_loc           = SetupBoundary_loc;
4491e5af17SJed Brown   problem->ics.qfunction                     = ICsEuler;
4591e5af17SJed Brown   problem->ics.qfunction_loc                 = ICsEuler_loc;
4691e5af17SJed Brown   problem->apply_vol_rhs.qfunction           = Euler;
4791e5af17SJed Brown   problem->apply_vol_rhs.qfunction_loc       = Euler_loc;
4891e5af17SJed Brown   problem->apply_vol_ifunction.qfunction     = IFunction_Euler;
4991e5af17SJed Brown   problem->apply_vol_ifunction.qfunction_loc = IFunction_Euler_loc;
5091e5af17SJed Brown   problem->apply_inflow.qfunction            = TravelingVortex_Inflow;
5191e5af17SJed Brown   problem->apply_inflow.qfunction_loc        = TravelingVortex_Inflow_loc;
5291e5af17SJed Brown   problem->apply_outflow.qfunction           = Euler_Outflow;
5391e5af17SJed Brown   problem->apply_outflow.qfunction_loc       = Euler_Outflow_loc;
5477841947SLeila Ghaffari   problem->non_zero_time                     = PETSC_TRUE;
5577841947SLeila Ghaffari   problem->print_info                        = PRINT_EULER_VORTEX;
5677841947SLeila Ghaffari 
5777841947SLeila Ghaffari   // ------------------------------------------------------
5877841947SLeila Ghaffari   //             Create the libCEED context
5977841947SLeila Ghaffari   // ------------------------------------------------------
6077841947SLeila Ghaffari   CeedScalar vortex_strength = 5.;   // -
61504dc8e0SLeila Ghaffari   CeedScalar c_tau           = 0.5;  // -
62932417b3SJed Brown   // c_tau = 0.5 is reported as "optimal" in Hughes et al 2010
63e535bb46SLeila Ghaffari   PetscReal center[3],                 // m
64e535bb46SLeila Ghaffari       mean_velocity[3] = {1., 1., 0};  // m/s
651864f1c2SLeila Ghaffari   PetscReal domain_min[3], domain_max[3], domain_size[3];
662b730f8bSJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
67ba6664aeSJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
6877841947SLeila Ghaffari 
6977841947SLeila Ghaffari   // ------------------------------------------------------
7077841947SLeila Ghaffari   //             Create the PETSc context
7177841947SLeila Ghaffari   // ------------------------------------------------------
7277841947SLeila Ghaffari   PetscScalar meter  = 1e-2;  // 1 meter in scaled length units
7377841947SLeila Ghaffari   PetscScalar second = 1e-2;  // 1 second in scaled time units
7477841947SLeila Ghaffari 
7577841947SLeila Ghaffari   // ------------------------------------------------------
7677841947SLeila Ghaffari   //              Command line Options
7777841947SLeila Ghaffari   // ------------------------------------------------------
7867490bc6SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for EULER_VORTEX problem", NULL);
7977841947SLeila Ghaffari   // -- Physics
802b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-vortex_strength", "Strength of Vortex", NULL, vortex_strength, &vortex_strength, NULL));
8177841947SLeila Ghaffari   PetscInt  n = problem->dim;
8277841947SLeila Ghaffari   PetscBool user_velocity;
832b730f8bSJeremy L Thompson   PetscCall(PetscOptionsRealArray("-mean_velocity", "Background velocity vector", NULL, mean_velocity, &n, &user_velocity));
84ba6664aeSJames Wright   for (PetscInt i = 0; i < 3; i++) center[i] = .5 * domain_size[i];
8577841947SLeila Ghaffari   n = problem->dim;
862b730f8bSJeremy L Thompson   PetscCall(PetscOptionsRealArray("-center", "Location of vortex center", NULL, center, &n, NULL));
872b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
882b730f8bSJeremy L Thompson   PetscCall(PetscOptionsEnum("-euler_test", "Euler test option", NULL, EulerTestTypes, (PetscEnum)(euler_test = EULER_TEST_ISENTROPIC_VORTEX),
892b730f8bSJeremy L Thompson                              (PetscEnum *)&euler_test, NULL));
902b730f8bSJeremy L Thompson   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
912b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL));
9277841947SLeila Ghaffari   // -- Units
932b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
9477841947SLeila Ghaffari   meter = fabs(meter);
952b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
9677841947SLeila Ghaffari   second = fabs(second);
9777841947SLeila Ghaffari 
9877841947SLeila Ghaffari   // -- Warnings
99e6225c47SLeila Ghaffari   if (stab == STAB_SUPG && !implicit) {
1002b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! Use -stab supg only with -implicit\n"));
101e6225c47SLeila Ghaffari   }
1022b730f8bSJeremy L Thompson   if (user_velocity && (euler_test == EULER_TEST_1 || euler_test == EULER_TEST_3)) {
1032b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! Background velocity vector for -euler_test t1 and -euler_test t3 is (0,0,0)\n"));
10477841947SLeila Ghaffari   }
10577841947SLeila Ghaffari 
10667490bc6SJeremy L Thompson   PetscOptionsEnd();
10777841947SLeila Ghaffari 
10877841947SLeila Ghaffari   // ------------------------------------------------------
10977841947SLeila Ghaffari   //           Set up the PETSc context
11077841947SLeila Ghaffari   // ------------------------------------------------------
11177841947SLeila Ghaffari   user->units->meter  = meter;
11277841947SLeila Ghaffari   user->units->second = second;
11377841947SLeila Ghaffari 
11477841947SLeila Ghaffari   // ------------------------------------------------------
11577841947SLeila Ghaffari   //           Set up the libCEED context
11677841947SLeila Ghaffari   // ------------------------------------------------------
11777841947SLeila Ghaffari   // -- Scale variables to desired units
118ba6664aeSJames Wright   for (PetscInt i = 0; i < 3; i++) {
119e535bb46SLeila Ghaffari     center[i] *= meter;
1201864f1c2SLeila Ghaffari     domain_size[i] *= meter;
1211864f1c2SLeila Ghaffari     mean_velocity[i] *= (meter / second);
122e535bb46SLeila Ghaffari   }
1231864f1c2SLeila Ghaffari   problem->dm_scale = meter;
12477841947SLeila Ghaffari 
12577841947SLeila Ghaffari   // -- QFunction Context
126e6225c47SLeila Ghaffari   user->phys->stab            = stab;
12777841947SLeila Ghaffari   user->phys->euler_test      = euler_test;
12877841947SLeila Ghaffari   user->phys->implicit        = implicit;
12977841947SLeila Ghaffari   user->phys->has_curr_time   = has_curr_time;
13077841947SLeila Ghaffari   user->phys->has_neumann     = has_neumann;
131841e4c73SJed Brown   euler_ctx->curr_time        = 0.;
132841e4c73SJed Brown   euler_ctx->implicit         = implicit;
133841e4c73SJed Brown   euler_ctx->euler_test       = euler_test;
134841e4c73SJed Brown   euler_ctx->center[0]        = center[0];
135841e4c73SJed Brown   euler_ctx->center[1]        = center[1];
136841e4c73SJed Brown   euler_ctx->center[2]        = center[2];
137841e4c73SJed Brown   euler_ctx->vortex_strength  = vortex_strength;
138841e4c73SJed Brown   euler_ctx->c_tau            = c_tau;
139841e4c73SJed Brown   euler_ctx->mean_velocity[0] = mean_velocity[0];
140841e4c73SJed Brown   euler_ctx->mean_velocity[1] = mean_velocity[1];
141841e4c73SJed Brown   euler_ctx->mean_velocity[2] = mean_velocity[2];
142841e4c73SJed Brown   euler_ctx->stabilization    = stab;
14377841947SLeila Ghaffari 
144*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &euler_context));
145*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(euler_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*euler_ctx), euler_ctx));
146*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(euler_context, CEED_MEM_HOST, FreeContextPetsc));
147*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(euler_context, "solution time", offsetof(struct EulerContext_, curr_time), 1,
148*a424bcd0SJames Wright                                                          "Physical time of the solution"));
149*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_context, &problem->ics.qfunction_context));
150*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_vol_rhs.qfunction_context));
151*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_vol_ifunction.qfunction_context));
152*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_inflow.qfunction_context));
153*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_outflow.qfunction_context));
154ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
15577841947SLeila Ghaffari }
15677841947SLeila Ghaffari 
157367c849eSJames Wright PetscErrorCode PRINT_EULER_VORTEX(User user, ProblemData *problem, AppCtx app_ctx) {
158367c849eSJames Wright   MPI_Comm     comm = user->comm;
159*a424bcd0SJames Wright   Ceed         ceed = user->ceed;
160841e4c73SJed Brown   EulerContext euler_ctx;
16177841947SLeila Ghaffari 
162841e4c73SJed Brown   PetscFunctionBeginUser;
163*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfunction_context, CEED_MEM_HOST, &euler_ctx));
1642b730f8bSJeremy L Thompson   PetscCall(PetscPrintf(comm,
16577841947SLeila Ghaffari                         "  Problem:\n"
16677841947SLeila Ghaffari                         "    Problem Name                       : %s\n"
16777841947SLeila Ghaffari                         "    Test Case                          : %s\n"
168e6225c47SLeila Ghaffari                         "    Background Velocity                : %f,%f,%f\n"
169e6225c47SLeila Ghaffari                         "    Stabilization                      : %s\n",
1702b730f8bSJeremy L Thompson                         app_ctx->problem_name, EulerTestTypes[euler_ctx->euler_test], euler_ctx->mean_velocity[0], euler_ctx->mean_velocity[1],
1712b730f8bSJeremy L Thompson                         euler_ctx->mean_velocity[2], StabilizationTypes[euler_ctx->stabilization]));
17277841947SLeila Ghaffari 
173*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfunction_context, &euler_ctx));
174ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
17577841947SLeila Ghaffari }
176