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