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