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 14991aef52SJames Wright PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 15a515125bSLeila Ghaffari EulerTestType euler_test; 16a515125bSLeila Ghaffari User user = *(User *)ctx; 17139613f2SLeila Ghaffari StabilizationType stab; 18b4c37c5cSJames Wright MPI_Comm comm = user->comm; 19b4c37c5cSJames Wright Ceed ceed = user->ceed; 20a515125bSLeila Ghaffari PetscBool implicit; 2115a3537eSJed Brown EulerContext euler_ctx; 22e07531f7SJames Wright CeedQFunctionContext euler_qfctx; 23*66d54740SJames Wright PetscInt dim; 24a515125bSLeila Ghaffari 2515a3537eSJed Brown PetscFunctionBeginUser; 262b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &euler_ctx)); 27a515125bSLeila Ghaffari 28a515125bSLeila Ghaffari // ------------------------------------------------------ 29a515125bSLeila Ghaffari // SET UP DENSITY_CURRENT 30a515125bSLeila Ghaffari // ------------------------------------------------------ 31e07531f7SJames Wright problem->ics.qf_func_ptr = ICsEuler; 32e07531f7SJames Wright problem->ics.qf_loc = ICsEuler_loc; 33e07531f7SJames Wright problem->apply_vol_rhs.qf_func_ptr = Euler; 34e07531f7SJames Wright problem->apply_vol_rhs.qf_loc = Euler_loc; 35e07531f7SJames Wright problem->apply_vol_ifunction.qf_func_ptr = IFunction_Euler; 36e07531f7SJames Wright problem->apply_vol_ifunction.qf_loc = IFunction_Euler_loc; 37e07531f7SJames Wright problem->apply_inflow.qf_func_ptr = TravelingVortex_Inflow; 38e07531f7SJames Wright problem->apply_inflow.qf_loc = TravelingVortex_Inflow_loc; 39e07531f7SJames Wright problem->apply_outflow.qf_func_ptr = Euler_Outflow; 40e07531f7SJames Wright problem->apply_outflow.qf_loc = Euler_Outflow_loc; 4158ce1233SJames Wright problem->compute_exact_solution_error = PETSC_TRUE; 42a515125bSLeila Ghaffari problem->print_info = PRINT_EULER_VORTEX; 43a515125bSLeila Ghaffari 44a515125bSLeila Ghaffari // ------------------------------------------------------ 45a515125bSLeila Ghaffari // Create the libCEED context 46a515125bSLeila Ghaffari // ------------------------------------------------------ 47a515125bSLeila Ghaffari CeedScalar vortex_strength = 5.; // - 48f821ee77SLeila Ghaffari CeedScalar c_tau = 0.5; // - 49d8a22b9eSJed Brown // c_tau = 0.5 is reported as "optimal" in Hughes et al 2010 50*66d54740SJames Wright PetscReal center[3] = {0.}, // m 513b451b28SLeila Ghaffari mean_velocity[3] = {1., 1., 0}; // m/s 52*66d54740SJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3] = {0.}; 532b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 54*66d54740SJames Wright PetscCall(DMGetDimension(dm, &dim)); 55*66d54740SJames Wright for (PetscInt i = 0; i < dim; i++) domain_size[i] = domain_max[i] - domain_min[i]; 56a515125bSLeila Ghaffari 57a515125bSLeila Ghaffari // ------------------------------------------------------ 58a515125bSLeila Ghaffari // Create the PETSc context 59a515125bSLeila Ghaffari // ------------------------------------------------------ 60a515125bSLeila Ghaffari PetscScalar meter = 1e-2; // 1 meter in scaled length units 61a515125bSLeila Ghaffari PetscScalar second = 1e-2; // 1 second in scaled time units 62a515125bSLeila Ghaffari 63a515125bSLeila Ghaffari // ------------------------------------------------------ 64a515125bSLeila Ghaffari // Command line Options 65a515125bSLeila Ghaffari // ------------------------------------------------------ 661485969bSJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for EULER_VORTEX problem", NULL); 67a515125bSLeila Ghaffari // -- Physics 682b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-vortex_strength", "Strength of Vortex", NULL, vortex_strength, &vortex_strength, NULL)); 69*66d54740SJames Wright PetscInt n = dim; 70a515125bSLeila Ghaffari PetscBool user_velocity; 712b916ea7SJeremy L Thompson PetscCall(PetscOptionsRealArray("-mean_velocity", "Background velocity vector", NULL, mean_velocity, &n, &user_velocity)); 72*66d54740SJames Wright for (PetscInt i = 0; i < dim; i++) center[i] = .5 * domain_size[i]; 73*66d54740SJames Wright n = dim; 742b916ea7SJeremy L Thompson PetscCall(PetscOptionsRealArray("-center", "Location of vortex center", NULL, center, &n, NULL)); 752b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 762b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-euler_test", "Euler test option", NULL, EulerTestTypes, (PetscEnum)(euler_test = EULER_TEST_ISENTROPIC_VORTEX), 772b916ea7SJeremy L Thompson (PetscEnum *)&euler_test, NULL)); 782b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 792b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 80a515125bSLeila Ghaffari // -- Units 812b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 82a515125bSLeila Ghaffari meter = fabs(meter); 832b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 84a515125bSLeila Ghaffari second = fabs(second); 85a515125bSLeila Ghaffari 86a515125bSLeila Ghaffari // -- Warnings 87139613f2SLeila Ghaffari if (stab == STAB_SUPG && !implicit) { 882b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Use -stab supg only with -implicit\n")); 89139613f2SLeila Ghaffari } 902b916ea7SJeremy L Thompson if (user_velocity && (euler_test == EULER_TEST_1 || euler_test == EULER_TEST_3)) { 912b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Background velocity vector for -euler_test t1 and -euler_test t3 is (0,0,0)\n")); 92a515125bSLeila Ghaffari } 93a515125bSLeila Ghaffari 941485969bSJeremy L Thompson PetscOptionsEnd(); 95a515125bSLeila Ghaffari 96a515125bSLeila Ghaffari // ------------------------------------------------------ 97a515125bSLeila Ghaffari // Set up the PETSc context 98a515125bSLeila Ghaffari // ------------------------------------------------------ 99a515125bSLeila Ghaffari user->units->meter = meter; 100a515125bSLeila Ghaffari user->units->second = second; 101a515125bSLeila Ghaffari 102a515125bSLeila Ghaffari // ------------------------------------------------------ 103a515125bSLeila Ghaffari // Set up the libCEED context 104a515125bSLeila Ghaffari // ------------------------------------------------------ 105a515125bSLeila Ghaffari // -- Scale variables to desired units 106493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) { 1073b451b28SLeila Ghaffari center[i] *= meter; 10805a512bdSLeila Ghaffari domain_size[i] *= meter; 10905a512bdSLeila Ghaffari mean_velocity[i] *= (meter / second); 1103b451b28SLeila Ghaffari } 11105a512bdSLeila Ghaffari problem->dm_scale = meter; 112a515125bSLeila Ghaffari 113a515125bSLeila Ghaffari // -- QFunction Context 114a515125bSLeila Ghaffari user->phys->implicit = implicit; 11515a3537eSJed Brown euler_ctx->curr_time = 0.; 11615a3537eSJed Brown euler_ctx->implicit = implicit; 11715a3537eSJed Brown euler_ctx->euler_test = euler_test; 11815a3537eSJed Brown euler_ctx->center[0] = center[0]; 11915a3537eSJed Brown euler_ctx->center[1] = center[1]; 12015a3537eSJed Brown euler_ctx->center[2] = center[2]; 12115a3537eSJed Brown euler_ctx->vortex_strength = vortex_strength; 12215a3537eSJed Brown euler_ctx->c_tau = c_tau; 12315a3537eSJed Brown euler_ctx->mean_velocity[0] = mean_velocity[0]; 12415a3537eSJed Brown euler_ctx->mean_velocity[1] = mean_velocity[1]; 12515a3537eSJed Brown euler_ctx->mean_velocity[2] = mean_velocity[2]; 12615a3537eSJed Brown euler_ctx->stabilization = stab; 127a515125bSLeila Ghaffari 128e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &euler_qfctx)); 129e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(euler_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*euler_ctx), euler_ctx)); 130e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(euler_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 131e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(euler_qfctx, "solution time", offsetof(struct EulerContext_, curr_time), 1, 132b4c37c5cSJames Wright "Physical time of the solution")); 133e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->ics.qfctx)); 134e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->apply_vol_rhs.qfctx)); 135e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->apply_vol_ifunction.qfctx)); 136e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->apply_inflow.qfctx)); 137e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->apply_outflow.qfctx)); 138e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&euler_qfctx)); 139d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 140a515125bSLeila Ghaffari } 141a515125bSLeila Ghaffari 142991aef52SJames Wright PetscErrorCode PRINT_EULER_VORTEX(User user, ProblemData problem, AppCtx app_ctx) { 1432d49c0afSJames Wright MPI_Comm comm = user->comm; 144b4c37c5cSJames Wright Ceed ceed = user->ceed; 14515a3537eSJed Brown EulerContext euler_ctx; 146a515125bSLeila Ghaffari 14715a3537eSJed Brown PetscFunctionBeginUser; 148e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfctx, CEED_MEM_HOST, &euler_ctx)); 1492b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, 150a515125bSLeila Ghaffari " Problem:\n" 151a515125bSLeila Ghaffari " Problem Name : %s\n" 152a515125bSLeila Ghaffari " Test Case : %s\n" 153139613f2SLeila Ghaffari " Background Velocity : %f,%f,%f\n" 154139613f2SLeila Ghaffari " Stabilization : %s\n", 1552b916ea7SJeremy L Thompson app_ctx->problem_name, EulerTestTypes[euler_ctx->euler_test], euler_ctx->mean_velocity[0], euler_ctx->mean_velocity[1], 1562b916ea7SJeremy L Thompson euler_ctx->mean_velocity[2], StabilizationTypes[euler_ctx->stabilization])); 157a515125bSLeila Ghaffari 158e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfctx, &euler_ctx)); 159d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 160a515125bSLeila Ghaffari } 161