1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 /// @file 9 /// Utility functions for setting up problems using the Newtonian Qfunction 10 11 #include "../qfunctions/newtonian.h" 12 13 #include "../navierstokes.h" 14 #include "../qfunctions/setupgeo.h" 15 16 // For use with PetscOptionsEnum 17 static const char *const StateVariables[] = {"CONSERVATIVE", "PRIMITIVE", "StateVariable", "STATEVAR_", NULL}; 18 19 // Compute relative error |a - b|/|s| 20 static PetscErrorCode CheckPrimitiveWithTolerance(StatePrimitive sY, StatePrimitive aY, StatePrimitive bY, const char *name, PetscReal rtol_pressure, 21 PetscReal rtol_velocity, PetscReal rtol_temperature) { 22 PetscFunctionBeginUser; 23 StatePrimitive eY; // relative error 24 eY.pressure = (aY.pressure - bY.pressure) / sY.pressure; 25 PetscScalar u = sqrt(Square(sY.velocity[0]) + Square(sY.velocity[1]) + Square(sY.velocity[2])); 26 for (int j = 0; j < 3; j++) eY.velocity[j] = (aY.velocity[j] - bY.velocity[j]) / u; 27 eY.temperature = (aY.temperature - bY.temperature) / sY.temperature; 28 if (fabs(eY.pressure) > rtol_pressure) printf("%s: pressure error %g\n", name, eY.pressure); 29 for (int j = 0; j < 3; j++) { 30 if (fabs(eY.velocity[j]) > rtol_velocity) printf("%s: velocity[%d] error %g\n", name, j, eY.velocity[j]); 31 } 32 if (fabs(eY.temperature) > rtol_temperature) printf("%s: temperature error %g\n", name, eY.temperature); 33 PetscFunctionReturn(0); 34 } 35 36 static PetscErrorCode UnitTests_Newtonian(User user, NewtonianIdealGasContext gas) { 37 Units units = user->units; 38 const CeedScalar eps = 1e-6; 39 const CeedScalar kg = units->kilogram, m = units->meter, sec = units->second, Pascal = units->Pascal; 40 PetscFunctionBeginUser; 41 const CeedScalar rho = 1.2 * kg / (m * m * m), u = 40 * m / sec; 42 CeedScalar U[5] = {rho, rho * u, rho * u * 1.1, rho * u * 1.2, 250e3 * Pascal + .5 * rho * u * u}; 43 const CeedScalar x[3] = {.1, .2, .3}; 44 State s = StateFromU(gas, U, x); 45 for (int i = 0; i < 8; i++) { 46 CeedScalar dU[5] = {0}, dx[3] = {0}; 47 if (i < 5) dU[i] = U[i]; 48 else dx[i - 5] = x[i - 5]; 49 State ds = StateFromU_fwd(gas, s, dU, x, dx); 50 for (int j = 0; j < 5; j++) dU[j] = (1 + eps * (i == j)) * U[j]; 51 for (int j = 0; j < 3; j++) dx[j] = (1 + eps * (i == 5 + j)) * x[j]; 52 State t = StateFromU(gas, dU, dx); 53 StatePrimitive dY; 54 dY.pressure = (t.Y.pressure - s.Y.pressure) / eps; 55 for (int j = 0; j < 3; j++) dY.velocity[j] = (t.Y.velocity[j] - s.Y.velocity[j]) / eps; 56 dY.temperature = (t.Y.temperature - s.Y.temperature) / eps; 57 char buf[128]; 58 snprintf(buf, sizeof buf, "StateFromU_fwd i=%d", i); 59 PetscCall(CheckPrimitiveWithTolerance(dY, ds.Y, dY, buf, 5e-6, 1e-6, 1e-6)); 60 } 61 PetscFunctionReturn(0); 62 } 63 64 PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) { 65 SetupContext setup_context; 66 User user = *(User *)ctx; 67 CeedInt degree = user->app_ctx->degree; 68 StabilizationType stab; 69 StateVariable state_var; 70 MPI_Comm comm = PETSC_COMM_WORLD; 71 PetscBool implicit; 72 PetscBool has_curr_time = PETSC_FALSE, unit_tests; 73 NewtonianIdealGasContext newtonian_ig_ctx; 74 CeedQFunctionContext newtonian_ig_context; 75 76 PetscFunctionBeginUser; 77 PetscCall(PetscCalloc1(1, &setup_context)); 78 PetscCall(PetscCalloc1(1, &newtonian_ig_ctx)); 79 80 // ------------------------------------------------------ 81 // Setup Generic Newtonian IG Problem 82 // ------------------------------------------------------ 83 problem->dim = 3; 84 problem->q_data_size_vol = 10; 85 problem->q_data_size_sur = 10; 86 problem->jac_data_size_sur = 11; 87 problem->setup_vol.qfunction = Setup; 88 problem->setup_vol.qfunction_loc = Setup_loc; 89 problem->setup_sur.qfunction = SetupBoundary; 90 problem->setup_sur.qfunction_loc = SetupBoundary_loc; 91 problem->bc = NULL; 92 problem->bc_ctx = setup_context; 93 problem->non_zero_time = PETSC_FALSE; 94 problem->print_info = PRINT_NEWTONIAN; 95 96 // ------------------------------------------------------ 97 // Create the libCEED context 98 // ------------------------------------------------------ 99 CeedScalar cv = 717.; // J/(kg K) 100 CeedScalar cp = 1004.; // J/(kg K) 101 CeedScalar g[3] = {0, 0, -9.81}; // m/s^2 102 CeedScalar lambda = -2. / 3.; // - 103 CeedScalar mu = 1.8e-5; // Pa s, dynamic viscosity 104 CeedScalar k = 0.02638; // W/(m K) 105 CeedScalar c_tau = 0.5 / degree; // - 106 CeedScalar Ctau_t = 1.0; // - 107 CeedScalar Ctau_v = 36.0; // TODO make function of degree 108 CeedScalar Ctau_C = 1.0; // TODO make function of degree 109 CeedScalar Ctau_M = 1.0; // TODO make function of degree 110 CeedScalar Ctau_E = 1.0; // TODO make function of degree 111 PetscReal domain_min[3], domain_max[3], domain_size[3]; 112 PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 113 for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 114 115 // ------------------------------------------------------ 116 // Create the PETSc context 117 // ------------------------------------------------------ 118 PetscScalar meter = 1; // 1 meter in scaled length units 119 PetscScalar kilogram = 1; // 1 kilogram in scaled mass units 120 PetscScalar second = 1; // 1 second in scaled time units 121 PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 122 PetscScalar W_per_m_K, Pascal, J_per_kg_K, m_per_squared_s; 123 124 // ------------------------------------------------------ 125 // Command line Options 126 // ------------------------------------------------------ 127 PetscOptionsBegin(comm, NULL, "Options for Newtonian Ideal Gas based problem", NULL); 128 // -- Conservative vs Primitive variables 129 PetscCall(PetscOptionsEnum("-state_var", "State variables used", NULL, StateVariables, (PetscEnum)(state_var = STATEVAR_CONSERVATIVE), 130 (PetscEnum *)&state_var, NULL)); 131 132 switch (state_var) { 133 case STATEVAR_CONSERVATIVE: 134 problem->ics.qfunction = ICsNewtonianIG; 135 problem->ics.qfunction_loc = ICsNewtonianIG_loc; 136 problem->apply_vol_rhs.qfunction = RHSFunction_Newtonian; 137 problem->apply_vol_rhs.qfunction_loc = RHSFunction_Newtonian_loc; 138 problem->apply_vol_ifunction.qfunction = IFunction_Newtonian_Conserv; 139 problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_Conserv_loc; 140 problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian_Conserv; 141 problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_Conserv_loc; 142 problem->apply_inflow.qfunction = BoundaryIntegral_Conserv; 143 problem->apply_inflow.qfunction_loc = BoundaryIntegral_Conserv_loc; 144 problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Conserv; 145 problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Conserv_loc; 146 problem->apply_outflow.qfunction = PressureOutflow_Conserv; 147 problem->apply_outflow.qfunction_loc = PressureOutflow_Conserv_loc; 148 problem->apply_outflow_jacobian.qfunction = PressureOutflow_Jacobian_Conserv; 149 problem->apply_outflow_jacobian.qfunction_loc = PressureOutflow_Jacobian_Conserv_loc; 150 break; 151 152 case STATEVAR_PRIMITIVE: 153 problem->ics.qfunction = ICsNewtonianIG_Prim; 154 problem->ics.qfunction_loc = ICsNewtonianIG_Prim_loc; 155 problem->apply_vol_ifunction.qfunction = IFunction_Newtonian_Prim; 156 problem->apply_vol_ifunction.qfunction_loc = IFunction_Newtonian_Prim_loc; 157 problem->apply_vol_ijacobian.qfunction = IJacobian_Newtonian_Prim; 158 problem->apply_vol_ijacobian.qfunction_loc = IJacobian_Newtonian_Prim_loc; 159 problem->apply_inflow.qfunction = BoundaryIntegral_Prim; 160 problem->apply_inflow.qfunction_loc = BoundaryIntegral_Prim_loc; 161 problem->apply_inflow_jacobian.qfunction = BoundaryIntegral_Jacobian_Prim; 162 problem->apply_inflow_jacobian.qfunction_loc = BoundaryIntegral_Jacobian_Prim_loc; 163 problem->apply_outflow.qfunction = PressureOutflow_Prim; 164 problem->apply_outflow.qfunction_loc = PressureOutflow_Prim_loc; 165 problem->apply_outflow_jacobian.qfunction = PressureOutflow_Jacobian_Prim; 166 problem->apply_outflow_jacobian.qfunction_loc = PressureOutflow_Jacobian_Prim_loc; 167 break; 168 } 169 170 // -- Physics 171 PetscCall(PetscOptionsScalar("-cv", "Heat capacity at constant volume", NULL, cv, &cv, NULL)); 172 PetscCall(PetscOptionsScalar("-cp", "Heat capacity at constant pressure", NULL, cp, &cp, NULL)); 173 PetscCall(PetscOptionsScalar("-lambda", "Stokes hypothesis second viscosity coefficient", NULL, lambda, &lambda, NULL)); 174 PetscCall(PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", NULL, mu, &mu, NULL)); 175 PetscCall(PetscOptionsScalar("-k", "Thermal conductivity", NULL, k, &k, NULL)); 176 177 PetscInt dim = problem->dim; 178 PetscCall(PetscOptionsRealArray("-g", "Gravitational acceleration", NULL, g, &dim, NULL)); 179 PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 180 PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 181 PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 182 PetscCall(PetscOptionsScalar("-Ctau_v", "Stabilization viscous constant", NULL, Ctau_v, &Ctau_v, NULL)); 183 PetscCall(PetscOptionsScalar("-Ctau_C", "Stabilization continuity constant", NULL, Ctau_C, &Ctau_C, NULL)); 184 PetscCall(PetscOptionsScalar("-Ctau_M", "Stabilization momentum constant", NULL, Ctau_M, &Ctau_M, NULL)); 185 PetscCall(PetscOptionsScalar("-Ctau_E", "Stabilization energy constant", NULL, Ctau_E, &Ctau_E, NULL)); 186 PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 187 PetscCall(PetscOptionsBool("-newtonian_unit_tests", "Run Newtonian unit tests", NULL, unit_tests = PETSC_FALSE, &unit_tests, NULL)); 188 189 // -- Units 190 PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 191 meter = fabs(meter); 192 PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL)); 193 kilogram = fabs(kilogram); 194 PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 195 second = fabs(second); 196 PetscCall(PetscOptionsScalar("-units_Kelvin", "1 Kelvin in scaled temperature units", NULL, Kelvin, &Kelvin, NULL)); 197 Kelvin = fabs(Kelvin); 198 199 // -- Warnings 200 if (stab == STAB_SUPG && !implicit) { 201 PetscCall(PetscPrintf(comm, "Warning! Use -stab supg only with -implicit\n")); 202 } 203 if (state_var == STATEVAR_PRIMITIVE && !implicit) { 204 SETERRQ(comm, PETSC_ERR_ARG_NULL, "RHSFunction is not provided for primitive variables (use -state_var primitive only with -implicit)\n"); 205 } 206 PetscOptionsEnd(); 207 208 // ------------------------------------------------------ 209 // Set up the PETSc context 210 // ------------------------------------------------------ 211 // -- Define derived units 212 Pascal = kilogram / (meter * PetscSqr(second)); 213 J_per_kg_K = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 214 m_per_squared_s = meter / PetscSqr(second); 215 W_per_m_K = kilogram * meter / (pow(second, 3) * Kelvin); 216 217 user->units->meter = meter; 218 user->units->kilogram = kilogram; 219 user->units->second = second; 220 user->units->Kelvin = Kelvin; 221 user->units->Pascal = Pascal; 222 user->units->J_per_kg_K = J_per_kg_K; 223 user->units->m_per_squared_s = m_per_squared_s; 224 user->units->W_per_m_K = W_per_m_K; 225 226 // ------------------------------------------------------ 227 // Set up the libCEED context 228 // ------------------------------------------------------ 229 // -- Scale variables to desired units 230 cv *= J_per_kg_K; 231 cp *= J_per_kg_K; 232 mu *= Pascal * second; 233 k *= W_per_m_K; 234 for (PetscInt i = 0; i < 3; i++) domain_size[i] *= meter; 235 for (PetscInt i = 0; i < 3; i++) g[i] *= m_per_squared_s; 236 problem->dm_scale = meter; 237 238 // -- Setup Context 239 setup_context->cv = cv; 240 setup_context->cp = cp; 241 setup_context->lx = domain_size[0]; 242 setup_context->ly = domain_size[1]; 243 setup_context->lz = domain_size[2]; 244 setup_context->time = 0; 245 PetscCall(PetscArraycpy(setup_context->g, g, 3)); 246 247 // -- Solver Settings 248 user->phys->stab = stab; 249 user->phys->implicit = implicit; 250 user->phys->state_var = state_var; 251 user->phys->has_curr_time = has_curr_time; 252 253 // -- QFunction Context 254 newtonian_ig_ctx->lambda = lambda; 255 newtonian_ig_ctx->mu = mu; 256 newtonian_ig_ctx->k = k; 257 newtonian_ig_ctx->cv = cv; 258 newtonian_ig_ctx->cp = cp; 259 newtonian_ig_ctx->c_tau = c_tau; 260 newtonian_ig_ctx->Ctau_t = Ctau_t; 261 newtonian_ig_ctx->Ctau_v = Ctau_v; 262 newtonian_ig_ctx->Ctau_C = Ctau_C; 263 newtonian_ig_ctx->Ctau_M = Ctau_M; 264 newtonian_ig_ctx->Ctau_E = Ctau_E; 265 newtonian_ig_ctx->stabilization = stab; 266 newtonian_ig_ctx->is_implicit = implicit; 267 newtonian_ig_ctx->state_var = state_var; 268 PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3)); 269 270 if (bc->num_freestream > 0) PetscCall(FreestreamBCSetup(problem, dm, ctx, newtonian_ig_ctx)); 271 272 CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context); 273 CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context); 274 CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc); 275 CeedQFunctionContextRegisterDouble(problem->ics.qfunction_context, "evaluation time", (char *)&setup_context->time - (char *)setup_context, 1, 276 "Time of evaluation"); 277 278 CeedQFunctionContextCreate(user->ceed, &newtonian_ig_context); 279 CeedQFunctionContextSetData(newtonian_ig_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*newtonian_ig_ctx), newtonian_ig_ctx); 280 CeedQFunctionContextSetDataDestroy(newtonian_ig_context, CEED_MEM_HOST, FreeContextPetsc); 281 CeedQFunctionContextRegisterDouble(newtonian_ig_context, "timestep size", offsetof(struct NewtonianIdealGasContext_, dt), 1, 282 "Size of timestep, delta t"); 283 CeedQFunctionContextRegisterDouble(newtonian_ig_context, "ijacobian time shift", offsetof(struct NewtonianIdealGasContext_, ijacobian_time_shift), 284 1, "Shift for mass matrix in IJacobian"); 285 problem->apply_vol_rhs.qfunction_context = newtonian_ig_context; 286 CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ifunction.qfunction_context); 287 CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_vol_ijacobian.qfunction_context); 288 CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow.qfunction_context); 289 CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_inflow_jacobian.qfunction_context); 290 CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_outflow.qfunction_context); 291 CeedQFunctionContextReferenceCopy(newtonian_ig_context, &problem->apply_outflow_jacobian.qfunction_context); 292 293 if (unit_tests) { 294 PetscCall(UnitTests_Newtonian(user, newtonian_ig_ctx)); 295 } 296 PetscFunctionReturn(0); 297 } 298 299 PetscErrorCode PRINT_NEWTONIAN(ProblemData *problem, AppCtx app_ctx) { 300 MPI_Comm comm = PETSC_COMM_WORLD; 301 NewtonianIdealGasContext newtonian_ctx; 302 303 PetscFunctionBeginUser; 304 CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ctx); 305 PetscCall(PetscPrintf(comm, 306 " Problem:\n" 307 " Problem Name : %s\n" 308 " Stabilization : %s\n", 309 app_ctx->problem_name, StabilizationTypes[newtonian_ctx->stabilization])); 310 CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ctx); 311 PetscFunctionReturn(0); 312 } 313