1 static char help[] = "Time-dependent Low Mach Flow in 2d and 3d channels with finite elements.\n\ 2 We solve the Low Mach flow problem for both conducting and non-conducting fluids,\n\ 3 using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n"; 4 5 /*F 6 The non-conducting Low Mach flow is time-dependent isoviscous Navier-Stokes flow. We discretize using the 7 finite element method on an unstructured mesh. The weak form equations are 8 9 \begin{align*} 10 < q, \nabla\cdot u > = 0 11 <v, du/dt> + <v, u \cdot \nabla u> + < \nabla v, \nu (\nabla u + {\nabla u}^T) > - < \nabla\cdot v, p > - < v, f > = 0 12 < w, u \cdot \nabla T > + < \nabla w, \alpha \nabla T > - < w, Q > = 0 13 \end{align*} 14 15 where $\nu$ is the kinematic viscosity and $\alpha$ is thermal diffusivity. 16 17 The conducting form is given in the ABLATE documentation [1,2] and derived in Principe and Codina [2]. 18 19 For visualization, use 20 21 -dm_view hdf5:$PWD/sol.h5 -sol_vec_view hdf5:$PWD/sol.h5::append -exact_vec_view hdf5:$PWD/sol.h5::append 22 23 [1] https://ubchrest.github.io/ablate/content/formulations/lowMachFlow/ 24 [2] https://github.com/UBCHREST/ablate/blob/main/ablateCore/flow/lowMachFlow.c 25 [3] J. Principe and R. Codina, "Mathematical models for thermally coupled low speed flows", Adv. in Theo. and App. Mech., 2(1), pp.93--112, 2009. 26 F*/ 27 28 #include <petscdmplex.h> 29 #include <petscsnes.h> 30 #include <petscts.h> 31 #include <petscds.h> 32 #include <petscbag.h> 33 34 typedef enum {MOD_INCOMPRESSIBLE, MOD_CONDUCTING, NUM_MOD_TYPES} ModType; 35 const char *modTypes[NUM_MOD_TYPES+1] = {"incompressible", "conducting", "unknown"}; 36 37 typedef enum {SOL_QUADRATIC, SOL_CUBIC, SOL_CUBIC_TRIG, SOL_TAYLOR_GREEN, SOL_PIPE, SOL_PIPE_WIGGLY, NUM_SOL_TYPES} SolType; 38 const char *solTypes[NUM_SOL_TYPES+1] = {"quadratic", "cubic", "cubic_trig", "taylor_green", "pipe", "pipe_wiggly", "unknown"}; 39 40 /* Fields */ 41 const PetscInt VEL = 0; 42 const PetscInt PRES = 1; 43 const PetscInt TEMP = 2; 44 /* Sources */ 45 const PetscInt MOMENTUM = 0; 46 const PetscInt MASS = 1; 47 const PetscInt ENERGY = 2; 48 /* Constants */ 49 const PetscInt STROUHAL = 0; 50 const PetscInt FROUDE = 1; 51 const PetscInt REYNOLDS = 2; 52 const PetscInt PECLET = 3; 53 const PetscInt P_TH = 4; 54 const PetscInt MU = 5; 55 const PetscInt NU = 6; 56 const PetscInt C_P = 7; 57 const PetscInt K = 8; 58 const PetscInt ALPHA = 9; 59 const PetscInt T_IN = 10; 60 const PetscInt G_DIR = 11; 61 const PetscInt EPSILON = 12; 62 63 typedef struct { 64 PetscReal Strouhal; /* Strouhal number */ 65 PetscReal Froude; /* Froude number */ 66 PetscReal Reynolds; /* Reynolds number */ 67 PetscReal Peclet; /* Peclet number */ 68 PetscReal p_th; /* Thermodynamic pressure */ 69 PetscReal mu; /* Dynamic viscosity */ 70 PetscReal nu; /* Kinematic viscosity */ 71 PetscReal c_p; /* Specific heat at constant pressure */ 72 PetscReal k; /* Thermal conductivity */ 73 PetscReal alpha; /* Thermal diffusivity */ 74 PetscReal T_in; /* Inlet temperature */ 75 PetscReal g_dir; /* Gravity direction */ 76 PetscReal epsilon; /* Strength of perturbation */ 77 } Parameter; 78 79 typedef struct { 80 /* Problem definition */ 81 PetscBag bag; /* Holds problem parameters */ 82 ModType modType; /* Model type */ 83 SolType solType; /* MMS solution type */ 84 PetscBool hasNullSpace; /* Problem has the constant null space for pressure */ 85 /* Flow diagnostics */ 86 DM dmCell; /* A DM with piecewise constant discretization */ 87 } AppCtx; 88 89 static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 90 { 91 PetscInt d; 92 for (d = 0; d < Nc; ++d) u[d] = 0.0; 93 return 0; 94 } 95 96 static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 97 { 98 PetscInt d; 99 for (d = 0; d < Nc; ++d) u[d] = 1.0; 100 return 0; 101 } 102 103 /* 104 CASE: quadratic 105 In 2D we use exact solution: 106 107 u = t + x^2 + y^2 108 v = t + 2x^2 - 2xy 109 p = x + y - 1 110 T = t + x + y + 1 111 f = <t (2x + 2y) + 2x^3 + 4x^2y - 2xy^2 -4\nu + 2, t (2x - 2y) + 4xy^2 + 2x^2y - 2y^3 -4\nu + 2> 112 Q = 1 + 2t + 3x^2 - 2xy + y^2 113 114 so that 115 116 \nabla \cdot u = 2x - 2x = 0 117 118 f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p 119 = <1, 1> + <t + x^2 + y^2, t + 2x^2 - 2xy> . <<2x, 4x - 2y>, <2y, -2x>> - \nu <4, 4> + <1, 1> 120 = <t (2x + 2y) + 2x^3 + 4x^2y - 2xy^2, t (2x - 2y) + 2x^2y + 4xy^2 - 2y^3> + <-4 \nu + 2, -4\nu + 2> 121 = <t (2x + 2y) + 2x^3 + 4x^2y - 2xy^2 - 4\nu + 2, t (2x - 2y) + 4xy^2 + 2x^2y - 2y^3 - 4\nu + 2> 122 123 Q = dT/dt + u \cdot \nabla T - \alpha \Delta T 124 = 1 + <t + x^2 + y^2, t + 2x^2 - 2xy> . <1, 1> - \alpha 0 125 = 1 + 2t + 3x^2 - 2xy + y^2 126 */ 127 128 static PetscErrorCode quadratic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 129 { 130 u[0] = time + X[0]*X[0] + X[1]*X[1]; 131 u[1] = time + 2.0*X[0]*X[0] - 2.0*X[0]*X[1]; 132 return 0; 133 } 134 static PetscErrorCode quadratic_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 135 { 136 u[0] = 1.0; 137 u[1] = 1.0; 138 return 0; 139 } 140 141 static PetscErrorCode quadratic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 142 { 143 p[0] = X[0] + X[1] - 1.0; 144 return 0; 145 } 146 147 static PetscErrorCode quadratic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 148 { 149 T[0] = time + X[0] + X[1] + 1.0; 150 return 0; 151 } 152 static PetscErrorCode quadratic_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 153 { 154 T[0] = 1.0; 155 return 0; 156 } 157 158 static void f0_quadratic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 159 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 160 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 161 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 162 { 163 const PetscReal nu = PetscRealPart(constants[NU]); 164 165 f0[0] -= t*(2*X[0] + 2*X[1]) + 2*X[0]*X[0]*X[0] + 4*X[0]*X[0]*X[1] - 2*X[0]*X[1]*X[1] - 4.0*nu + 2; 166 f0[1] -= t*(2*X[0] - 2*X[1]) + 4*X[0]*X[1]*X[1] + 2*X[0]*X[0]*X[1] - 2*X[1]*X[1]*X[1] - 4.0*nu + 2; 167 } 168 169 static void f0_quadratic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 170 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 171 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 172 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 173 { 174 f0[0] -= 2*t + 1 + 3*X[0]*X[0] - 2*X[0]*X[1] + X[1]*X[1]; 175 } 176 177 /* 178 CASE: quadratic 179 In 2D we use exact solution: 180 181 u = t + x^2 + y^2 182 v = t + 2x^2 - 2xy 183 p = x + y - 1 184 T = t + x + y + 1 185 rho = p^{th} / T 186 187 so that 188 189 \nabla \cdot u = 2x - 2x = 0 190 grad u = <<2 x, 4x - 2y>, <2 y, -2x>> 191 epsilon(u) = 1/2 (grad u + grad u^T) = <<2x, 2x>, <2x, -2x>> 192 epsilon'(u) = epsilon(u) - 1/3 (div u) I = epsilon(u) 193 div epsilon'(u) = <2, 2> 194 195 f = rho S du/dt + rho u \cdot \nabla u - 2\mu/Re div \epsilon'(u) + \nabla p + rho / F^2 \hat y 196 = rho S <1, 1> + rho <t + x^2 + y^2, t + 2x^2 - 2xy> . <<2x, 4x - 2y>, <2y, -2x>> - 2\mu/Re <2, 2> + <1, 1> + rho/F^2 <0, 1> 197 = rho S <1, 1> + rho <t (2x + 2y) + 2x^3 + 4x^2y - 2xy^2, t (2x - 2y) + 2x^2y + 4xy^2 - 2y^3> - mu/Re <4, 4> + <1, 1> + rho/F^2 <0, 1> 198 199 g = S rho_t + div (rho u) 200 = -S pth T_t/T^2 + rho div (u) + u . grad rho 201 = -S pth 1/T^2 - pth u . grad T / T^2 202 = -pth / T^2 (S + 2t + 3 x^2 - 2xy + y^2) 203 204 Q = rho c_p S dT/dt + rho c_p u . grad T - 1/Pe div k grad T 205 = c_p S pth / T + c_p pth (2t + 3 x^2 - 2xy + y^2) / T - k/Pe 0 206 = c_p pth / T (S + 2t + 3 x^2 - 2xy + y^2) 207 */ 208 static void f0_conduct_quadratic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 209 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 210 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 211 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 212 { 213 const PetscReal S = PetscRealPart(constants[STROUHAL]); 214 const PetscReal F = PetscRealPart(constants[FROUDE]); 215 const PetscReal Re = PetscRealPart(constants[REYNOLDS]); 216 const PetscReal mu = PetscRealPart(constants[MU]); 217 const PetscReal p_th = PetscRealPart(constants[P_TH]); 218 const PetscReal rho = p_th / (t + X[0] + X[1] + 1.); 219 const PetscInt gd = (PetscInt) PetscRealPart(constants[G_DIR]); 220 221 f0[0] -= rho * S + rho * (2.*t*(X[0] + X[1]) + 2.*X[0]*X[0]*X[0] + 4.*X[0]*X[0]*X[1] - 2.*X[0]*X[1]*X[1]) - 4.*mu/Re + 1.; 222 f0[1] -= rho * S + rho * (2.*t*(X[0] - X[1]) + 2.*X[0]*X[0]*X[1] + 4.*X[0]*X[1]*X[1] - 2.*X[1]*X[1]*X[1]) - 4.*mu/Re + 1.; 223 f0[gd] -= rho/PetscSqr(F); 224 } 225 226 static void f0_conduct_quadratic_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 227 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 228 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 229 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 230 { 231 const PetscReal S = PetscRealPart(constants[STROUHAL]); 232 const PetscReal p_th = PetscRealPart(constants[P_TH]); 233 234 f0[0] += p_th * (S + 2.*t + 3.*X[0]*X[0] - 2.*X[0]*X[1] + X[1]*X[1]) / PetscSqr(t + X[0] + X[1] + 1.); 235 } 236 237 static void f0_conduct_quadratic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 238 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 239 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 240 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 241 { 242 const PetscReal S = PetscRealPart(constants[STROUHAL]); 243 const PetscReal c_p = PetscRealPart(constants[C_P]); 244 const PetscReal p_th = PetscRealPart(constants[P_TH]); 245 246 f0[0] -= c_p*p_th * (S + 2.*t + 3.*X[0]*X[0] - 2.*X[0]*X[1] + X[1]*X[1]) / (t + X[0] + X[1] + 1.); 247 } 248 249 /* 250 CASE: cubic 251 In 2D we use exact solution: 252 253 u = t + x^3 + y^3 254 v = t + 2x^3 - 3x^2y 255 p = 3/2 x^2 + 3/2 y^2 - 1 256 T = t + 1/2 x^2 + 1/2 y^2 257 f = < t(3x^2 + 3y^2) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y) + 3x + 1, 258 t(3x^2 - 6xy) + 6x^2y^3 + 3x^4y - 6xy^4 - \nu(12x - 6y) + 3y + 1> 259 Q = x^4 + xy^3 + 2x^3y - 3x^2y^2 + xt + yt - 2\alpha + 1 260 261 so that 262 263 \nabla \cdot u = 3x^2 - 3x^2 = 0 264 265 du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p - f 266 = <1,1> + <t(3x^2 + 3y^2) + 3x^5 + 6x^3y^2 - 6x^2y^3, t(3x^2 - 6xy) + 6x^2y^3 + 3x^4y - 6xy^4> - \nu<6x + 6y, 12x - 6y> + <3x, 3y> - <t(3x^2 + 3y^2) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y) + 3x + 1, t(3x^2 - 6xy) + 6x^2y^3 + 3x^4y - 6xy^4 - \nu(12x - 6y) + 3y + 1> = 0 267 268 dT/dt + u \cdot \nabla T - \alpha \Delta T - Q = 1 + (x^3 + y^3) x + (2x^3 - 3x^2y) y - 2*\alpha - (x^4 + xy^3 + 2x^3y - 3x^2y^2 - 2*\alpha +1) = 0 269 */ 270 static PetscErrorCode cubic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 271 { 272 u[0] = time + X[0]*X[0]*X[0] + X[1]*X[1]*X[1]; 273 u[1] = time + 2.0*X[0]*X[0]*X[0] - 3.0*X[0]*X[0]*X[1]; 274 return 0; 275 } 276 static PetscErrorCode cubic_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 277 { 278 u[0] = 1.0; 279 u[1] = 1.0; 280 return 0; 281 } 282 283 static PetscErrorCode cubic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 284 { 285 p[0] = 3.0*X[0]*X[0]/2.0 + 3.0*X[1]*X[1]/2.0 - 1.0; 286 return 0; 287 } 288 289 static PetscErrorCode cubic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 290 { 291 T[0] = time + X[0]*X[0]/2.0 + X[1]*X[1]/2.0; 292 return 0; 293 } 294 static PetscErrorCode cubic_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 295 { 296 T[0] = 1.0; 297 return 0; 298 } 299 300 static void f0_cubic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 301 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 302 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 303 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 304 { 305 const PetscReal nu = PetscRealPart(constants[NU]); 306 307 f0[0] -= (t*(3*X[0]*X[0] + 3*X[1]*X[1]) + 3*X[0]*X[0]*X[0]*X[0]*X[0] + 6*X[0]*X[0]*X[0]*X[1]*X[1] - 6*X[0]*X[0]*X[1]*X[1]*X[1] - ( 6*X[0] + 6*X[1])*nu + 3*X[0] + 1); 308 f0[1] -= (t*(3*X[0]*X[0] - 6*X[0]*X[1]) + 3*X[0]*X[0]*X[0]*X[0]*X[1] + 6*X[0]*X[0]*X[1]*X[1]*X[1] - 6*X[0]*X[1]*X[1]*X[1]*X[1] - (12*X[0] - 6*X[1])*nu + 3*X[1] + 1); 309 } 310 311 static void f0_cubic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 312 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 313 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 314 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 315 { 316 const PetscReal alpha = PetscRealPart(constants[ALPHA]); 317 318 f0[0] -= X[0]*X[0]*X[0]*X[0] + 2.0*X[0]*X[0]*X[0]*X[1] - 3.0*X[0]*X[0]*X[1]*X[1] + X[0]*X[1]*X[1]*X[1] + X[0]*t + X[1]*t - 2.0*alpha + 1; 319 } 320 321 /* 322 CASE: cubic-trigonometric 323 In 2D we use exact solution: 324 325 u = beta cos t + x^3 + y^3 326 v = beta sin t + 2x^3 - 3x^2y 327 p = 3/2 x^2 + 3/2 y^2 - 1 328 T = 20 cos t + 1/2 x^2 + 1/2 y^2 329 f = < beta cos t 3x^2 + beta sin t (3y^2 - 1) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y) + 3x, 330 beta cos t (6x^2 - 6xy) - beta sin t (3x^2) + 3x^4y + 6x^2y^3 - 6xy^4 - \nu(12x - 6y) + 3y> 331 Q = beta cos t x + beta sin t (y - 1) + x^4 + 2x^3y - 3x^2y^2 + xy^3 - 2\alpha 332 333 so that 334 335 \nabla \cdot u = 3x^2 - 3x^2 = 0 336 337 f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p 338 = <-sin t, cos t> + <cos t + x^3 + y^3, sin t + 2x^3 - 3x^2y> <<3x^2, 6x^2 - 6xy>, <3y^2, -3x^2>> - \nu <6x + 6y, 12x - 6y> + <3x, 3y> 339 = <-sin t, cos t> + <cos t 3x^2 + 3x^5 + 3x^2y^3 + sin t 3y^2 + 6x^3y^2 - 9x^2y^3, cos t (6x^2 - 6xy) + 6x^5 - 6x^4y + 6x^2y^3 - 6xy^4 + sin t (-3x^2) - 6x^5 + 9x^4y> - \nu <6x + 6y, 12x - 6y> + <3x, 3y> 340 = <cos t (3x^2) + sin t (3y^2 - 1) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu (6x + 6y) + 3x, 341 cos t (6x^2 - 6xy) - sin t (3x^2) + 3x^4y + 6x^2y^3 - 6xy^4 - \nu (12x - 6y) + 3y> 342 343 Q = dT/dt + u \cdot \nabla T - \alpha \Delta T 344 = -sin t + <cos t + x^3 + y^3, sin t + 2x^3 - 3x^2y> . <x, y> - 2 \alpha 345 = -sin t + cos t (x) + x^4 + xy^3 + sin t (y) + 2x^3y - 3x^2y^2 - 2 \alpha 346 = cos t x + sin t (y - 1) + (x^4 + 2x^3y - 3x^2y^2 + xy^3 - 2 \alpha) 347 */ 348 static PetscErrorCode cubic_trig_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 349 { 350 u[0] = 100.*PetscCosReal(time) + X[0]*X[0]*X[0] + X[1]*X[1]*X[1]; 351 u[1] = 100.*PetscSinReal(time) + 2.0*X[0]*X[0]*X[0] - 3.0*X[0]*X[0]*X[1]; 352 return 0; 353 } 354 static PetscErrorCode cubic_trig_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 355 { 356 u[0] = -100.*PetscSinReal(time); 357 u[1] = 100.*PetscCosReal(time); 358 return 0; 359 } 360 361 static PetscErrorCode cubic_trig_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 362 { 363 p[0] = 3.0*X[0]*X[0]/2.0 + 3.0*X[1]*X[1]/2.0 - 1.0; 364 return 0; 365 } 366 367 static PetscErrorCode cubic_trig_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 368 { 369 T[0] = 100.*PetscCosReal(time) + X[0]*X[0]/2.0 + X[1]*X[1]/2.0; 370 return 0; 371 } 372 static PetscErrorCode cubic_trig_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 373 { 374 T[0] = -100.*PetscSinReal(time); 375 return 0; 376 } 377 378 static void f0_cubic_trig_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 379 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 380 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 381 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 382 { 383 const PetscReal nu = PetscRealPart(constants[NU]); 384 385 f0[0] -= 100.*PetscCosReal(t)*(3*X[0]*X[0]) + 100.*PetscSinReal(t)*(3*X[1]*X[1] - 1.) + 3*X[0]*X[0]*X[0]*X[0]*X[0] + 6*X[0]*X[0]*X[0]*X[1]*X[1] - 6*X[0]*X[0]*X[1]*X[1]*X[1] - ( 6*X[0] + 6*X[1])*nu + 3*X[0]; 386 f0[1] -= 100.*PetscCosReal(t)*(6*X[0]*X[0] - 6*X[0]*X[1]) - 100.*PetscSinReal(t)*(3*X[0]*X[0]) + 3*X[0]*X[0]*X[0]*X[0]*X[1] + 6*X[0]*X[0]*X[1]*X[1]*X[1] - 6*X[0]*X[1]*X[1]*X[1]*X[1] - (12*X[0] - 6*X[1])*nu + 3*X[1]; 387 } 388 389 static void f0_cubic_trig_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 390 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 391 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 392 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 393 { 394 const PetscReal alpha = PetscRealPart(constants[ALPHA]); 395 396 f0[0] -= 100.*PetscCosReal(t)*X[0] + 100.*PetscSinReal(t)*(X[1] - 1.) + X[0]*X[0]*X[0]*X[0] + 2.0*X[0]*X[0]*X[0]*X[1] - 3.0*X[0]*X[0]*X[1]*X[1] + X[0]*X[1]*X[1]*X[1] - 2.0*alpha; 397 } 398 399 /* 400 CASE: Taylor-Green vortex 401 In 2D we use exact solution: 402 403 u = 1 - cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t) 404 v = 1 + sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t) 405 p = -1/4 [cos(2 \pi(x - t)) + cos(2 \pi(y - t))] exp(-4 \pi^2 \nu t) 406 T = t + x + y 407 f = <\nu \pi^2 exp(-2\nu \pi^2 t) cos(\pi(x-t)) sin(\pi(y-t)), -\nu \pi^2 exp(-2\nu \pi^2 t) sin(\pi(x-t)) cos(\pi(y-t)) > 408 Q = 3 + sin(\pi(x-y)) exp(-2\nu \pi^2 t) 409 410 so that 411 412 \nabla \cdot u = \pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t) - \pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t) = 0 413 414 f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p 415 = <-\pi (sin(\pi(x - t)) sin(\pi(y - t)) - cos(\pi(x - t)) cos(\pi(y - t)) - 2\pi cos(\pi(x - t)) sin(\pi(y - t))) exp(-2 \pi^2 \nu t), 416 \pi (sin(\pi(x - t)) sin(\pi(y - t)) - cos(\pi(x - t)) cos(\pi(y - t)) - 2\pi sin(\pi(x - t)) cos(\pi(y - t))) exp(-2 \pi^2 \nu t)> 417 + < \pi (1 - cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)) sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t), 418 \pi (1 - cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)) cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)> 419 + <-\pi (1 + sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)) cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t), 420 -\pi (1 + sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)) sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)> 421 + <-2\pi^2 cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t), 422 2\pi^2 sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)> 423 + < \pi/2 sin(2\pi(x - t)) exp(-4 \pi^2 \nu t), 424 \pi/2 sin(2\pi(y - t)) exp(-4 \pi^2 \nu t)> 425 = <-\pi (sin(\pi(x - t)) sin(\pi(y - t)) - cos(\pi(x - t)) cos(\pi(y - t)) - 2\pi cos(\pi(x - t)) sin(\pi(y - t))) exp(-2 \pi^2 \nu t), 426 \pi (sin(\pi(x - t)) sin(\pi(y - t)) - cos(\pi(x - t)) cos(\pi(y - t)) - 2\pi sin(\pi(x - t)) cos(\pi(y - t))) exp(-2 \pi^2 \nu t)> 427 + < \pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t), 428 \pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)> 429 + <-\pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t), 430 -\pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)> 431 + <-\pi/2 sin(2\pi(x - t)) exp(-4 \pi^2 \nu t), 432 -\pi/2 sin(2\pi(y - t)) exp(-4 \pi^2 \nu t)> 433 + <-2\pi^2 cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t), 434 2\pi^2 sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)> 435 + < \pi/2 sin(2\pi(x - t)) exp(-4 \pi^2 \nu t), 436 \pi/2 sin(2\pi(y - t)) exp(-4 \pi^2 \nu t)> 437 = <-\pi (sin(\pi(x - t)) sin(\pi(y - t)) - cos(\pi(x - t)) cos(\pi(y - t))) exp(-2 \pi^2 \nu t), 438 \pi (sin(\pi(x - t)) sin(\pi(y - t)) - cos(\pi(x - t)) cos(\pi(y - t))) exp(-2 \pi^2 \nu t)> 439 + < \pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t), 440 \pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)> 441 + <-\pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t), 442 -\pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)> 443 = < \pi cos(\pi(x - t)) cos(\pi(y - t)), 444 \pi sin(\pi(x - t)) sin(\pi(y - t))> 445 + <-\pi cos(\pi(x - t)) cos(\pi(y - t)), 446 -\pi sin(\pi(x - t)) sin(\pi(y - t))> = 0 447 Q = dT/dt + u \cdot \nabla T - \alpha \Delta T 448 = 1 + u \cdot <1, 1> - 0 449 = 1 + u + v 450 */ 451 452 static PetscErrorCode taylor_green_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 453 { 454 u[0] = 1 - PetscCosReal(PETSC_PI*(X[0]-time))*PetscSinReal(PETSC_PI*(X[1]-time))*PetscExpReal(-2*PETSC_PI*PETSC_PI*time); 455 u[1] = 1 + PetscSinReal(PETSC_PI*(X[0]-time))*PetscCosReal(PETSC_PI*(X[1]-time))*PetscExpReal(-2*PETSC_PI*PETSC_PI*time); 456 return 0; 457 } 458 static PetscErrorCode taylor_green_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 459 { 460 u[0] = -PETSC_PI*(PetscSinReal(PETSC_PI*(X[0]-time))*PetscSinReal(PETSC_PI*(X[1]-time)) 461 - PetscCosReal(PETSC_PI*(X[0]-time))*PetscCosReal(PETSC_PI*(X[1]-time)) 462 - 2*PETSC_PI*PetscCosReal(PETSC_PI*(X[0]-time))*PetscSinReal(PETSC_PI*(X[1]-time)))*PetscExpReal(-2*PETSC_PI*PETSC_PI*time); 463 u[1] = PETSC_PI*(PetscSinReal(PETSC_PI*(X[0]-time))*PetscSinReal(PETSC_PI*(X[1]-time)) 464 - PetscCosReal(PETSC_PI*(X[0]-time))*PetscCosReal(PETSC_PI*(X[1]-time)) 465 - 2*PETSC_PI*PetscSinReal(PETSC_PI*(X[0]-time))*PetscCosReal(PETSC_PI*(X[1]-time)))*PetscExpReal(-2*PETSC_PI*PETSC_PI*time); 466 return 0; 467 } 468 469 static PetscErrorCode taylor_green_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 470 { 471 p[0] = -0.25*(PetscCosReal(2*PETSC_PI*(X[0]-time)) + PetscCosReal(2*PETSC_PI*(X[1]-time)))*PetscExpReal(-4*PETSC_PI*PETSC_PI*time); 472 return 0; 473 } 474 475 static PetscErrorCode taylor_green_p_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 476 { 477 p[0] = PETSC_PI*(0.5*(PetscSinReal(2*PETSC_PI*(X[0]-time)) + PetscSinReal(2*PETSC_PI*(X[1]-time))) 478 + PETSC_PI*(PetscCosReal(2*PETSC_PI*(X[0]-time)) + PetscCosReal(2*PETSC_PI*(X[1]-time))))*PetscExpReal(-4*PETSC_PI*PETSC_PI*time); 479 return 0; 480 } 481 482 static PetscErrorCode taylor_green_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 483 { 484 T[0] = time + X[0] + X[1]; 485 return 0; 486 } 487 static PetscErrorCode taylor_green_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 488 { 489 T[0] = 1.0; 490 return 0; 491 } 492 493 static void f0_taylor_green_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 494 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 495 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 496 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 497 { 498 PetscScalar vel[2]; 499 500 taylor_green_u(dim, t, X, Nf, vel, NULL); 501 f0[0] -= 1.0 + vel[0] + vel[1]; 502 } 503 504 /* 505 CASE: Pipe flow 506 Poiseuille flow, with the incoming fluid having a parabolic temperature profile and the side walls being held at T_in 507 508 u = \Delta Re/(2 mu) y (1 - y) 509 v = 0 510 p = -\Delta x 511 T = y (1 - y) + T_in 512 rho = p^{th} / T 513 514 so that 515 516 \nabla \cdot u = 0 - 0 = 0 517 grad u = \Delta Re/(2 mu) <<0, 0>, <1 - 2y, 0>> 518 epsilon(u) = 1/2 (grad u + grad u^T) = \Delta Re/(4 mu) <<0, 1 - 2y>, <<1 - 2y, 0>> 519 epsilon'(u) = epsilon(u) - 1/3 (div u) I = epsilon(u) 520 div epsilon'(u) = -\Delta Re/(2 mu) <1, 0> 521 522 f = rho S du/dt + rho u \cdot \nabla u - 2\mu/Re div \epsilon'(u) + \nabla p + rho / F^2 \hat y 523 = 0 + 0 - div (2\mu/Re \epsilon'(u) - pI) + rho / F^2 \hat y 524 = -\Delta div <<x, (1 - 2y)/2>, <<(1 - 2y)/2, x>> + rho / F^2 \hat y 525 = \Delta <1, 0> - \Delta <1, 0> + rho/F^2 <0, 1> 526 = rho/F^2 <0, 1> 527 528 g = S rho_t + div (rho u) 529 = 0 + rho div (u) + u . grad rho 530 = 0 + 0 - pth u . grad T / T^2 531 = 0 532 533 Q = rho c_p S dT/dt + rho c_p u . grad T - 1/Pe div k grad T 534 = 0 + c_p pth / T 0 + 2 k/Pe 535 = 2 k/Pe 536 537 The boundary conditions on the top and bottom are zero velocity and T_in temperature. The boundary term is 538 539 (2\mu/Re \epsilon'(u) - p I) . n = \Delta <<x, (1 - 2y)/2>, <<(1 - 2y)/2, x>> . n 540 541 so that 542 543 x = 0: \Delta <<0, (1 - 2y)/2>, <<(1 - 2y)/2, 0>> . <-1, 0> = <0, (2y - 1)/2> 544 x = 1: \Delta <<1, (1 - 2y)/2>, <<(1 - 2y)/2, 1>> . <1, 0> = <1, (1 - 2y)/2> 545 */ 546 static PetscErrorCode pipe_u(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 547 { 548 Parameter *param = (Parameter *) ctx; 549 550 u[0] = (0.5*param->Reynolds / param->mu) * X[1]*(1.0 - X[1]); 551 u[1] = 0.0; 552 return 0; 553 } 554 static PetscErrorCode pipe_u_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 555 { 556 u[0] = 0.0; 557 u[1] = 0.0; 558 return 0; 559 } 560 561 static PetscErrorCode pipe_p(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 562 { 563 p[0] = -X[0]; 564 return 0; 565 } 566 static PetscErrorCode pipe_p_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 567 { 568 p[0] = 0.0; 569 return 0; 570 } 571 572 static PetscErrorCode pipe_T(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 573 { 574 Parameter *param = (Parameter *) ctx; 575 576 T[0] = X[1]*(1.0 - X[1]) + param->T_in; 577 return 0; 578 } 579 static PetscErrorCode pipe_T_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 580 { 581 T[0] = 0.0; 582 return 0; 583 } 584 585 static void f0_conduct_pipe_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 586 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 587 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 588 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 589 { 590 const PetscReal F = PetscRealPart(constants[FROUDE]); 591 const PetscReal p_th = PetscRealPart(constants[P_TH]); 592 const PetscReal T_in = PetscRealPart(constants[T_IN]); 593 const PetscReal rho = p_th / (X[1]*(1. - X[1]) + T_in); 594 const PetscInt gd = (PetscInt) PetscRealPart(constants[G_DIR]); 595 596 f0[gd] -= rho/PetscSqr(F); 597 } 598 599 static void f0_conduct_bd_pipe_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 600 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 601 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 602 PetscReal t, const PetscReal X[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 603 { 604 PetscReal sigma[4] = {X[0], 0.5*(1. - 2.*X[1]), 0.5*(1. - 2.*X[1]), X[0]}; 605 PetscInt d, e; 606 607 for (d = 0; d < dim; ++d) { 608 for (e = 0; e < dim; ++e) { 609 f0[d] -= sigma[d*dim + e] * n[e]; 610 } 611 } 612 } 613 614 static void f0_conduct_pipe_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 615 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 616 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 617 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 618 { 619 f0[0] += 0.0; 620 } 621 622 static void f0_conduct_pipe_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 623 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 624 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 625 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 626 { 627 const PetscReal k = PetscRealPart(constants[K]); 628 const PetscReal Pe = PetscRealPart(constants[PECLET]); 629 630 f0[0] -= 2*k/Pe; 631 } 632 633 /* 634 CASE: Wiggly pipe flow 635 Perturbed Poiseuille flow, with the incoming fluid having a perturbed parabolic temperature profile and the side walls being held at T_in 636 637 u = \Delta Re/(2 mu) [y (1 - y) + a sin(pi y)] 638 v = 0 639 p = -\Delta x 640 T = y (1 - y) + a sin(pi y) + T_in 641 rho = p^{th} / T 642 643 so that 644 645 \nabla \cdot u = 0 - 0 = 0 646 grad u = \Delta Re/(2 mu) <<0, 0>, <1 - 2y + a pi cos(pi y), 0>> 647 epsilon(u) = 1/2 (grad u + grad u^T) = \Delta Re/(4 mu) <<0, 1 - 2y + a pi cos(pi y)>, <<1 - 2y + a pi cos(pi y), 0>> 648 epsilon'(u) = epsilon(u) - 1/3 (div u) I = epsilon(u) 649 div epsilon'(u) = -\Delta Re/(2 mu) <1 + a pi^2/2 sin(pi y), 0> 650 651 f = rho S du/dt + rho u \cdot \nabla u - 2\mu/Re div \epsilon'(u) + \nabla p + rho / F^2 \hat y 652 = 0 + 0 - div (2\mu/Re \epsilon'(u) - pI) + rho / F^2 \hat y 653 = -\Delta div <<x, (1 - 2y)/2 + a pi/2 cos(pi y)>, <<(1 - 2y)/2 + a pi/2 cos(pi y), x>> + rho / F^2 \hat y 654 = -\Delta <1 - 1 - a pi^2/2 sin(pi y), 0> + rho/F^2 <0, 1> 655 = a \Delta pi^2/2 sin(pi y) <1, 0> + rho/F^2 <0, 1> 656 657 g = S rho_t + div (rho u) 658 = 0 + rho div (u) + u . grad rho 659 = 0 + 0 - pth u . grad T / T^2 660 = 0 661 662 Q = rho c_p S dT/dt + rho c_p u . grad T - 1/Pe div k grad T 663 = 0 + c_p pth / T 0 - k/Pe div <0, 1 - 2y + a pi cos(pi y)> 664 = - k/Pe (-2 - a pi^2 sin(pi y)) 665 = 2 k/Pe (1 + a pi^2/2 sin(pi y)) 666 667 The boundary conditions on the top and bottom are zero velocity and T_in temperature. The boundary term is 668 669 (2\mu/Re \epsilon'(u) - p I) . n = \Delta <<x, (1 - 2y)/2 + a pi/2 cos(pi y)>, <<(1 - 2y)/2 + a pi/2 cos(pi y), x>> . n 670 671 so that 672 673 x = 0: \Delta <<0, (1 - 2y)/2>, <<(1 - 2y)/2, 0>> . <-1, 0> = <0, (2y - 1)/2 - a pi/2 cos(pi y)> 674 x = 1: \Delta <<1, (1 - 2y)/2>, <<(1 - 2y)/2, 1>> . < 1, 0> = <1, (1 - 2y)/2 + a pi/2 cos(pi y)> 675 */ 676 static PetscErrorCode pipe_wiggly_u(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 677 { 678 Parameter *param = (Parameter *) ctx; 679 680 u[0] = (0.5*param->Reynolds / param->mu) * (X[1]*(1.0 - X[1]) + param->epsilon*PetscSinReal(PETSC_PI*X[1])); 681 u[1] = 0.0; 682 return 0; 683 } 684 static PetscErrorCode pipe_wiggly_u_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 685 { 686 u[0] = 0.0; 687 u[1] = 0.0; 688 return 0; 689 } 690 691 static PetscErrorCode pipe_wiggly_p(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 692 { 693 p[0] = -X[0]; 694 return 0; 695 } 696 static PetscErrorCode pipe_wiggly_p_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 697 { 698 p[0] = 0.0; 699 return 0; 700 } 701 702 static PetscErrorCode pipe_wiggly_T(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 703 { 704 Parameter *param = (Parameter *) ctx; 705 706 T[0] = X[1]*(1.0 - X[1]) + param->epsilon*PetscSinReal(PETSC_PI*X[1]) + param->T_in; 707 return 0; 708 } 709 static PetscErrorCode pipe_wiggly_T_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 710 { 711 T[0] = 0.0; 712 return 0; 713 } 714 715 static void f0_conduct_pipe_wiggly_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 716 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 717 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 718 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 719 { 720 const PetscReal F = PetscRealPart(constants[FROUDE]); 721 const PetscReal p_th = PetscRealPart(constants[P_TH]); 722 const PetscReal T_in = PetscRealPart(constants[T_IN]); 723 const PetscReal eps = PetscRealPart(constants[EPSILON]); 724 const PetscReal rho = p_th / (X[1]*(1. - X[1]) + T_in); 725 const PetscInt gd = (PetscInt) PetscRealPart(constants[G_DIR]); 726 727 f0[0] -= eps*0.5*PetscSqr(PETSC_PI)*PetscSinReal(PETSC_PI*X[1]); 728 f0[gd] -= rho/PetscSqr(F); 729 } 730 731 static void f0_conduct_bd_pipe_wiggly_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 732 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 733 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 734 PetscReal t, const PetscReal X[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 735 { 736 const PetscReal eps = PetscRealPart(constants[EPSILON]); 737 PetscReal sigma[4] = {X[0], 0.5*(1. - 2.*X[1]) + eps*0.5*PETSC_PI*PetscCosReal(PETSC_PI*X[1]), 0.5*(1. - 2.*X[1]) + eps*0.5*PETSC_PI*PetscCosReal(PETSC_PI*X[1]), X[0]}; 738 PetscInt d, e; 739 740 for (d = 0; d < dim; ++d) { 741 for (e = 0; e < dim; ++e) { 742 f0[d] -= sigma[d*dim + e] * n[e]; 743 } 744 } 745 } 746 747 static void f0_conduct_pipe_wiggly_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 748 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 749 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 750 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 751 { 752 f0[0] += 0.0; 753 } 754 755 static void f0_conduct_pipe_wiggly_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 756 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 757 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 758 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 759 { 760 const PetscReal k = PetscRealPart(constants[K]); 761 const PetscReal Pe = PetscRealPart(constants[PECLET]); 762 const PetscReal eps = PetscRealPart(constants[EPSILON]); 763 764 f0[0] -= 2*k/Pe*(1.0 + eps*0.5*PetscSqr(PETSC_PI)*PetscSinReal(PETSC_PI*X[1])); 765 } 766 767 /* Physics Kernels */ 768 769 static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 770 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 771 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 772 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 773 { 774 PetscInt d; 775 for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d]; 776 } 777 778 /* -\frac{Sp^{th}}{T^2} \frac{\partial T}{\partial t} + \frac{p^{th}}{T} \nabla \cdot \vb{u} - \frac{p^{th}}{T^2} \vb{u} \cdot \nabla T */ 779 static void f0_conduct_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 780 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 781 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 782 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 783 { 784 const PetscReal S = PetscRealPart(constants[STROUHAL]); 785 const PetscReal p_th = PetscRealPart(constants[P_TH]); 786 PetscInt d; 787 788 // -\frac{S p^{th}}{T^2} \frac{\partial T}{\partial t} 789 f0[0] += -u_t[uOff[TEMP]] * S * p_th / PetscSqr(u[uOff[TEMP]]); 790 791 // \frac{p^{th}}{T} \nabla \cdot \vb{u} 792 for (d = 0; d < dim; ++d) { 793 f0[0] += p_th / u[uOff[TEMP]] * u_x[uOff_x[VEL] + d*dim + d]; 794 } 795 796 // - \frac{p^{th}}{T^2} \vb{u} \cdot \nabla T 797 for (d = 0; d < dim; ++d) { 798 f0[0] -= p_th / (u[uOff[TEMP]] * u[uOff[TEMP]]) * u[uOff[VEL] + d] * u_x[uOff_x[TEMP] + d]; 799 } 800 801 // Add in any fixed source term 802 if (NfAux > 0) { 803 f0[0] += a[aOff[MASS]]; 804 } 805 } 806 807 /* \vb{u}_t + \vb{u} \cdot \nabla\vb{u} */ 808 static void f0_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 809 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 810 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 811 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 812 { 813 const PetscInt Nc = dim; 814 PetscInt c, d; 815 816 for (c = 0; c < Nc; ++c) { 817 /* \vb{u}_t */ 818 f0[c] += u_t[uOff[VEL] + c]; 819 /* \vb{u} \cdot \nabla\vb{u} */ 820 for (d = 0; d < dim; ++d) f0[c] += u[uOff[VEL] + d]*u_x[uOff_x[VEL] + c*dim + d]; 821 } 822 } 823 824 /* \rho S \frac{\partial \vb{u}}{\partial t} + \rho \vb{u} \cdot \nabla \vb{u} + \rho \frac{\hat{\vb{z}}}{F^2} */ 825 static void f0_conduct_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 826 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 827 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 828 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 829 { 830 const PetscReal S = PetscRealPart(constants[STROUHAL]); 831 const PetscReal F = PetscRealPart(constants[FROUDE]); 832 const PetscReal p_th = PetscRealPart(constants[P_TH]); 833 const PetscReal rho = p_th / PetscRealPart(u[uOff[TEMP]]); 834 const PetscInt gdir = (PetscInt) PetscRealPart(constants[G_DIR]); 835 PetscInt Nc = dim; 836 PetscInt c, d; 837 838 // \rho S \frac{\partial \vb{u}}{\partial t} 839 for (d = 0; d < dim; ++d) { 840 f0[d] = rho * S * u_t[uOff[VEL] + d]; 841 } 842 843 // \rho \vb{u} \cdot \nabla \vb{u} 844 for (c = 0; c < Nc; ++c) { 845 for (d = 0; d < dim; ++d) { 846 f0[c] += rho * u[uOff[VEL] + d] * u_x[uOff_x[VEL] + c*dim + d]; 847 } 848 } 849 850 // rho \hat{z}/F^2 851 f0[gdir] += rho / (F*F); 852 853 // Add in any fixed source term 854 if (NfAux > 0) { 855 for (d = 0; d < dim; ++d) { 856 f0[d] += a[aOff[MOMENTUM] + d]; 857 } 858 } 859 } 860 861 /*f1_v = \nu[grad(u) + grad(u)^T] - pI */ 862 static void f1_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 863 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 864 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 865 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 866 { 867 const PetscReal nu = PetscRealPart(constants[NU]); 868 const PetscInt Nc = dim; 869 PetscInt c, d; 870 871 for (c = 0; c < Nc; ++c) { 872 for (d = 0; d < dim; ++d) { 873 f1[c*dim+d] = nu*(u_x[c*dim+d] + u_x[d*dim+c]); 874 } 875 f1[c*dim+c] -= u[uOff[1]]; 876 } 877 } 878 879 /* 2 \mu/Re (1/2 (\nabla \vb{u} + \nabla \vb{u}^T) - 1/3 (\nabla \cdot \vb{u}) I) - p I */ 880 static void f1_conduct_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 881 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 882 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 883 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 884 { 885 const PetscReal Re = PetscRealPart(constants[REYNOLDS]); 886 const PetscReal mu = PetscRealPart(constants[MU]); 887 const PetscReal coef = mu / Re; 888 PetscReal u_div = 0.0; 889 const PetscInt Nc = dim; 890 PetscInt c, d; 891 892 for (c = 0; c < Nc; ++c) { 893 u_div += PetscRealPart(u_x[uOff_x[VEL] + c*dim + c]); 894 } 895 896 for (c = 0; c < Nc; ++c) { 897 // 2 \mu/Re 1/2 (\nabla \vb{u} + \nabla \vb{u}^T 898 for (d = 0; d < dim; ++d) { 899 f1[c*dim + d] += coef * (u_x[uOff_x[VEL] + c*dim + d] + u_x[uOff_x[VEL] + d*dim + c]); 900 } 901 // -2/3 \mu/Re (\nabla \cdot \vb{u}) I 902 f1[c * dim + c] -= 2.0 * coef / 3.0 * u_div; 903 } 904 905 // -p I 906 for (c = 0; c < Nc; ++c) { 907 f1[c*dim + c] -= u[uOff[PRES]]; 908 } 909 } 910 911 /* T_t + \vb{u} \cdot \nabla T */ 912 static void f0_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 913 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 914 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 915 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 916 { 917 PetscInt d; 918 919 /* T_t */ 920 f0[0] += u_t[uOff[TEMP]]; 921 /* \vb{u} \cdot \nabla T */ 922 for (d = 0; d < dim; ++d) f0[0] += u[uOff[VEL] + d] * u_x[uOff_x[TEMP] + d]; 923 } 924 925 /* \frac{C_p S p^{th}}{T} \frac{\partial T}{\partial t} + \frac{C_p p^{th}}{T} \vb{u} \cdot \nabla T */ 926 static void f0_conduct_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 927 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 928 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 929 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 930 { 931 const PetscReal S = PetscRealPart(constants[STROUHAL]); 932 const PetscReal c_p = PetscRealPart(constants[C_P]); 933 const PetscReal p_th = PetscRealPart(constants[P_TH]); 934 PetscInt d; 935 936 // \frac{C_p S p^{th}}{T} \frac{\partial T}{\partial t} 937 f0[0] = c_p * S * p_th / u[uOff[TEMP]] * u_t[uOff[TEMP]]; 938 939 // \frac{C_p p^{th}}{T} \vb{u} \cdot \nabla T 940 for (d = 0; d < dim; ++d) { 941 f0[0] += c_p * p_th / u[uOff[TEMP]] * u[uOff[VEL] + d] * u_x[uOff_x[TEMP] + d]; 942 } 943 944 // Add in any fixed source term 945 if (NfAux > 0) { 946 f0[0] += a[aOff[ENERGY]]; 947 } 948 } 949 950 static void f1_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 951 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 952 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 953 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 954 { 955 const PetscReal alpha = PetscRealPart(constants[ALPHA]); 956 PetscInt d; 957 958 for (d = 0; d < dim; ++d) f1[d] = alpha*u_x[uOff_x[2]+d]; 959 } 960 961 /* \frac{k}{Pe} \nabla T */ 962 static void f1_conduct_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 963 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 964 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 965 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 966 { 967 const PetscReal Pe = PetscRealPart(constants[PECLET]); 968 const PetscReal k = PetscRealPart(constants[K]); 969 PetscInt d; 970 971 // \frac{k}{Pe} \nabla T 972 for (d = 0; d < dim; ++d) { 973 f1[d] = k / Pe * u_x[uOff_x[TEMP] + d]; 974 } 975 } 976 977 static void g1_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 978 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 979 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 980 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 981 { 982 PetscInt d; 983 for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; 984 } 985 986 static void g0_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 987 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 988 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 989 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 990 { 991 PetscInt c, d; 992 const PetscInt Nc = dim; 993 994 for (d = 0; d < dim; ++d) g0[d*dim+d] = u_tShift; 995 996 for (c = 0; c < Nc; ++c) { 997 for (d = 0; d < dim; ++d) { 998 g0[c*Nc+d] += u_x[ c*Nc+d]; 999 } 1000 } 1001 } 1002 1003 static void g1_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1004 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1005 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1006 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 1007 { 1008 PetscInt NcI = dim; 1009 PetscInt NcJ = dim; 1010 PetscInt c, d, e; 1011 1012 for (c = 0; c < NcI; ++c) { 1013 for (d = 0; d < NcJ; ++d) { 1014 for (e = 0; e < dim; ++e) { 1015 if (c == d) { 1016 g1[(c*NcJ+d)*dim+e] += u[e]; 1017 } 1018 } 1019 } 1020 } 1021 } 1022 1023 static void g0_conduct_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1024 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1025 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1026 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 1027 { 1028 const PetscReal p_th = PetscRealPart(constants[P_TH]); 1029 PetscInt d; 1030 1031 // - \phi_i \frac{p^{th}}{T^2} \frac{\partial T}{\partial x_c} \psi_{j, u_c} 1032 for (d = 0; d < dim; ++d) { 1033 g0[d] = -p_th / PetscSqr(u[uOff[TEMP]]) * u_x[uOff_x[TEMP] + d]; 1034 } 1035 } 1036 1037 static void g1_conduct_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1038 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1039 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1040 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 1041 { 1042 const PetscReal p_th = PetscRealPart(constants[P_TH]); 1043 PetscInt d; 1044 1045 // \phi_i \frac{p^{th}}{T} \frac{\partial \psi_{u_c,j}}{\partial x_c} 1046 for (d = 0; d < dim; ++d) { 1047 g1[d * dim + d] = p_th / u[uOff[TEMP]]; 1048 } 1049 } 1050 1051 static void g0_conduct_qT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1052 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1053 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1054 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 1055 { 1056 const PetscReal S = PetscRealPart(constants[STROUHAL]); 1057 const PetscReal p_th = PetscRealPart(constants[P_TH]); 1058 PetscInt d; 1059 1060 // - \phi_i \frac{S p^{th}}{T^2} \psi_j 1061 g0[0] -= S * p_th / PetscSqr(u[uOff[TEMP]]) * u_tShift; 1062 // \phi_i 2 \frac{S p^{th}}{T^3} T_t \psi_j 1063 g0[0] += 2.0 * S * p_th / PetscPowScalarInt(u[uOff[TEMP]], 3) * u_t[uOff[TEMP]]; 1064 // \phi_i \frac{p^{th}}{T^2} \left( - \nabla \cdot \vb{u} \psi_j + \frac{2}{T} \vb{u} \cdot \nabla T \psi_j \right) 1065 for (d = 0; d < dim; ++d) { 1066 g0[0] += p_th / PetscSqr(u[uOff[TEMP]]) * (-u_x[uOff_x[VEL] + d * dim + d] + 2.0 / u[uOff[TEMP]] * u[uOff[VEL] + d] * u_x[uOff_x[TEMP] + d]); 1067 } 1068 } 1069 1070 static void g1_conduct_qT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1071 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1072 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1073 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 1074 { 1075 const PetscReal p_th = PetscRealPart(constants[P_TH]); 1076 PetscInt d; 1077 1078 // - \phi_i \frac{p^{th}}{T^2} \vb{u} \cdot \nabla \psi_j 1079 for (d = 0; d < dim; ++d) { 1080 g1[d] = -p_th / PetscSqr(u[uOff[TEMP]]) * u[uOff[VEL] + d]; 1081 } 1082 } 1083 1084 static void g2_vp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1085 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1086 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1087 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 1088 { 1089 PetscInt d; 1090 for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; 1091 } 1092 1093 static void g3_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1094 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1095 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1096 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 1097 { 1098 const PetscReal nu = PetscRealPart(constants[NU]); 1099 const PetscInt Nc = dim; 1100 PetscInt c, d; 1101 1102 for (c = 0; c < Nc; ++c) { 1103 for (d = 0; d < dim; ++d) { 1104 g3[((c*Nc+c)*dim+d)*dim+d] += nu; 1105 g3[((c*Nc+d)*dim+d)*dim+c] += nu; 1106 } 1107 } 1108 } 1109 1110 static void g0_conduct_vT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1111 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1112 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1113 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 1114 { 1115 const PetscReal S = PetscRealPart(constants[STROUHAL]); 1116 const PetscReal F = PetscRealPart(constants[FROUDE]); 1117 const PetscReal p_th = PetscRealPart(constants[P_TH]); 1118 const PetscInt gdir = (PetscInt) PetscRealPart(constants[G_DIR]); 1119 const PetscInt Nc = dim; 1120 PetscInt c, d; 1121 1122 // - \vb{\phi}_i \cdot \vb{u}_t \frac{p^{th} S}{T^2} \psi_j 1123 for (d = 0; d < dim; ++d) { 1124 g0[d] -= p_th * S / PetscSqr(u[uOff[TEMP]]) * u_t[uOff[VEL] + d]; 1125 } 1126 1127 // - \vb{\phi}_i \cdot \vb{u} \cdot \nabla \vb{u} \frac{p^{th}}{T^2} \psi_j 1128 for (c = 0; c < Nc; ++c) { 1129 for (d = 0; d < dim; ++d) { 1130 g0[c] -= p_th / PetscSqr(u[uOff[TEMP]]) * u[uOff[VEL] + d] * u_x[uOff_x[VEL] + c * dim + d]; 1131 } 1132 } 1133 1134 // - \vb{\phi}_i \cdot \vu{z} \frac{p^{th}}{T^2 F^2} \psi_j 1135 g0[gdir] -= p_th / PetscSqr(u[uOff[TEMP]] * F); 1136 } 1137 1138 static void g0_conduct_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1139 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1140 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1141 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 1142 { 1143 const PetscReal S = PetscRealPart(constants[STROUHAL]); 1144 const PetscReal p_th = PetscRealPart(constants[P_TH]); 1145 const PetscInt Nc = dim; 1146 PetscInt c, d; 1147 1148 // \vb{\phi}_i \cdot S \rho \psi_j 1149 for (d = 0; d < dim; ++d) { 1150 g0[d * dim + d] = S * p_th / u[uOff[TEMP]] * u_tShift; 1151 } 1152 1153 // \phi^c_i \cdot \rho \frac{\partial u^c}{\partial x^d} \psi^d_j 1154 for (c = 0; c < Nc; ++c) { 1155 for (d = 0; d < dim; ++d) { 1156 g0[c * Nc + d] += p_th / u[uOff[TEMP]] * u_x[uOff_x[VEL] + c * Nc + d]; 1157 } 1158 } 1159 } 1160 1161 static void g1_conduct_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1162 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1163 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1164 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 1165 { 1166 const PetscReal p_th = PetscRealPart(constants[P_TH]); 1167 const PetscInt NcI = dim; 1168 const PetscInt NcJ = dim; 1169 PetscInt c, d, e; 1170 1171 // \phi^c_i \rho u^e \frac{\partial \psi^d_j}{\partial x^e} 1172 for (c = 0; c < NcI; ++c) { 1173 for (d = 0; d < NcJ; ++d) { 1174 for (e = 0; e < dim; ++e) { 1175 if (c == d) { 1176 g1[(c * NcJ + d) * dim + e] += p_th / u[uOff[TEMP]] * u[uOff[VEL] + e]; 1177 } 1178 } 1179 } 1180 } 1181 } 1182 1183 static void g3_conduct_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1184 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1185 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1186 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 1187 { 1188 const PetscReal Re = PetscRealPart(constants[REYNOLDS]); 1189 const PetscReal mu = PetscRealPart(constants[MU]); 1190 const PetscInt Nc = dim; 1191 PetscInt c, d; 1192 1193 for (c = 0; c < Nc; ++c) { 1194 for (d = 0; d < dim; ++d) { 1195 // \frac{\partial \phi^c_i}{\partial x^d} \mu/Re \frac{\partial \psi^c_i}{\partial x^d} 1196 g3[((c * Nc + c) * dim + d) * dim + d] += mu / Re; // gradU 1197 // \frac{\partial \phi^c_i}{\partial x^d} \mu/Re \frac{\partial \psi^d_i}{\partial x^c} 1198 g3[((c * Nc + d) * dim + d) * dim + c] += mu / Re; // gradU transpose 1199 // \frac{\partial \phi^c_i}{\partial x^d} -2/3 \mu/Re \frac{\partial \psi^d_i}{\partial x^c} 1200 g3[((c * Nc + d) * dim + c) * dim + d] -= 2.0 / 3.0 * mu / Re; 1201 } 1202 } 1203 } 1204 1205 static void g2_conduct_vp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1206 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1207 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1208 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 1209 { 1210 PetscInt d; 1211 for (d = 0; d < dim; ++d) { 1212 g2[d * dim + d] = -1.0; 1213 } 1214 } 1215 1216 static void g0_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1217 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1218 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1219 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 1220 { 1221 g0[0] = u_tShift; 1222 } 1223 1224 static void g0_wu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1225 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1226 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1227 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 1228 { 1229 PetscInt d; 1230 for (d = 0; d < dim; ++d) g0[d] = u_x[uOff_x[2]+d]; 1231 } 1232 1233 static void g1_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1234 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1235 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1236 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 1237 { 1238 PetscInt d; 1239 for (d = 0; d < dim; ++d) g1[d] = u[uOff[0]+d]; 1240 } 1241 1242 static void g3_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1243 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1244 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1245 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 1246 { 1247 const PetscReal alpha = PetscRealPart(constants[ALPHA]); 1248 PetscInt d; 1249 1250 for (d = 0; d < dim; ++d) g3[d*dim+d] = alpha; 1251 } 1252 1253 static void g0_conduct_wu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1254 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1255 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1256 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 1257 { 1258 const PetscReal p_th = PetscRealPart(constants[P_TH]); 1259 const PetscReal c_p = PetscRealPart(constants[C_P]); 1260 PetscInt d; 1261 1262 // \phi_i \frac{C_p p^{th}}{T} \nabla T \cdot \psi_j 1263 for (d = 0; d < dim; ++d) { 1264 g0[d] = c_p * p_th / u[uOff[TEMP]] * u_x[uOff_x[TEMP] + d]; 1265 } 1266 } 1267 1268 static void g0_conduct_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1269 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1270 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1271 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 1272 { 1273 const PetscReal S = PetscRealPart(constants[STROUHAL]); 1274 const PetscReal p_th = PetscRealPart(constants[P_TH]); 1275 const PetscReal c_p = PetscRealPart(constants[C_P]); 1276 PetscInt d; 1277 1278 // \psi_i C_p S p^{th}\T \psi_{j} 1279 g0[0] += c_p * S * p_th / u[uOff[TEMP]] * u_tShift; 1280 // - \phi_i C_p S p^{th}/T^2 T_t \psi_j 1281 g0[0] -= c_p * S * p_th / PetscSqr(u[uOff[TEMP]]) * u_t[uOff[TEMP]]; 1282 // - \phi_i C_p p^{th}/T^2 \vb{u} \cdot \nabla T \psi_j 1283 for (d = 0; d < dim; ++d) { 1284 g0[0] -= c_p * p_th / PetscSqr(u[uOff[TEMP]]) * u[uOff[VEL] + d] * u_x[uOff_x[TEMP] + d]; 1285 } 1286 } 1287 1288 static void g1_conduct_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1289 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1290 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1291 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 1292 { 1293 const PetscReal p_th = PetscRealPart(constants[P_TH]); 1294 const PetscReal c_p = PetscRealPart(constants[C_P]); 1295 PetscInt d; 1296 1297 // \phi_i C_p p^{th}/T \vb{u} \cdot \nabla \psi_j 1298 for (d = 0; d < dim; ++d) { 1299 g1[d] += c_p * p_th / u[uOff[TEMP]] * u[uOff[VEL] + d]; 1300 } 1301 } 1302 1303 static void g3_conduct_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1304 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1305 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1306 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 1307 { 1308 const PetscReal Pe = PetscRealPart(constants[PECLET]); 1309 const PetscReal k = PetscRealPart(constants[K]); 1310 PetscInt d; 1311 1312 // \nabla \phi_i \frac{k}{Pe} \nabla \phi_j 1313 for (d = 0; d < dim; ++d) { 1314 g3[d * dim + d] = k / Pe; 1315 } 1316 } 1317 1318 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 1319 { 1320 PetscInt mod, sol; 1321 PetscErrorCode ierr; 1322 1323 PetscFunctionBeginUser; 1324 options->modType = MOD_INCOMPRESSIBLE; 1325 options->solType = SOL_QUADRATIC; 1326 options->hasNullSpace = PETSC_TRUE; 1327 options->dmCell = NULL; 1328 1329 ierr = PetscOptionsBegin(comm, "", "Low Mach flow Problem Options", "DMPLEX");CHKERRQ(ierr); 1330 mod = options->modType; 1331 ierr = PetscOptionsEList("-mod_type", "The model type", "ex76.c", modTypes, NUM_MOD_TYPES, modTypes[options->modType], &mod, NULL);CHKERRQ(ierr); 1332 options->modType = (ModType) mod; 1333 sol = options->solType; 1334 ierr = PetscOptionsEList("-sol_type", "The solution type", "ex76.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL);CHKERRQ(ierr); 1335 options->solType = (SolType) sol; 1336 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1337 PetscFunctionReturn(0); 1338 } 1339 1340 static PetscErrorCode SetupParameters(DM dm, AppCtx *user) 1341 { 1342 PetscBag bag; 1343 Parameter *p; 1344 PetscReal dir; 1345 PetscInt dim; 1346 PetscErrorCode ierr; 1347 1348 PetscFunctionBeginUser; 1349 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1350 dir = (PetscReal) (dim-1); 1351 /* setup PETSc parameter bag */ 1352 ierr = PetscBagGetData(user->bag, (void **) &p);CHKERRQ(ierr); 1353 ierr = PetscBagSetName(user->bag, "par", "Low Mach flow parameters");CHKERRQ(ierr); 1354 bag = user->bag; 1355 ierr = PetscBagRegisterReal(bag, &p->Strouhal, 1.0, "S", "Strouhal number");CHKERRQ(ierr); 1356 ierr = PetscBagRegisterReal(bag, &p->Froude, 1.0, "Fr", "Froude number");CHKERRQ(ierr); 1357 ierr = PetscBagRegisterReal(bag, &p->Reynolds, 1.0, "Re", "Reynolds number");CHKERRQ(ierr); 1358 ierr = PetscBagRegisterReal(bag, &p->Peclet, 1.0, "Pe", "Peclet number");CHKERRQ(ierr); 1359 ierr = PetscBagRegisterReal(bag, &p->p_th, 1.0, "p_th", "Thermodynamic pressure");CHKERRQ(ierr); 1360 ierr = PetscBagRegisterReal(bag, &p->mu, 1.0, "mu", "Dynamic viscosity");CHKERRQ(ierr); 1361 ierr = PetscBagRegisterReal(bag, &p->nu, 1.0, "nu", "Kinematic viscosity");CHKERRQ(ierr); 1362 ierr = PetscBagRegisterReal(bag, &p->c_p, 1.0, "c_p", "Specific heat at constant pressure");CHKERRQ(ierr); 1363 ierr = PetscBagRegisterReal(bag, &p->k, 1.0, "k", "Thermal conductivity");CHKERRQ(ierr); 1364 ierr = PetscBagRegisterReal(bag, &p->alpha, 1.0, "alpha", "Thermal diffusivity");CHKERRQ(ierr); 1365 ierr = PetscBagRegisterReal(bag, &p->T_in, 1.0, "T_in", "Inlet temperature");CHKERRQ(ierr); 1366 ierr = PetscBagRegisterReal(bag, &p->g_dir, dir, "g_dir", "Gravity direction");CHKERRQ(ierr); 1367 ierr = PetscBagRegisterReal(bag, &p->epsilon, 1.0, "epsilon", "Perturbation strength");CHKERRQ(ierr); 1368 PetscFunctionReturn(0); 1369 } 1370 1371 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 1372 { 1373 PetscErrorCode ierr; 1374 1375 PetscFunctionBeginUser; 1376 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1377 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1378 ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 1379 ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 1380 PetscFunctionReturn(0); 1381 } 1382 1383 static PetscErrorCode UniformBoundaryConditions(DM dm, DMLabel label, PetscSimplePointFunc exactFuncs[], PetscSimplePointFunc exactFuncs_t[], AppCtx *user) 1384 { 1385 PetscDS ds; 1386 PetscInt id; 1387 void *ctx; 1388 PetscErrorCode ierr; 1389 1390 PetscFunctionBeginUser; 1391 ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 1392 ierr = PetscBagGetData(user->bag, &ctx);CHKERRQ(ierr); 1393 id = 3; 1394 ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "top wall velocity", label, 1, &id, VEL, 0, NULL, (void (*)(void)) exactFuncs[VEL], (void (*)(void)) exactFuncs_t[VEL], ctx, NULL);CHKERRQ(ierr); 1395 id = 1; 1396 ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "bottom wall velocity", label, 1, &id, VEL, 0, NULL, (void (*)(void)) exactFuncs[VEL], (void (*)(void)) exactFuncs_t[VEL], ctx, NULL);CHKERRQ(ierr); 1397 id = 2; 1398 ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "right wall velocity", label, 1, &id, VEL, 0, NULL, (void (*)(void)) exactFuncs[VEL], (void (*)(void)) exactFuncs_t[VEL], ctx, NULL);CHKERRQ(ierr); 1399 id = 4; 1400 ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "left wall velocity", label, 1, &id, VEL, 0, NULL, (void (*)(void)) exactFuncs[VEL], (void (*)(void)) exactFuncs_t[VEL], ctx, NULL);CHKERRQ(ierr); 1401 id = 3; 1402 ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "top wall temp", label, 1, &id, TEMP, 0, NULL, (void (*)(void)) exactFuncs[TEMP], (void (*)(void)) exactFuncs_t[TEMP], ctx, NULL);CHKERRQ(ierr); 1403 id = 1; 1404 ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "bottom wall temp", label, 1, &id, TEMP, 0, NULL, (void (*)(void)) exactFuncs[TEMP], (void (*)(void)) exactFuncs_t[TEMP], ctx, NULL);CHKERRQ(ierr); 1405 id = 2; 1406 ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "right wall temp", label, 1, &id, TEMP, 0, NULL, (void (*)(void)) exactFuncs[TEMP], (void (*)(void)) exactFuncs_t[TEMP], ctx, NULL);CHKERRQ(ierr); 1407 id = 4; 1408 ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "left wall temp", label, 1, &id, TEMP, 0, NULL, (void (*)(void)) exactFuncs[TEMP], (void (*)(void)) exactFuncs_t[TEMP], ctx, NULL);CHKERRQ(ierr); 1409 PetscFunctionReturn(0); 1410 } 1411 1412 static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 1413 { 1414 PetscSimplePointFunc exactFuncs[3]; 1415 PetscSimplePointFunc exactFuncs_t[3]; 1416 PetscDS ds; 1417 PetscWeakForm wf; 1418 DMLabel label; 1419 Parameter *ctx; 1420 PetscInt id, bd; 1421 PetscErrorCode ierr; 1422 1423 PetscFunctionBeginUser; 1424 ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr); 1425 ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 1426 ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr); 1427 1428 switch(user->modType) { 1429 case MOD_INCOMPRESSIBLE: 1430 ierr = PetscDSSetResidual(ds, VEL, f0_v, f1_v);CHKERRQ(ierr); 1431 ierr = PetscDSSetResidual(ds, PRES, f0_q, NULL);CHKERRQ(ierr); 1432 ierr = PetscDSSetResidual(ds, TEMP, f0_w, f1_w);CHKERRQ(ierr); 1433 1434 ierr = PetscDSSetJacobian(ds, VEL, VEL, g0_vu, g1_vu, NULL, g3_vu);CHKERRQ(ierr); 1435 ierr = PetscDSSetJacobian(ds, VEL, PRES, NULL, NULL, g2_vp, NULL);CHKERRQ(ierr); 1436 ierr = PetscDSSetJacobian(ds, PRES, VEL, NULL, g1_qu, NULL, NULL);CHKERRQ(ierr); 1437 ierr = PetscDSSetJacobian(ds, TEMP, VEL, g0_wu, NULL, NULL, NULL);CHKERRQ(ierr); 1438 ierr = PetscDSSetJacobian(ds, TEMP, TEMP, g0_wT, g1_wT, NULL, g3_wT);CHKERRQ(ierr); 1439 1440 switch(user->solType) { 1441 case SOL_QUADRATIC: 1442 ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, VEL, 0, 1, f0_quadratic_v, 0, NULL);CHKERRQ(ierr); 1443 ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, TEMP, 0, 1, f0_quadratic_w, 0, NULL);CHKERRQ(ierr); 1444 1445 exactFuncs[VEL] = quadratic_u; 1446 exactFuncs[PRES] = quadratic_p; 1447 exactFuncs[TEMP] = quadratic_T; 1448 exactFuncs_t[VEL] = quadratic_u_t; 1449 exactFuncs_t[PRES] = NULL; 1450 exactFuncs_t[TEMP] = quadratic_T_t; 1451 1452 ierr = UniformBoundaryConditions(dm, label, exactFuncs, exactFuncs_t, user);CHKERRQ(ierr); 1453 break; 1454 case SOL_CUBIC: 1455 ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, VEL, 0, 1, f0_cubic_v, 0, NULL);CHKERRQ(ierr); 1456 ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, TEMP, 0, 1, f0_cubic_w, 0, NULL);CHKERRQ(ierr); 1457 1458 exactFuncs[VEL] = cubic_u; 1459 exactFuncs[PRES] = cubic_p; 1460 exactFuncs[TEMP] = cubic_T; 1461 exactFuncs_t[VEL] = cubic_u_t; 1462 exactFuncs_t[PRES] = NULL; 1463 exactFuncs_t[TEMP] = cubic_T_t; 1464 1465 ierr = UniformBoundaryConditions(dm, label, exactFuncs, exactFuncs_t, user);CHKERRQ(ierr); 1466 break; 1467 case SOL_CUBIC_TRIG: 1468 ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, VEL, 0, 1, f0_cubic_trig_v, 0, NULL);CHKERRQ(ierr); 1469 ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, TEMP, 0, 1, f0_cubic_trig_w, 0, NULL);CHKERRQ(ierr); 1470 1471 exactFuncs[VEL] = cubic_trig_u; 1472 exactFuncs[PRES] = cubic_trig_p; 1473 exactFuncs[TEMP] = cubic_trig_T; 1474 exactFuncs_t[VEL] = cubic_trig_u_t; 1475 exactFuncs_t[PRES] = NULL; 1476 exactFuncs_t[TEMP] = cubic_trig_T_t; 1477 1478 ierr = UniformBoundaryConditions(dm, label, exactFuncs, exactFuncs_t, user);CHKERRQ(ierr); 1479 break; 1480 case SOL_TAYLOR_GREEN: 1481 ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, TEMP, 0, 1, f0_taylor_green_w, 0, NULL);CHKERRQ(ierr); 1482 1483 exactFuncs[VEL] = taylor_green_u; 1484 exactFuncs[PRES] = taylor_green_p; 1485 exactFuncs[TEMP] = taylor_green_T; 1486 exactFuncs_t[VEL] = taylor_green_u_t; 1487 exactFuncs_t[PRES] = taylor_green_p_t; 1488 exactFuncs_t[TEMP] = taylor_green_T_t; 1489 1490 ierr = UniformBoundaryConditions(dm, label, exactFuncs, exactFuncs_t, user);CHKERRQ(ierr); 1491 break; 1492 default: SETERRQ2(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType); 1493 } 1494 break; 1495 case MOD_CONDUCTING: 1496 ierr = PetscDSSetResidual(ds, VEL, f0_conduct_v, f1_conduct_v);CHKERRQ(ierr); 1497 ierr = PetscDSSetResidual(ds, PRES, f0_conduct_q, NULL);CHKERRQ(ierr); 1498 ierr = PetscDSSetResidual(ds, TEMP, f0_conduct_w, f1_conduct_w);CHKERRQ(ierr); 1499 1500 ierr = PetscDSSetJacobian(ds, VEL, VEL, g0_conduct_vu, g1_conduct_vu, NULL, g3_conduct_vu);CHKERRQ(ierr); 1501 ierr = PetscDSSetJacobian(ds, VEL, PRES, NULL, NULL, g2_conduct_vp, NULL);CHKERRQ(ierr); 1502 ierr = PetscDSSetJacobian(ds, VEL, TEMP, g0_conduct_vT, NULL, NULL, NULL);CHKERRQ(ierr); 1503 ierr = PetscDSSetJacobian(ds, PRES, VEL, g0_conduct_qu, g1_conduct_qu, NULL, NULL);CHKERRQ(ierr); 1504 ierr = PetscDSSetJacobian(ds, PRES, TEMP, g0_conduct_qT, g1_conduct_qT, NULL, NULL);CHKERRQ(ierr); 1505 ierr = PetscDSSetJacobian(ds, TEMP, VEL, g0_conduct_wu, NULL, NULL, NULL);CHKERRQ(ierr); 1506 ierr = PetscDSSetJacobian(ds, TEMP, TEMP, g0_conduct_wT, g1_conduct_wT, NULL, g3_conduct_wT);CHKERRQ(ierr); 1507 1508 switch(user->solType) { 1509 case SOL_QUADRATIC: 1510 ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, VEL, 0, 1, f0_conduct_quadratic_v, 0, NULL);CHKERRQ(ierr); 1511 ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, PRES, 0, 1, f0_conduct_quadratic_q, 0, NULL);CHKERRQ(ierr); 1512 ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, TEMP, 0, 1, f0_conduct_quadratic_w, 0, NULL);CHKERRQ(ierr); 1513 1514 exactFuncs[VEL] = quadratic_u; 1515 exactFuncs[PRES] = quadratic_p; 1516 exactFuncs[TEMP] = quadratic_T; 1517 exactFuncs_t[VEL] = quadratic_u_t; 1518 exactFuncs_t[PRES] = NULL; 1519 exactFuncs_t[TEMP] = quadratic_T_t; 1520 1521 ierr = UniformBoundaryConditions(dm, label, exactFuncs, exactFuncs_t, user);CHKERRQ(ierr); 1522 break; 1523 case SOL_PIPE: 1524 user->hasNullSpace = PETSC_FALSE; 1525 ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, VEL, 0, 1, f0_conduct_pipe_v, 0, NULL);CHKERRQ(ierr); 1526 ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, PRES, 0, 1, f0_conduct_pipe_q, 0, NULL);CHKERRQ(ierr); 1527 ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, TEMP, 0, 1, f0_conduct_pipe_w, 0, NULL);CHKERRQ(ierr); 1528 1529 exactFuncs[VEL] = pipe_u; 1530 exactFuncs[PRES] = pipe_p; 1531 exactFuncs[TEMP] = pipe_T; 1532 exactFuncs_t[VEL] = pipe_u_t; 1533 exactFuncs_t[PRES] = pipe_p_t; 1534 exactFuncs_t[TEMP] = pipe_T_t; 1535 1536 ierr = PetscBagGetData(user->bag, (void **) &ctx);CHKERRQ(ierr); 1537 id = 2; 1538 ierr = DMAddBoundary(dm, DM_BC_NATURAL, "right wall", label, 1, &id, 0, 0, NULL, NULL, NULL, ctx, &bd);CHKERRQ(ierr); 1539 ierr = PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 1540 ierr = PetscWeakFormSetIndexBdResidual(wf, label, id, VEL, 0, 0, f0_conduct_bd_pipe_v, 0, NULL);CHKERRQ(ierr); 1541 id = 4; 1542 ierr = DMAddBoundary(dm, DM_BC_NATURAL, "left wall", label, 1, &id, 0, 0, NULL, NULL, NULL, ctx, &bd);CHKERRQ(ierr); 1543 ierr = PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 1544 ierr = PetscWeakFormSetIndexBdResidual(wf, label, id, VEL, 0, 0, f0_conduct_bd_pipe_v, 0, NULL);CHKERRQ(ierr); 1545 id = 4; 1546 ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "left wall temperature", label, 1, &id, TEMP, 0, NULL, (void (*)(void)) exactFuncs[TEMP], (void (*)(void)) exactFuncs_t[TEMP], ctx, NULL);CHKERRQ(ierr); 1547 id = 3; 1548 ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "top wall velocity", label, 1, &id, VEL, 0, NULL, (void (*)(void)) exactFuncs[VEL], (void (*)(void)) exactFuncs_t[VEL], ctx, NULL);CHKERRQ(ierr); 1549 ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "top wall temperature", label, 1, &id, TEMP, 0, NULL, (void (*)(void)) exactFuncs[TEMP], (void (*)(void)) exactFuncs_t[TEMP], ctx, NULL);CHKERRQ(ierr); 1550 id = 1; 1551 ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "bottom wall velocity", label, 1, &id, VEL, 0, NULL, (void (*)(void)) exactFuncs[VEL], (void (*)(void)) exactFuncs_t[VEL], ctx, NULL);CHKERRQ(ierr); 1552 ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "bottom wall temperature", label, 1, &id, TEMP, 0, NULL, (void (*)(void)) exactFuncs[TEMP], (void (*)(void)) exactFuncs_t[TEMP], ctx, NULL);CHKERRQ(ierr); 1553 break; 1554 case SOL_PIPE_WIGGLY: 1555 user->hasNullSpace = PETSC_FALSE; 1556 ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, VEL, 0, 1, f0_conduct_pipe_wiggly_v, 0, NULL);CHKERRQ(ierr); 1557 ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, PRES, 0, 1, f0_conduct_pipe_wiggly_q, 0, NULL);CHKERRQ(ierr); 1558 ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, TEMP, 0, 1, f0_conduct_pipe_wiggly_w, 0, NULL);CHKERRQ(ierr); 1559 1560 exactFuncs[VEL] = pipe_wiggly_u; 1561 exactFuncs[PRES] = pipe_wiggly_p; 1562 exactFuncs[TEMP] = pipe_wiggly_T; 1563 exactFuncs_t[VEL] = pipe_wiggly_u_t; 1564 exactFuncs_t[PRES] = pipe_wiggly_p_t; 1565 exactFuncs_t[TEMP] = pipe_wiggly_T_t; 1566 1567 ierr = PetscBagGetData(user->bag, (void **) &ctx);CHKERRQ(ierr); 1568 id = 2; 1569 ierr = DMAddBoundary(dm, DM_BC_NATURAL, "right wall", label, 1, &id, 0, 0, NULL, NULL, NULL, ctx, &bd);CHKERRQ(ierr); 1570 ierr = PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 1571 ierr = PetscWeakFormSetIndexBdResidual(wf, label, id, VEL, 0, 0, f0_conduct_bd_pipe_wiggly_v, 0, NULL);CHKERRQ(ierr); 1572 id = 4; 1573 ierr = DMAddBoundary(dm, DM_BC_NATURAL, "left wall", label, 1, &id, 0, 0, NULL, NULL, NULL, ctx, &bd);CHKERRQ(ierr); 1574 ierr = PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 1575 ierr = PetscWeakFormSetIndexBdResidual(wf, label, id, VEL, 0, 0, f0_conduct_bd_pipe_wiggly_v, 0, NULL);CHKERRQ(ierr); 1576 id = 4; 1577 ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "left wall temperature", label, 1, &id, TEMP, 0, NULL, (void (*)(void)) exactFuncs[TEMP], (void (*)(void)) exactFuncs_t[TEMP], ctx, NULL);CHKERRQ(ierr); 1578 id = 3; 1579 ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "top wall velocity", label, 1, &id, VEL, 0, NULL, (void (*)(void)) exactFuncs[VEL], (void (*)(void)) exactFuncs_t[VEL], ctx, NULL);CHKERRQ(ierr); 1580 ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "top wall temperature", label, 1, &id, TEMP, 0, NULL, (void (*)(void)) exactFuncs[TEMP], (void (*)(void)) exactFuncs_t[TEMP], ctx, NULL);CHKERRQ(ierr); 1581 id = 1; 1582 ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "bottom wall velocity", label, 1, &id, VEL, 0, NULL, (void (*)(void)) exactFuncs[VEL], (void (*)(void)) exactFuncs_t[VEL], ctx, NULL);CHKERRQ(ierr); 1583 ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "bottom wall temperature", label, 1, &id, TEMP, 0, NULL, (void (*)(void)) exactFuncs[TEMP], (void (*)(void)) exactFuncs_t[TEMP], ctx, NULL);CHKERRQ(ierr); 1584 break; 1585 default: SETERRQ2(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType); 1586 } 1587 break; 1588 default: SETERRQ2(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Unsupported model type: %s (%D)", solTypes[PetscMin(user->modType, NUM_MOD_TYPES)], user->modType); 1589 } 1590 /* Setup constants */ 1591 { 1592 Parameter *param; 1593 PetscScalar constants[13]; 1594 1595 ierr = PetscBagGetData(user->bag, (void **) ¶m);CHKERRQ(ierr); 1596 1597 constants[STROUHAL] = param->Strouhal; 1598 constants[FROUDE] = param->Froude; 1599 constants[REYNOLDS] = param->Reynolds; 1600 constants[PECLET] = param->Peclet; 1601 constants[P_TH] = param->p_th; 1602 constants[MU] = param->mu; 1603 constants[NU] = param->nu; 1604 constants[C_P] = param->c_p; 1605 constants[K] = param->k; 1606 constants[ALPHA] = param->alpha; 1607 constants[T_IN] = param->T_in; 1608 constants[G_DIR] = param->g_dir; 1609 constants[EPSILON] = param->epsilon; 1610 ierr = PetscDSSetConstants(ds, 13, constants);CHKERRQ(ierr); 1611 } 1612 1613 ierr = PetscBagGetData(user->bag, (void **) &ctx);CHKERRQ(ierr); 1614 ierr = PetscDSSetExactSolution(ds, VEL, exactFuncs[VEL], ctx);CHKERRQ(ierr); 1615 ierr = PetscDSSetExactSolution(ds, PRES, exactFuncs[PRES], ctx);CHKERRQ(ierr); 1616 ierr = PetscDSSetExactSolution(ds, TEMP, exactFuncs[TEMP], ctx);CHKERRQ(ierr); 1617 ierr = PetscDSSetExactSolutionTimeDerivative(ds, VEL, exactFuncs_t[VEL], ctx);CHKERRQ(ierr); 1618 ierr = PetscDSSetExactSolutionTimeDerivative(ds, PRES, exactFuncs_t[PRES], ctx);CHKERRQ(ierr); 1619 ierr = PetscDSSetExactSolutionTimeDerivative(ds, TEMP, exactFuncs_t[TEMP], ctx);CHKERRQ(ierr); 1620 PetscFunctionReturn(0); 1621 } 1622 1623 static PetscErrorCode CreateCellDM(DM dm, AppCtx *user) 1624 { 1625 PetscFE fe, fediv; 1626 DMPolytopeType ct; 1627 PetscInt dim, cStart; 1628 PetscBool simplex; 1629 PetscErrorCode ierr; 1630 1631 PetscFunctionBeginUser; 1632 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1633 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr); 1634 ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr); 1635 simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 1636 1637 ierr = DMGetField(dm, VEL, NULL, (PetscObject *) &fe);CHKERRQ(ierr); 1638 ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "div_", PETSC_DEFAULT, &fediv);CHKERRQ(ierr); 1639 ierr = PetscFECopyQuadrature(fe, fediv);CHKERRQ(ierr); 1640 ierr = PetscObjectSetName((PetscObject) fediv, "divergence");CHKERRQ(ierr); 1641 1642 ierr = DMDestroy(&user->dmCell);CHKERRQ(ierr); 1643 ierr = DMClone(dm, &user->dmCell);CHKERRQ(ierr); 1644 ierr = DMSetField(user->dmCell, 0, NULL, (PetscObject) fediv);CHKERRQ(ierr); 1645 ierr = DMCreateDS(user->dmCell);CHKERRQ(ierr); 1646 ierr = PetscFEDestroy(&fediv);CHKERRQ(ierr); 1647 PetscFunctionReturn(0); 1648 } 1649 1650 static PetscErrorCode GetCellDM(DM dm, AppCtx *user, DM *dmCell) 1651 { 1652 PetscInt cStart, cEnd, cellStart = -1, cellEnd = -1; 1653 PetscErrorCode ierr; 1654 1655 PetscFunctionBeginUser; 1656 ierr = DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1657 if (user->dmCell) {ierr = DMPlexGetSimplexOrBoxCells(user->dmCell, 0, &cellStart, &cellEnd);CHKERRQ(ierr);} 1658 if (cStart != cellStart || cEnd != cellEnd) {ierr = CreateCellDM(dm, user);CHKERRQ(ierr);} 1659 *dmCell = user->dmCell; 1660 PetscFunctionReturn(0); 1661 } 1662 1663 static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 1664 { 1665 DM cdm = dm; 1666 PetscFE fe[3]; 1667 Parameter *param; 1668 DMPolytopeType ct; 1669 PetscInt dim, cStart; 1670 PetscBool simplex; 1671 PetscErrorCode ierr; 1672 1673 PetscFunctionBeginUser; 1674 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1675 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr); 1676 ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr); 1677 simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 1678 /* Create finite element */ 1679 ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr); 1680 ierr = PetscObjectSetName((PetscObject) fe[0], "velocity");CHKERRQ(ierr); 1681 1682 ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr); 1683 ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr); 1684 ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr); 1685 1686 ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "temp_", PETSC_DEFAULT, &fe[2]);CHKERRQ(ierr); 1687 ierr = PetscFECopyQuadrature(fe[0], fe[2]);CHKERRQ(ierr); 1688 ierr = PetscObjectSetName((PetscObject) fe[2], "temperature");CHKERRQ(ierr); 1689 1690 /* Set discretization and boundary conditions for each mesh */ 1691 ierr = DMSetField(dm, VEL, NULL, (PetscObject) fe[VEL]);CHKERRQ(ierr); 1692 ierr = DMSetField(dm, PRES, NULL, (PetscObject) fe[PRES]);CHKERRQ(ierr); 1693 ierr = DMSetField(dm, TEMP, NULL, (PetscObject) fe[TEMP]);CHKERRQ(ierr); 1694 ierr = DMCreateDS(dm);CHKERRQ(ierr); 1695 ierr = SetupProblem(dm, user);CHKERRQ(ierr); 1696 ierr = PetscBagGetData(user->bag, (void **) ¶m);CHKERRQ(ierr); 1697 while (cdm) { 1698 ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr); 1699 ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 1700 } 1701 ierr = PetscFEDestroy(&fe[VEL]);CHKERRQ(ierr); 1702 ierr = PetscFEDestroy(&fe[PRES]);CHKERRQ(ierr); 1703 ierr = PetscFEDestroy(&fe[TEMP]);CHKERRQ(ierr); 1704 1705 if (user->hasNullSpace) { 1706 PetscObject pressure; 1707 MatNullSpace nullspacePres; 1708 1709 ierr = DMGetField(dm, PRES, NULL, &pressure);CHKERRQ(ierr); 1710 ierr = MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres);CHKERRQ(ierr); 1711 ierr = PetscObjectCompose(pressure, "nullspace", (PetscObject) nullspacePres);CHKERRQ(ierr); 1712 ierr = MatNullSpaceDestroy(&nullspacePres);CHKERRQ(ierr); 1713 } 1714 PetscFunctionReturn(0); 1715 } 1716 1717 static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace) 1718 { 1719 Vec vec; 1720 PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero}; 1721 PetscErrorCode ierr; 1722 1723 PetscFunctionBeginUser; 1724 if (ofield != PRES) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Nullspace must be for pressure field at index %D, not %D", PRES, ofield); 1725 funcs[nfield] = constant; 1726 ierr = DMCreateGlobalVector(dm, &vec);CHKERRQ(ierr); 1727 ierr = DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec);CHKERRQ(ierr); 1728 ierr = VecNormalize(vec, NULL);CHKERRQ(ierr); 1729 ierr = PetscObjectSetName((PetscObject) vec, "Pressure Null Space");CHKERRQ(ierr); 1730 ierr = VecViewFromOptions(vec, NULL, "-pressure_nullspace_view");CHKERRQ(ierr); 1731 ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_FALSE, 1, &vec, nullSpace);CHKERRQ(ierr); 1732 ierr = VecDestroy(&vec);CHKERRQ(ierr); 1733 PetscFunctionReturn(0); 1734 } 1735 1736 static PetscErrorCode RemoveDiscretePressureNullspace_Private(TS ts, Vec u) 1737 { 1738 DM dm; 1739 AppCtx *user; 1740 MatNullSpace nullsp; 1741 PetscErrorCode ierr; 1742 1743 PetscFunctionBegin; 1744 ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 1745 ierr = DMGetApplicationContext(dm, &user);CHKERRQ(ierr); 1746 if (!user->hasNullSpace) PetscFunctionReturn(0); 1747 ierr = CreatePressureNullSpace(dm, 1, 1, &nullsp);CHKERRQ(ierr); 1748 ierr = MatNullSpaceRemove(nullsp, u);CHKERRQ(ierr); 1749 ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr); 1750 PetscFunctionReturn(0); 1751 } 1752 1753 /* Make the discrete pressure discretely divergence free */ 1754 static PetscErrorCode RemoveDiscretePressureNullspace(TS ts) 1755 { 1756 Vec u; 1757 PetscErrorCode ierr; 1758 1759 PetscFunctionBegin; 1760 ierr = TSGetSolution(ts, &u);CHKERRQ(ierr); 1761 ierr = RemoveDiscretePressureNullspace_Private(ts, u);CHKERRQ(ierr); 1762 PetscFunctionReturn(0); 1763 } 1764 1765 static void divergence(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1766 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1767 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1768 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar divu[]) 1769 { 1770 PetscInt d; 1771 1772 divu[0] = 0.; 1773 for (d = 0; d < dim; ++d) divu[0] += u_x[d*dim+d]; 1774 } 1775 1776 static PetscErrorCode SetInitialConditions(TS ts, Vec u) 1777 { 1778 AppCtx *user; 1779 DM dm; 1780 PetscReal t; 1781 PetscErrorCode ierr; 1782 1783 PetscFunctionBegin; 1784 ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 1785 ierr = TSGetTime(ts, &t);CHKERRQ(ierr); 1786 ierr = DMComputeExactSolution(dm, t, u, NULL);CHKERRQ(ierr); 1787 ierr = DMGetApplicationContext(dm, &user);CHKERRQ(ierr); 1788 ierr = RemoveDiscretePressureNullspace_Private(ts, u);CHKERRQ(ierr); 1789 PetscFunctionReturn(0); 1790 } 1791 1792 static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx) 1793 { 1794 PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 1795 void *ctxs[3]; 1796 PetscPointFunc diagnostics[1] = {divergence}; 1797 DM dm, dmCell = NULL; 1798 PetscDS ds; 1799 Vec v, divu; 1800 PetscReal ferrors[3], massFlux; 1801 PetscInt f; 1802 PetscErrorCode ierr; 1803 1804 PetscFunctionBeginUser; 1805 ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 1806 ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 1807 1808 for (f = 0; f < 3; ++f) {ierr = PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f]);CHKERRQ(ierr);} 1809 ierr = DMComputeL2FieldDiff(dm, crtime, exactFuncs, ctxs, u, ferrors);CHKERRQ(ierr); 1810 ierr = GetCellDM(dm, (AppCtx *) ctx, &dmCell);CHKERRQ(ierr); 1811 ierr = DMGetGlobalVector(dmCell, &divu);CHKERRQ(ierr); 1812 ierr = DMProjectField(dmCell, crtime, u, diagnostics, INSERT_VALUES, divu);CHKERRQ(ierr); 1813 ierr = VecViewFromOptions(divu, NULL, "-divu_vec_view");CHKERRQ(ierr); 1814 ierr = VecNorm(divu, NORM_2, &massFlux);CHKERRQ(ierr); 1815 ierr = PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: [%2.3g, %2.3g, %2.3g] ||div u||: %2.3g\n", (int) step, (double) crtime, (double) ferrors[0], (double) ferrors[1], (double) ferrors[2], (double) massFlux);CHKERRQ(ierr); 1816 1817 ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr); 1818 1819 ierr = DMGetGlobalVector(dm, &v);CHKERRQ(ierr); 1820 ierr = DMProjectFunction(dm, crtime, exactFuncs, ctxs, INSERT_ALL_VALUES, v);CHKERRQ(ierr); 1821 ierr = PetscObjectSetName((PetscObject) v, "Exact Solution");CHKERRQ(ierr); 1822 ierr = VecViewFromOptions(v, NULL, "-exact_vec_view");CHKERRQ(ierr); 1823 ierr = DMRestoreGlobalVector(dm, &v);CHKERRQ(ierr); 1824 1825 ierr = VecViewFromOptions(divu, NULL, "-div_vec_view");CHKERRQ(ierr); 1826 ierr = DMRestoreGlobalVector(dmCell, &divu);CHKERRQ(ierr); 1827 1828 PetscFunctionReturn(0); 1829 } 1830 1831 int main(int argc, char **argv) 1832 { 1833 DM dm; /* problem definition */ 1834 TS ts; /* timestepper */ 1835 Vec u; /* solution */ 1836 AppCtx user; /* user-defined work context */ 1837 PetscReal t; 1838 PetscErrorCode ierr; 1839 1840 ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 1841 ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 1842 ierr = PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag);CHKERRQ(ierr); 1843 ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 1844 ierr = SetupParameters(dm, &user);CHKERRQ(ierr); 1845 ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr); 1846 ierr = TSSetDM(ts, dm);CHKERRQ(ierr); 1847 ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr); 1848 /* Setup problem */ 1849 ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr); 1850 ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr); 1851 1852 ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 1853 ierr = PetscObjectSetName((PetscObject) u, "Numerical Solution");CHKERRQ(ierr); 1854 if (user.hasNullSpace) {ierr = DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace);CHKERRQ(ierr);} 1855 1856 ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user);CHKERRQ(ierr); 1857 ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user);CHKERRQ(ierr); 1858 ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user);CHKERRQ(ierr); 1859 ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 1860 ierr = TSSetPreStep(ts, RemoveDiscretePressureNullspace);CHKERRQ(ierr); 1861 ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 1862 1863 ierr = TSSetComputeInitialCondition(ts, SetInitialConditions);CHKERRQ(ierr); /* Must come after SetFromOptions() */ 1864 ierr = SetInitialConditions(ts, u);CHKERRQ(ierr); 1865 ierr = TSGetTime(ts, &t);CHKERRQ(ierr); 1866 ierr = DMSetOutputSequenceNumber(dm, 0, t);CHKERRQ(ierr); 1867 ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr); 1868 ierr = TSMonitorSet(ts, MonitorError, &user, NULL);CHKERRQ(ierr);CHKERRQ(ierr); 1869 1870 ierr = TSSolve(ts, u);CHKERRQ(ierr); 1871 ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr); 1872 1873 ierr = VecDestroy(&u);CHKERRQ(ierr); 1874 ierr = DMDestroy(&user.dmCell);CHKERRQ(ierr); 1875 ierr = DMDestroy(&dm);CHKERRQ(ierr); 1876 ierr = TSDestroy(&ts);CHKERRQ(ierr); 1877 ierr = PetscBagDestroy(&user.bag);CHKERRQ(ierr); 1878 ierr = PetscFinalize(); 1879 return ierr; 1880 } 1881 1882 /*TEST 1883 1884 testset: 1885 requires: triangle !single 1886 args: -dm_plex_separate_marker \ 1887 -div_petscdualspace_lagrange_use_moments -div_petscdualspace_lagrange_moment_order 3 \ 1888 -snes_error_if_not_converged -snes_convergence_test correct_pressure \ 1889 -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 1890 -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 \ 1891 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 1892 -fieldsplit_0_pc_type lu \ 1893 -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 1894 1895 test: 1896 suffix: 2d_tri_p2_p1_p1 1897 args: -sol_type quadratic \ 1898 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 1899 -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1 1900 1901 test: 1902 # Using -dm_refine 5 -convest_num_refine 2 gives L_2 convergence rate: [0.89, 0.011, 1.0] 1903 suffix: 2d_tri_p2_p1_p1_tconv 1904 args: -sol_type cubic_trig \ 1905 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 1906 -ts_max_steps 4 -ts_dt 0.1 -ts_convergence_estimate -convest_num_refine 1 1907 1908 test: 1909 # Using -dm_refine 3 -convest_num_refine 3 gives L_2 convergence rate: [3.0, 2.5, 1.9] 1910 suffix: 2d_tri_p2_p1_p1_sconv 1911 args: -sol_type cubic \ 1912 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 1913 -ts_max_steps 1 -ts_dt 1e-4 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 1914 1915 test: 1916 suffix: 2d_tri_p3_p2_p2 1917 args: -sol_type cubic \ 1918 -vel_petscspace_degree 3 -pres_petscspace_degree 2 -temp_petscspace_degree 2 \ 1919 -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1 1920 1921 test: 1922 # Using -dm_refine 3 -convest_num_refine 3 gives L_2 convergence rate: [3.0, 2.1, 3.1] 1923 suffix: 2d_tri_p2_p1_p1_tg_sconv 1924 args: -sol_type taylor_green \ 1925 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 1926 -ts_max_steps 1 -ts_dt 1e-8 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 1927 1928 test: 1929 # Using -dm_refine 3 -convest_num_refine 2 gives L_2 convergence rate: [1.2, 1.5, 1.2] 1930 suffix: 2d_tri_p2_p1_p1_tg_tconv 1931 args: -sol_type taylor_green \ 1932 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 1933 -ts_max_steps 4 -ts_dt 0.1 -ts_convergence_estimate -convest_num_refine 1 1934 1935 testset: 1936 requires: triangle !single 1937 args: -dm_plex_separate_marker -mod_type conducting \ 1938 -div_petscdualspace_lagrange_use_moments -div_petscdualspace_lagrange_moment_order 3 \ 1939 -snes_error_if_not_converged -snes_max_linear_solve_fail 5 \ 1940 -ksp_type fgmres -ksp_max_it 2 -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 1941 -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 \ 1942 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 1943 -fieldsplit_0_pc_type lu \ 1944 -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 1945 1946 test: 1947 # At this resolution, the rhs is inconsistent on some Newton steps 1948 suffix: 2d_tri_p2_p1_p1_conduct 1949 args: -sol_type quadratic \ 1950 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 1951 -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1 \ 1952 -pc_fieldsplit_schur_precondition full \ 1953 -fieldsplit_pressure_ksp_max_it 2 -fieldsplit_pressure_pc_type svd 1954 1955 test: 1956 suffix: 2d_tri_p2_p1_p2_conduct_pipe 1957 args: -sol_type pipe \ 1958 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 2 \ 1959 -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1 1960 1961 test: 1962 suffix: 2d_tri_p2_p1_p2_conduct_pipe_wiggly_sconv 1963 args: -sol_type pipe_wiggly -Fr 1e10 \ 1964 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 2 \ 1965 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 1966 -ts_max_steps 1 -ts_dt 1e10 \ 1967 -ksp_atol 1e-12 -ksp_max_it 300 \ 1968 -fieldsplit_pressure_ksp_atol 1e-14 1969 1970 TEST*/ 1971