1dc936754SJeremy L Thompson // Copyright (c) 2017-2024, 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 13e419654dSJeremy L Thompson #include <ceed.h> 14e419654dSJeremy L Thompson #include <petscdm.h> 15e419654dSJeremy L Thompson 16149fb536SJames Wright #include <navierstokes.h> 176dd99bedSLeila Ghaffari 18991aef52SJames Wright PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 19a515125bSLeila Ghaffari EulerTestType euler_test; 20a515125bSLeila Ghaffari User user = *(User *)ctx; 21139613f2SLeila Ghaffari StabilizationType stab; 22b4c37c5cSJames Wright MPI_Comm comm = user->comm; 23b4c37c5cSJames Wright Ceed ceed = user->ceed; 24a515125bSLeila Ghaffari PetscBool implicit; 2515a3537eSJed Brown EulerContext euler_ctx; 2615a3537eSJed Brown CeedQFunctionContext euler_context; 27a515125bSLeila Ghaffari 2815a3537eSJed Brown PetscFunctionBeginUser; 292b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &euler_ctx)); 30a515125bSLeila Ghaffari 31a515125bSLeila Ghaffari // ------------------------------------------------------ 32a515125bSLeila Ghaffari // SET UP DENSITY_CURRENT 33a515125bSLeila Ghaffari // ------------------------------------------------------ 34a515125bSLeila Ghaffari problem->dim = 3; 359785fe93SJed Brown problem->ics.qfunction = ICsEuler; 369785fe93SJed Brown problem->ics.qfunction_loc = ICsEuler_loc; 379785fe93SJed Brown problem->apply_vol_rhs.qfunction = Euler; 389785fe93SJed Brown problem->apply_vol_rhs.qfunction_loc = Euler_loc; 399785fe93SJed Brown problem->apply_vol_ifunction.qfunction = IFunction_Euler; 409785fe93SJed Brown problem->apply_vol_ifunction.qfunction_loc = IFunction_Euler_loc; 419785fe93SJed Brown problem->apply_inflow.qfunction = TravelingVortex_Inflow; 429785fe93SJed Brown problem->apply_inflow.qfunction_loc = TravelingVortex_Inflow_loc; 439785fe93SJed Brown problem->apply_outflow.qfunction = Euler_Outflow; 449785fe93SJed Brown problem->apply_outflow.qfunction_loc = Euler_Outflow_loc; 4558ce1233SJames Wright problem->compute_exact_solution_error = PETSC_TRUE; 46a515125bSLeila Ghaffari problem->print_info = PRINT_EULER_VORTEX; 47a515125bSLeila Ghaffari 48a515125bSLeila Ghaffari // ------------------------------------------------------ 49a515125bSLeila Ghaffari // Create the libCEED context 50a515125bSLeila Ghaffari // ------------------------------------------------------ 51a515125bSLeila Ghaffari CeedScalar vortex_strength = 5.; // - 52f821ee77SLeila Ghaffari CeedScalar c_tau = 0.5; // - 53d8a22b9eSJed Brown // c_tau = 0.5 is reported as "optimal" in Hughes et al 2010 543b451b28SLeila Ghaffari PetscReal center[3], // m 553b451b28SLeila Ghaffari mean_velocity[3] = {1., 1., 0}; // m/s 5605a512bdSLeila Ghaffari PetscReal domain_min[3], domain_max[3], domain_size[3]; 572b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 58493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 59a515125bSLeila Ghaffari 60a515125bSLeila Ghaffari // ------------------------------------------------------ 61a515125bSLeila Ghaffari // Create the PETSc context 62a515125bSLeila Ghaffari // ------------------------------------------------------ 63a515125bSLeila Ghaffari PetscScalar meter = 1e-2; // 1 meter in scaled length units 64a515125bSLeila Ghaffari PetscScalar second = 1e-2; // 1 second in scaled time units 65a515125bSLeila Ghaffari 66a515125bSLeila Ghaffari // ------------------------------------------------------ 67a515125bSLeila Ghaffari // Command line Options 68a515125bSLeila Ghaffari // ------------------------------------------------------ 691485969bSJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for EULER_VORTEX problem", NULL); 70a515125bSLeila Ghaffari // -- Physics 712b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-vortex_strength", "Strength of Vortex", NULL, vortex_strength, &vortex_strength, NULL)); 72a515125bSLeila Ghaffari PetscInt n = problem->dim; 73a515125bSLeila Ghaffari PetscBool user_velocity; 742b916ea7SJeremy L Thompson PetscCall(PetscOptionsRealArray("-mean_velocity", "Background velocity vector", NULL, mean_velocity, &n, &user_velocity)); 75493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) center[i] = .5 * domain_size[i]; 76a515125bSLeila Ghaffari n = problem->dim; 772b916ea7SJeremy L Thompson PetscCall(PetscOptionsRealArray("-center", "Location of vortex center", NULL, center, &n, NULL)); 782b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 792b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-euler_test", "Euler test option", NULL, EulerTestTypes, (PetscEnum)(euler_test = EULER_TEST_ISENTROPIC_VORTEX), 802b916ea7SJeremy L Thompson (PetscEnum *)&euler_test, NULL)); 812b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 822b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 83a515125bSLeila Ghaffari // -- Units 842b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 85a515125bSLeila Ghaffari meter = fabs(meter); 862b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 87a515125bSLeila Ghaffari second = fabs(second); 88a515125bSLeila Ghaffari 89a515125bSLeila Ghaffari // -- Warnings 90139613f2SLeila Ghaffari if (stab == STAB_SUPG && !implicit) { 912b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Use -stab supg only with -implicit\n")); 92139613f2SLeila Ghaffari } 932b916ea7SJeremy L Thompson if (user_velocity && (euler_test == EULER_TEST_1 || euler_test == EULER_TEST_3)) { 942b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Background velocity vector for -euler_test t1 and -euler_test t3 is (0,0,0)\n")); 95a515125bSLeila Ghaffari } 96a515125bSLeila Ghaffari 971485969bSJeremy L Thompson PetscOptionsEnd(); 98a515125bSLeila Ghaffari 99a515125bSLeila Ghaffari // ------------------------------------------------------ 100a515125bSLeila Ghaffari // Set up the PETSc context 101a515125bSLeila Ghaffari // ------------------------------------------------------ 102a515125bSLeila Ghaffari user->units->meter = meter; 103a515125bSLeila Ghaffari user->units->second = second; 104a515125bSLeila Ghaffari 105a515125bSLeila Ghaffari // ------------------------------------------------------ 106a515125bSLeila Ghaffari // Set up the libCEED context 107a515125bSLeila Ghaffari // ------------------------------------------------------ 108a515125bSLeila Ghaffari // -- Scale variables to desired units 109493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) { 1103b451b28SLeila Ghaffari center[i] *= meter; 11105a512bdSLeila Ghaffari domain_size[i] *= meter; 11205a512bdSLeila Ghaffari mean_velocity[i] *= (meter / second); 1133b451b28SLeila Ghaffari } 11405a512bdSLeila Ghaffari problem->dm_scale = meter; 115a515125bSLeila Ghaffari 116a515125bSLeila Ghaffari // -- QFunction Context 117a515125bSLeila Ghaffari user->phys->implicit = implicit; 11815a3537eSJed Brown euler_ctx->curr_time = 0.; 11915a3537eSJed Brown euler_ctx->implicit = implicit; 12015a3537eSJed Brown euler_ctx->euler_test = euler_test; 12115a3537eSJed Brown euler_ctx->center[0] = center[0]; 12215a3537eSJed Brown euler_ctx->center[1] = center[1]; 12315a3537eSJed Brown euler_ctx->center[2] = center[2]; 12415a3537eSJed Brown euler_ctx->vortex_strength = vortex_strength; 12515a3537eSJed Brown euler_ctx->c_tau = c_tau; 12615a3537eSJed Brown euler_ctx->mean_velocity[0] = mean_velocity[0]; 12715a3537eSJed Brown euler_ctx->mean_velocity[1] = mean_velocity[1]; 12815a3537eSJed Brown euler_ctx->mean_velocity[2] = mean_velocity[2]; 12915a3537eSJed Brown euler_ctx->stabilization = stab; 130a515125bSLeila Ghaffari 131b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &euler_context)); 132b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(euler_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*euler_ctx), euler_ctx)); 133b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(euler_context, CEED_MEM_HOST, FreeContextPetsc)); 134b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(euler_context, "solution time", offsetof(struct EulerContext_, curr_time), 1, 135b4c37c5cSJames Wright "Physical time of the solution")); 136b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_context, &problem->ics.qfunction_context)); 137b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_vol_rhs.qfunction_context)); 138b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_vol_ifunction.qfunction_context)); 139b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_inflow.qfunction_context)); 140b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_context, &problem->apply_outflow.qfunction_context)); 141*19e95c73SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&euler_context)); 142d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 143a515125bSLeila Ghaffari } 144a515125bSLeila Ghaffari 145991aef52SJames Wright PetscErrorCode PRINT_EULER_VORTEX(User user, ProblemData problem, AppCtx app_ctx) { 1462d49c0afSJames Wright MPI_Comm comm = user->comm; 147b4c37c5cSJames Wright Ceed ceed = user->ceed; 14815a3537eSJed Brown EulerContext euler_ctx; 149a515125bSLeila Ghaffari 15015a3537eSJed Brown PetscFunctionBeginUser; 151b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfunction_context, CEED_MEM_HOST, &euler_ctx)); 1522b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, 153a515125bSLeila Ghaffari " Problem:\n" 154a515125bSLeila Ghaffari " Problem Name : %s\n" 155a515125bSLeila Ghaffari " Test Case : %s\n" 156139613f2SLeila Ghaffari " Background Velocity : %f,%f,%f\n" 157139613f2SLeila Ghaffari " Stabilization : %s\n", 1582b916ea7SJeremy L Thompson app_ctx->problem_name, EulerTestTypes[euler_ctx->euler_test], euler_ctx->mean_velocity[0], euler_ctx->mean_velocity[1], 1592b916ea7SJeremy L Thompson euler_ctx->mean_velocity[2], StabilizationTypes[euler_ctx->stabilization])); 160a515125bSLeila Ghaffari 161b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfunction_context, &euler_ctx)); 162d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 163a515125bSLeila Ghaffari } 164