xref: /libCEED/examples/fluids/problems/eulervortex.c (revision de327db462e8afadb04b2fd99bdb39cfb29a1d98)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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"
176ef2784eSLeila Ghaffari 
18731c13d7SJames Wright PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx, SimpleBC bc) {
1977841947SLeila Ghaffari   EulerTestType        euler_test;
2077841947SLeila Ghaffari   User                 user = *(User *)ctx;
21e6225c47SLeila Ghaffari   StabilizationType    stab;
22a424bcd0SJames Wright   MPI_Comm             comm = user->comm;
23a424bcd0SJames Wright   Ceed                 ceed = user->ceed;
2477841947SLeila Ghaffari   PetscBool            implicit;
25841e4c73SJed Brown   EulerContext         euler_ctx;
26841e4c73SJed Brown   CeedQFunctionContext euler_context;
2777841947SLeila Ghaffari 
28841e4c73SJed Brown   PetscFunctionBeginUser;
292b730f8bSJeremy L Thompson   PetscCall(PetscCalloc1(1, &euler_ctx));
3077841947SLeila Ghaffari 
3177841947SLeila Ghaffari   // ------------------------------------------------------
3277841947SLeila Ghaffari   //               SET UP DENSITY_CURRENT
3377841947SLeila Ghaffari   // ------------------------------------------------------
3477841947SLeila Ghaffari   problem->dim                               = 3;
3591e5af17SJed Brown   problem->ics.qfunction                     = ICsEuler;
3691e5af17SJed Brown   problem->ics.qfunction_loc                 = ICsEuler_loc;
3791e5af17SJed Brown   problem->apply_vol_rhs.qfunction           = Euler;
3891e5af17SJed Brown   problem->apply_vol_rhs.qfunction_loc       = Euler_loc;
3991e5af17SJed Brown   problem->apply_vol_ifunction.qfunction     = IFunction_Euler;
4091e5af17SJed Brown   problem->apply_vol_ifunction.qfunction_loc = IFunction_Euler_loc;
4191e5af17SJed Brown   problem->apply_inflow.qfunction            = TravelingVortex_Inflow;
4291e5af17SJed Brown   problem->apply_inflow.qfunction_loc        = TravelingVortex_Inflow_loc;
4391e5af17SJed Brown   problem->apply_outflow.qfunction           = Euler_Outflow;
4491e5af17SJed Brown   problem->apply_outflow.qfunction_loc       = Euler_Outflow_loc;
45*de327db4SJames Wright   problem->compute_exact_solution_error      = PETSC_TRUE;
4677841947SLeila Ghaffari   problem->print_info                        = PRINT_EULER_VORTEX;
4777841947SLeila Ghaffari 
4877841947SLeila Ghaffari   // ------------------------------------------------------
4977841947SLeila Ghaffari   //             Create the libCEED context
5077841947SLeila Ghaffari   // ------------------------------------------------------
5177841947SLeila Ghaffari   CeedScalar vortex_strength = 5.;   // -
52504dc8e0SLeila Ghaffari   CeedScalar c_tau           = 0.5;  // -
53932417b3SJed Brown   // c_tau = 0.5 is reported as "optimal" in Hughes et al 2010
54e535bb46SLeila Ghaffari   PetscReal center[3],                 // m
55e535bb46SLeila Ghaffari       mean_velocity[3] = {1., 1., 0};  // m/s
561864f1c2SLeila Ghaffari   PetscReal domain_min[3], domain_max[3], domain_size[3];
572b730f8bSJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
58ba6664aeSJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
5977841947SLeila Ghaffari 
6077841947SLeila Ghaffari   // ------------------------------------------------------
6177841947SLeila Ghaffari   //             Create the PETSc context
6277841947SLeila Ghaffari   // ------------------------------------------------------
6377841947SLeila Ghaffari   PetscScalar meter  = 1e-2;  // 1 meter in scaled length units
6477841947SLeila Ghaffari   PetscScalar second = 1e-2;  // 1 second in scaled time units
6577841947SLeila Ghaffari 
6677841947SLeila Ghaffari   // ------------------------------------------------------
6777841947SLeila Ghaffari   //              Command line Options
6877841947SLeila Ghaffari   // ------------------------------------------------------
6967490bc6SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for EULER_VORTEX problem", NULL);
7077841947SLeila Ghaffari   // -- Physics
712b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-vortex_strength", "Strength of Vortex", NULL, vortex_strength, &vortex_strength, NULL));
7277841947SLeila Ghaffari   PetscInt  n = problem->dim;
7377841947SLeila Ghaffari   PetscBool user_velocity;
742b730f8bSJeremy L Thompson   PetscCall(PetscOptionsRealArray("-mean_velocity", "Background velocity vector", NULL, mean_velocity, &n, &user_velocity));
75ba6664aeSJames Wright   for (PetscInt i = 0; i < 3; i++) center[i] = .5 * domain_size[i];
7677841947SLeila Ghaffari   n = problem->dim;
772b730f8bSJeremy L Thompson   PetscCall(PetscOptionsRealArray("-center", "Location of vortex center", NULL, center, &n, NULL));
782b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
792b730f8bSJeremy L Thompson   PetscCall(PetscOptionsEnum("-euler_test", "Euler test option", NULL, EulerTestTypes, (PetscEnum)(euler_test = EULER_TEST_ISENTROPIC_VORTEX),
802b730f8bSJeremy L Thompson                              (PetscEnum *)&euler_test, NULL));
812b730f8bSJeremy L Thompson   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
822b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL));
8377841947SLeila Ghaffari   // -- Units
842b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
8577841947SLeila Ghaffari   meter = fabs(meter);
862b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
8777841947SLeila Ghaffari   second = fabs(second);
8877841947SLeila Ghaffari 
8977841947SLeila Ghaffari   // -- Warnings
90e6225c47SLeila Ghaffari   if (stab == STAB_SUPG && !implicit) {
912b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! Use -stab supg only with -implicit\n"));
92e6225c47SLeila Ghaffari   }
932b730f8bSJeremy L Thompson   if (user_velocity && (euler_test == EULER_TEST_1 || euler_test == EULER_TEST_3)) {
942b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! Background velocity vector for -euler_test t1 and -euler_test t3 is (0,0,0)\n"));
9577841947SLeila Ghaffari   }
9677841947SLeila Ghaffari 
9767490bc6SJeremy L Thompson   PetscOptionsEnd();
9877841947SLeila Ghaffari 
9977841947SLeila Ghaffari   // ------------------------------------------------------
10077841947SLeila Ghaffari   //           Set up the PETSc context
10177841947SLeila Ghaffari   // ------------------------------------------------------
10277841947SLeila Ghaffari   user->units->meter  = meter;
10377841947SLeila Ghaffari   user->units->second = second;
10477841947SLeila Ghaffari 
10577841947SLeila Ghaffari   // ------------------------------------------------------
10677841947SLeila Ghaffari   //           Set up the libCEED context
10777841947SLeila Ghaffari   // ------------------------------------------------------
10877841947SLeila Ghaffari   // -- Scale variables to desired units
109ba6664aeSJames Wright   for (PetscInt i = 0; i < 3; i++) {
110e535bb46SLeila Ghaffari     center[i] *= meter;
1111864f1c2SLeila Ghaffari     domain_size[i] *= meter;
1121864f1c2SLeila Ghaffari     mean_velocity[i] *= (meter / second);
113e535bb46SLeila Ghaffari   }
1141864f1c2SLeila Ghaffari   problem->dm_scale = meter;
11577841947SLeila Ghaffari 
11677841947SLeila Ghaffari   // -- QFunction Context
11777841947SLeila Ghaffari   user->phys->implicit        = implicit;
118841e4c73SJed Brown   euler_ctx->curr_time        = 0.;
119841e4c73SJed Brown   euler_ctx->implicit         = implicit;
120841e4c73SJed Brown   euler_ctx->euler_test       = euler_test;
121841e4c73SJed Brown   euler_ctx->center[0]        = center[0];
122841e4c73SJed Brown   euler_ctx->center[1]        = center[1];
123841e4c73SJed Brown   euler_ctx->center[2]        = center[2];
124841e4c73SJed Brown   euler_ctx->vortex_strength  = vortex_strength;
125841e4c73SJed Brown   euler_ctx->c_tau            = c_tau;
126841e4c73SJed Brown   euler_ctx->mean_velocity[0] = mean_velocity[0];
127841e4c73SJed Brown   euler_ctx->mean_velocity[1] = mean_velocity[1];
128841e4c73SJed Brown   euler_ctx->mean_velocity[2] = mean_velocity[2];
129841e4c73SJed Brown   euler_ctx->stabilization    = stab;
13077841947SLeila Ghaffari 
131a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &euler_context));
132a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(euler_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*euler_ctx), euler_ctx));
133a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(euler_context, CEED_MEM_HOST, FreeContextPetsc));
134a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(euler_context, "solution time", offsetof(struct EulerContext_, curr_time), 1,
135a424bcd0SJames Wright                                                          "Physical time of the solution"));
136a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_context, &problem->ics.qfunction_context));
137a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_vol_rhs.qfunction_context));
138a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_vol_ifunction.qfunction_context));
139a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_inflow.qfunction_context));
140a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_outflow.qfunction_context));
141ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
14277841947SLeila Ghaffari }
14377841947SLeila Ghaffari 
144731c13d7SJames Wright PetscErrorCode PRINT_EULER_VORTEX(User user, ProblemData problem, AppCtx app_ctx) {
145367c849eSJames Wright   MPI_Comm     comm = user->comm;
146a424bcd0SJames Wright   Ceed         ceed = user->ceed;
147841e4c73SJed Brown   EulerContext euler_ctx;
14877841947SLeila Ghaffari 
149841e4c73SJed Brown   PetscFunctionBeginUser;
150a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfunction_context, CEED_MEM_HOST, &euler_ctx));
1512b730f8bSJeremy L Thompson   PetscCall(PetscPrintf(comm,
15277841947SLeila Ghaffari                         "  Problem:\n"
15377841947SLeila Ghaffari                         "    Problem Name                       : %s\n"
15477841947SLeila Ghaffari                         "    Test Case                          : %s\n"
155e6225c47SLeila Ghaffari                         "    Background Velocity                : %f,%f,%f\n"
156e6225c47SLeila Ghaffari                         "    Stabilization                      : %s\n",
1572b730f8bSJeremy L Thompson                         app_ctx->problem_name, EulerTestTypes[euler_ctx->euler_test], euler_ctx->mean_velocity[0], euler_ctx->mean_velocity[1],
1582b730f8bSJeremy L Thompson                         euler_ctx->mean_velocity[2], StabilizationTypes[euler_ctx->stabilization]));
15977841947SLeila Ghaffari 
160a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfunction_context, &euler_ctx));
161ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
16277841947SLeila Ghaffari }
163