1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3 4 /// @file 5 /// Utility functions for setting up EULER_VORTEX 6 7 #include "../qfunctions/eulervortex.h" 8 9 #include <ceed.h> 10 #include <petscdm.h> 11 12 #include <navierstokes.h> 13 14 static PetscErrorCode EulerVortexOutflowBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) { 15 HoneeBCStruct honee_bc; 16 17 PetscFunctionBeginUser; 18 PetscCall(BCDefinitionGetContext(bc_def, &honee_bc)); 19 PetscCall(HoneeBCCreateIFunctionQF(bc_def, Euler_Outflow, Euler_Outflow_loc, honee_bc->qfctx, qf)); 20 PetscFunctionReturn(PETSC_SUCCESS); 21 } 22 23 static PetscErrorCode EulerVortexInflowBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) { 24 HoneeBCStruct honee_bc; 25 26 PetscFunctionBeginUser; 27 PetscCall(BCDefinitionGetContext(bc_def, &honee_bc)); 28 PetscCall(HoneeBCCreateIFunctionQF(bc_def, TravelingVortex_Inflow, TravelingVortex_Inflow_loc, honee_bc->qfctx, qf)); 29 PetscFunctionReturn(PETSC_SUCCESS); 30 } 31 32 PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx) { 33 EulerTestType euler_test; 34 Honee honee = *(Honee *)ctx; 35 StabilizationType stab; 36 MPI_Comm comm = honee->comm; 37 Ceed ceed = honee->ceed; 38 PetscBool implicit; 39 EulerContext euler_ctx; 40 CeedQFunctionContext euler_qfctx; 41 PetscInt dim; 42 43 PetscFunctionBeginUser; 44 PetscCall(PetscCalloc1(1, &euler_ctx)); 45 46 // ------------------------------------------------------ 47 // SET UP EULER VORTEX 48 // ------------------------------------------------------ 49 problem->ics.qf_func_ptr = ICsEuler; 50 problem->ics.qf_loc = ICsEuler_loc; 51 problem->apply_vol_rhs.qf_func_ptr = Euler; 52 problem->apply_vol_rhs.qf_loc = Euler_loc; 53 problem->apply_vol_ifunction.qf_func_ptr = IFunction_Euler; 54 problem->apply_vol_ifunction.qf_loc = IFunction_Euler_loc; 55 problem->num_comps_jac_data = 0; 56 problem->compute_exact_solution_error = PETSC_TRUE; 57 problem->print_info = PRINT_EULER_VORTEX; 58 59 // ------------------------------------------------------ 60 // Create the libCEED context 61 // ------------------------------------------------------ 62 CeedScalar vortex_strength = 5.; // - 63 CeedScalar c_tau = 0.5; // - 64 // c_tau = 0.5 is reported as "optimal" in Hughes et al 2010 65 PetscReal center[3] = {0.}, // m 66 mean_velocity[3] = {1., 1., 0}; // m/s 67 PetscReal domain_min[3], domain_max[3], domain_size[3] = {0.}; 68 PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 69 PetscCall(DMGetDimension(dm, &dim)); 70 for (PetscInt i = 0; i < dim; i++) domain_size[i] = domain_max[i] - domain_min[i]; 71 72 // ------------------------------------------------------ 73 // Create the PETSc context 74 // ------------------------------------------------------ 75 PetscScalar meter = 1e-2; // 1 meter in scaled length units 76 PetscScalar second = 1e-2; // 1 second in scaled time units 77 78 // ------------------------------------------------------ 79 // Command line Options 80 // ------------------------------------------------------ 81 PetscOptionsBegin(comm, NULL, "Options for EULER_VORTEX problem", NULL); 82 // -- Physics 83 PetscCall(PetscOptionsScalar("-vortex_strength", "Strength of Vortex", NULL, vortex_strength, &vortex_strength, NULL)); 84 PetscInt n = dim; 85 PetscBool user_velocity; 86 PetscCall(PetscOptionsRealArray("-mean_velocity", "Background velocity vector", NULL, mean_velocity, &n, &user_velocity)); 87 for (PetscInt i = 0; i < dim; i++) center[i] = .5 * domain_size[i]; 88 n = dim; 89 PetscCall(PetscOptionsRealArray("-center", "Location of vortex center", NULL, center, &n, NULL)); 90 PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 91 PetscCall(PetscOptionsEnum("-euler_test", "Euler test option", NULL, EulerTestTypes, (PetscEnum)(euler_test = EULER_TEST_ISENTROPIC_VORTEX), 92 (PetscEnum *)&euler_test, NULL)); 93 PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 94 PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 95 // -- Units 96 PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 97 meter = fabs(meter); 98 PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 99 second = fabs(second); 100 101 // -- Warnings 102 if (stab == STAB_SUPG && !implicit) { 103 PetscCall(PetscPrintf(comm, "Warning! Use -stab supg only with -implicit\n")); 104 } 105 if (user_velocity && (euler_test == EULER_TEST_1 || euler_test == EULER_TEST_3)) { 106 PetscCall(PetscPrintf(comm, "Warning! Background velocity vector for -euler_test t1 and -euler_test t3 is (0,0,0)\n")); 107 } 108 109 PetscOptionsEnd(); 110 111 // ------------------------------------------------------ 112 // Set up the PETSc context 113 // ------------------------------------------------------ 114 honee->units->meter = meter; 115 honee->units->second = second; 116 117 // ------------------------------------------------------ 118 // Set up the libCEED context 119 // ------------------------------------------------------ 120 // -- Scale variables to desired units 121 for (PetscInt i = 0; i < 3; i++) { 122 center[i] *= meter; 123 domain_size[i] *= meter; 124 mean_velocity[i] *= (meter / second); 125 } 126 127 // -- QFunction Context 128 honee->phys->implicit = implicit; 129 euler_ctx->curr_time = 0.; 130 euler_ctx->implicit = implicit; 131 euler_ctx->euler_test = euler_test; 132 euler_ctx->center[0] = center[0]; 133 euler_ctx->center[1] = center[1]; 134 euler_ctx->center[2] = center[2]; 135 euler_ctx->vortex_strength = vortex_strength; 136 euler_ctx->c_tau = c_tau; 137 euler_ctx->mean_velocity[0] = mean_velocity[0]; 138 euler_ctx->mean_velocity[1] = mean_velocity[1]; 139 euler_ctx->mean_velocity[2] = mean_velocity[2]; 140 euler_ctx->stabilization = stab; 141 142 PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &euler_qfctx)); 143 PetscCallCeed(ceed, CeedQFunctionContextSetData(euler_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*euler_ctx), euler_ctx)); 144 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(euler_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 145 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(euler_qfctx, "solution time", offsetof(struct EulerContext_, curr_time), 1, 146 "Physical time of the solution")); 147 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->ics.qfctx)); 148 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->apply_vol_rhs.qfctx)); 149 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &problem->apply_vol_ifunction.qfctx)); 150 151 for (PetscCount b = 0; b < problem->num_bc_defs; b++) { 152 BCDefinition bc_def = problem->bc_defs[b]; 153 const char *name; 154 155 PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL)); 156 if (!strcmp(name, "outflow")) { 157 HoneeBCStruct honee_bc; 158 159 PetscCall(PetscNew(&honee_bc)); 160 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &honee_bc->qfctx)); 161 honee_bc->honee = honee; 162 honee_bc->num_comps_jac_data = 0; 163 PetscCall(BCDefinitionSetContext(bc_def, HoneeBCDestroy, honee_bc)); 164 165 PetscCall(BCDefinitionSetIFunction(bc_def, EulerVortexOutflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp)); 166 PetscCall(BCDefinitionSetIJacobian(bc_def, NULL, NULL)); 167 } else if (!strcmp(name, "inflow")) { 168 HoneeBCStruct honee_bc; 169 170 PetscCall(PetscNew(&honee_bc)); 171 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(euler_qfctx, &honee_bc->qfctx)); 172 honee_bc->honee = honee; 173 honee_bc->num_comps_jac_data = 0; 174 PetscCall(BCDefinitionSetContext(bc_def, HoneeBCDestroy, honee_bc)); 175 176 PetscCall(BCDefinitionSetIFunction(bc_def, EulerVortexInflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp)); 177 PetscCall(BCDefinitionSetIJacobian(bc_def, NULL, NULL)); 178 } 179 } 180 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&euler_qfctx)); 181 PetscFunctionReturn(PETSC_SUCCESS); 182 } 183 184 PetscErrorCode PRINT_EULER_VORTEX(Honee honee, ProblemData problem, AppCtx app_ctx) { 185 MPI_Comm comm = honee->comm; 186 Ceed ceed = honee->ceed; 187 EulerContext euler_ctx; 188 189 PetscFunctionBeginUser; 190 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfctx, CEED_MEM_HOST, &euler_ctx)); 191 PetscCall(PetscPrintf(comm, 192 " Problem:\n" 193 " Problem Name : %s\n" 194 " Test Case : %s\n" 195 " Background Velocity : %f,%f,%f\n" 196 " Stabilization : %s\n", 197 app_ctx->problem_name, EulerTestTypes[euler_ctx->euler_test], euler_ctx->mean_velocity[0], euler_ctx->mean_velocity[1], 198 euler_ctx->mean_velocity[2], StabilizationTypes[euler_ctx->stabilization])); 199 200 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfctx, &euler_ctx)); 201 PetscFunctionReturn(PETSC_SUCCESS); 202 } 203