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 14*f978755dSJames Wright static PetscErrorCode EulerVortexOutflowBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) { 15*f978755dSJames Wright HoneeBCStruct honee_bc; 16*f978755dSJames Wright 17*f978755dSJames Wright PetscFunctionBeginUser; 18*f978755dSJames Wright PetscCall(BCDefinitionGetContext(bc_def, &honee_bc)); 19*f978755dSJames Wright PetscCall(HoneeBCCreateIFunctionQF(bc_def, Euler_Outflow, Euler_Outflow_loc, honee_bc->qfctx, qf)); 20*f978755dSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 21*f978755dSJames Wright } 22*f978755dSJames Wright 23991aef52SJames Wright PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 24a515125bSLeila Ghaffari EulerTestType euler_test; 250c373b74SJames Wright Honee honee = *(Honee *)ctx; 26139613f2SLeila Ghaffari StabilizationType stab; 270c373b74SJames Wright MPI_Comm comm = honee->comm; 280c373b74SJames Wright Ceed ceed = honee->ceed; 29a515125bSLeila Ghaffari PetscBool implicit; 3015a3537eSJed Brown EulerContext euler_ctx; 31e07531f7SJames Wright CeedQFunctionContext euler_qfctx; 3266d54740SJames Wright PetscInt dim; 33a515125bSLeila Ghaffari 3415a3537eSJed Brown PetscFunctionBeginUser; 352b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &euler_ctx)); 36a515125bSLeila Ghaffari 37a515125bSLeila Ghaffari // ------------------------------------------------------ 38*f978755dSJames Wright // SET UP EULER VORTEX 39a515125bSLeila Ghaffari // ------------------------------------------------------ 40e07531f7SJames Wright problem->ics.qf_func_ptr = ICsEuler; 41e07531f7SJames Wright problem->ics.qf_loc = ICsEuler_loc; 42e07531f7SJames Wright problem->apply_vol_rhs.qf_func_ptr = Euler; 43e07531f7SJames Wright problem->apply_vol_rhs.qf_loc = Euler_loc; 44e07531f7SJames Wright problem->apply_vol_ifunction.qf_func_ptr = IFunction_Euler; 45e07531f7SJames Wright problem->apply_vol_ifunction.qf_loc = IFunction_Euler_loc; 46e07531f7SJames Wright problem->apply_inflow.qf_func_ptr = TravelingVortex_Inflow; 47e07531f7SJames Wright problem->apply_inflow.qf_loc = TravelingVortex_Inflow_loc; 4828160fc2SJames Wright problem->jac_data_size_vol = 0; 4928160fc2SJames Wright problem->jac_data_size_sur = 0; 5058ce1233SJames Wright problem->compute_exact_solution_error = PETSC_TRUE; 51a515125bSLeila Ghaffari problem->print_info = PRINT_EULER_VORTEX; 52a515125bSLeila Ghaffari 53a515125bSLeila Ghaffari // ------------------------------------------------------ 54a515125bSLeila Ghaffari // Create the libCEED context 55a515125bSLeila Ghaffari // ------------------------------------------------------ 56a515125bSLeila Ghaffari CeedScalar vortex_strength = 5.; // - 57f821ee77SLeila Ghaffari CeedScalar c_tau = 0.5; // - 58d8a22b9eSJed Brown // c_tau = 0.5 is reported as "optimal" in Hughes et al 2010 5966d54740SJames Wright PetscReal center[3] = {0.}, // m 603b451b28SLeila Ghaffari mean_velocity[3] = {1., 1., 0}; // m/s 6166d54740SJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3] = {0.}; 622b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 6366d54740SJames Wright PetscCall(DMGetDimension(dm, &dim)); 6466d54740SJames Wright for (PetscInt i = 0; i < dim; i++) domain_size[i] = domain_max[i] - domain_min[i]; 65a515125bSLeila Ghaffari 66a515125bSLeila Ghaffari // ------------------------------------------------------ 67a515125bSLeila Ghaffari // Create the PETSc context 68a515125bSLeila Ghaffari // ------------------------------------------------------ 69a515125bSLeila Ghaffari PetscScalar meter = 1e-2; // 1 meter in scaled length units 70a515125bSLeila Ghaffari PetscScalar second = 1e-2; // 1 second in scaled time units 71a515125bSLeila Ghaffari 72a515125bSLeila Ghaffari // ------------------------------------------------------ 73a515125bSLeila Ghaffari // Command line Options 74a515125bSLeila Ghaffari // ------------------------------------------------------ 751485969bSJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for EULER_VORTEX problem", NULL); 76a515125bSLeila Ghaffari // -- Physics 772b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-vortex_strength", "Strength of Vortex", NULL, vortex_strength, &vortex_strength, NULL)); 7866d54740SJames Wright PetscInt n = dim; 79a515125bSLeila Ghaffari PetscBool user_velocity; 802b916ea7SJeremy L Thompson PetscCall(PetscOptionsRealArray("-mean_velocity", "Background velocity vector", NULL, mean_velocity, &n, &user_velocity)); 8166d54740SJames Wright for (PetscInt i = 0; i < dim; i++) center[i] = .5 * domain_size[i]; 8266d54740SJames Wright n = dim; 832b916ea7SJeremy L Thompson PetscCall(PetscOptionsRealArray("-center", "Location of vortex center", NULL, center, &n, NULL)); 842b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 852b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-euler_test", "Euler test option", NULL, EulerTestTypes, (PetscEnum)(euler_test = EULER_TEST_ISENTROPIC_VORTEX), 862b916ea7SJeremy L Thompson (PetscEnum *)&euler_test, NULL)); 872b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 882b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 89a515125bSLeila Ghaffari // -- Units 902b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 91a515125bSLeila Ghaffari meter = fabs(meter); 922b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 93a515125bSLeila Ghaffari second = fabs(second); 94a515125bSLeila Ghaffari 95a515125bSLeila Ghaffari // -- Warnings 96139613f2SLeila Ghaffari if (stab == STAB_SUPG && !implicit) { 972b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Use -stab supg only with -implicit\n")); 98139613f2SLeila Ghaffari } 992b916ea7SJeremy L Thompson if (user_velocity && (euler_test == EULER_TEST_1 || euler_test == EULER_TEST_3)) { 1002b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Background velocity vector for -euler_test t1 and -euler_test t3 is (0,0,0)\n")); 101a515125bSLeila Ghaffari } 102a515125bSLeila Ghaffari 1031485969bSJeremy L Thompson PetscOptionsEnd(); 104a515125bSLeila Ghaffari 105a515125bSLeila Ghaffari // ------------------------------------------------------ 106a515125bSLeila Ghaffari // Set up the PETSc context 107a515125bSLeila Ghaffari // ------------------------------------------------------ 1080c373b74SJames Wright honee->units->meter = meter; 1090c373b74SJames Wright honee->units->second = second; 110a515125bSLeila Ghaffari 111a515125bSLeila Ghaffari // ------------------------------------------------------ 112a515125bSLeila Ghaffari // Set up the libCEED context 113a515125bSLeila Ghaffari // ------------------------------------------------------ 114a515125bSLeila Ghaffari // -- Scale variables to desired units 115493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) { 1163b451b28SLeila Ghaffari center[i] *= meter; 11705a512bdSLeila Ghaffari domain_size[i] *= meter; 11805a512bdSLeila Ghaffari mean_velocity[i] *= (meter / second); 1193b451b28SLeila Ghaffari } 120a515125bSLeila Ghaffari 121a515125bSLeila Ghaffari // -- QFunction Context 1220c373b74SJames Wright honee->phys->implicit = implicit; 12315a3537eSJed Brown euler_ctx->curr_time = 0.; 12415a3537eSJed Brown euler_ctx->implicit = implicit; 12515a3537eSJed Brown euler_ctx->euler_test = euler_test; 12615a3537eSJed Brown euler_ctx->center[0] = center[0]; 12715a3537eSJed Brown euler_ctx->center[1] = center[1]; 12815a3537eSJed Brown euler_ctx->center[2] = center[2]; 12915a3537eSJed Brown euler_ctx->vortex_strength = vortex_strength; 13015a3537eSJed Brown euler_ctx->c_tau = c_tau; 13115a3537eSJed Brown euler_ctx->mean_velocity[0] = mean_velocity[0]; 13215a3537eSJed Brown euler_ctx->mean_velocity[1] = mean_velocity[1]; 13315a3537eSJed Brown euler_ctx->mean_velocity[2] = mean_velocity[2]; 13415a3537eSJed Brown euler_ctx->stabilization = stab; 135a515125bSLeila Ghaffari 1360c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &euler_qfctx)); 137e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(euler_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*euler_ctx), euler_ctx)); 138e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(euler_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 139e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(euler_qfctx, "solution time", offsetof(struct EulerContext_, curr_time), 1, 140b4c37c5cSJames Wright "Physical time of the solution")); 141e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->ics.qfctx)); 142e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->apply_vol_rhs.qfctx)); 143e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->apply_vol_ifunction.qfctx)); 144e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->apply_inflow.qfctx)); 145*f978755dSJames Wright 146*f978755dSJames Wright for (PetscCount b = 0; b < problem->num_bc_defs; b++) { 147*f978755dSJames Wright BCDefinition bc_def = problem->bc_defs[b]; 148*f978755dSJames Wright const char *name; 149*f978755dSJames Wright 150*f978755dSJames Wright PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL)); 151*f978755dSJames Wright if (!strcmp(name, "outflow")) { 152*f978755dSJames Wright HoneeBCStruct honee_bc; 153*f978755dSJames Wright 154*f978755dSJames Wright PetscCall(PetscNew(&honee_bc)); 155*f978755dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &honee_bc->qfctx)); 156*f978755dSJames Wright honee_bc->honee = honee; 157*f978755dSJames Wright honee_bc->jac_data_size_sur = honee->phys->implicit ? problem->jac_data_size_sur : 0; 158*f978755dSJames Wright PetscCall(BCDefinitionSetContext(bc_def, HoneeBCDestroy, honee_bc)); 159*f978755dSJames Wright 160*f978755dSJames Wright PetscCall(BCDefinitionSetIFunction(bc_def, EulerVortexOutflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp)); 161*f978755dSJames Wright PetscCall(BCDefinitionSetIJacobian(bc_def, NULL, NULL)); 162*f978755dSJames Wright } 163*f978755dSJames Wright } 164e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&euler_qfctx)); 165d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 166a515125bSLeila Ghaffari } 167a515125bSLeila Ghaffari 1680c373b74SJames Wright PetscErrorCode PRINT_EULER_VORTEX(Honee honee, ProblemData problem, AppCtx app_ctx) { 1690c373b74SJames Wright MPI_Comm comm = honee->comm; 1700c373b74SJames Wright Ceed ceed = honee->ceed; 17115a3537eSJed Brown EulerContext euler_ctx; 172a515125bSLeila Ghaffari 17315a3537eSJed Brown PetscFunctionBeginUser; 174e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfctx, CEED_MEM_HOST, &euler_ctx)); 1752b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, 176a515125bSLeila Ghaffari " Problem:\n" 177a515125bSLeila Ghaffari " Problem Name : %s\n" 178a515125bSLeila Ghaffari " Test Case : %s\n" 179139613f2SLeila Ghaffari " Background Velocity : %f,%f,%f\n" 180139613f2SLeila Ghaffari " Stabilization : %s\n", 1812b916ea7SJeremy L Thompson app_ctx->problem_name, EulerTestTypes[euler_ctx->euler_test], euler_ctx->mean_velocity[0], euler_ctx->mean_velocity[1], 1822b916ea7SJeremy L Thompson euler_ctx->mean_velocity[2], StabilizationTypes[euler_ctx->stabilization])); 183a515125bSLeila Ghaffari 184e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfctx, &euler_ctx)); 185d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 186a515125bSLeila Ghaffari } 187