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 in a rectangular\n\ 3 domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n"; 4 5 /*F 6 This 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 For visualization, use 18 19 -dm_view hdf5:$PWD/sol.h5 -sol_vec_view hdf5:$PWD/sol.h5::append -exact_vec_view hdf5:$PWD/sol.h5::append 20 F*/ 21 22 #include <petscdmplex.h> 23 #include <petscsnes.h> 24 #include <petscts.h> 25 #include <petscds.h> 26 #include <petscbag.h> 27 28 typedef enum {SOL_QUADRATIC, SOL_CUBIC, SOL_CUBIC_TRIG, NUM_SOL_TYPES} SolType; 29 const char *solTypes[NUM_SOL_TYPES+1] = {"quadratic", "cubic", "cubic_trig", "unknown"}; 30 31 typedef struct { 32 PetscReal nu; /* Kinematic viscosity */ 33 PetscReal alpha; /* Thermal diffusivity */ 34 PetscReal T_in; /* Inlet temperature*/ 35 } Parameter; 36 37 typedef struct { 38 /* Problem definition */ 39 PetscBag bag; /* Holds problem parameters */ 40 SolType solType; /* MMS solution type */ 41 } AppCtx; 42 43 static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 44 { 45 PetscInt d; 46 for (d = 0; d < Nc; ++d) u[d] = 0.0; 47 return 0; 48 } 49 50 static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 51 { 52 PetscInt d; 53 for (d = 0; d < Nc; ++d) u[d] = 1.0; 54 return 0; 55 } 56 57 /* 58 CASE: quadratic 59 In 2D we use exact solution: 60 61 u = t + x^2 + y^2 62 v = t + 2x^2 - 2xy 63 p = x + y - 1 64 T = t + x + y 65 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> 66 Q = 1 + 2t + 3x^2 - 2xy + y^2 67 68 so that 69 70 \nabla \cdot u = 2x - 2x = 0 71 72 f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p 73 = <1, 1> + <t + x^2 + y^2, t + 2x^2 - 2xy> . <<2x, 4x - 2y>, <2y, -2x>> - \nu <4, 4> + <1, 1> 74 = <t (2x + 2y) + 2x^3 + 4x^2y - 2xy^2, t (2x - 2y) + 2x^2y + 4xy^2 - 2y^3> + <-4 \nu + 2, -4\nu + 2> 75 = <t (2x + 2y) + 2x^3 + 4x^2y - 2xy^2 - 4\nu + 2, t (2x - 2y) + 4xy^2 + 2x^2y - 2y^3 - 4\nu + 2> 76 77 Q = dT/dt + u \cdot \nabla T - \alpha \Delta T 78 = 1 + <t + x^2 + y^2, t + 2x^2 - 2xy> . <1, 1> - \alpha 0 79 = 1 + 2t + 3x^2 - 2xy + y^2 80 */ 81 82 static PetscErrorCode quadratic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 83 { 84 u[0] = time + X[0]*X[0] + X[1]*X[1]; 85 u[1] = time + 2.0*X[0]*X[0] - 2.0*X[0]*X[1]; 86 return 0; 87 } 88 static PetscErrorCode quadratic_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 89 { 90 u[0] = 1.0; 91 u[1] = 1.0; 92 return 0; 93 } 94 95 static PetscErrorCode quadratic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 96 { 97 p[0] = X[0] + X[1] - 1.0; 98 return 0; 99 } 100 101 static PetscErrorCode quadratic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 102 { 103 T[0] = time + X[0] + X[1]; 104 return 0; 105 } 106 static PetscErrorCode quadratic_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 107 { 108 T[0] = 1.0; 109 return 0; 110 } 111 112 /* f0_v = du/dt - f */ 113 static void f0_quadratic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 114 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 115 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 116 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 117 { 118 const PetscReal nu = PetscRealPart(constants[0]); 119 PetscInt Nc = dim; 120 PetscInt c, d; 121 122 for (d = 0; d<dim; ++d) f0[d] = u_t[uOff[0]+d]; 123 124 for (c = 0; c < Nc; ++c) { 125 for (d = 0; d < dim; ++d) f0[c] += u[d]*u_x[c*dim+d]; 126 } 127 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); 128 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); 129 } 130 131 /* f0_w = dT/dt + u.grad(T) - Q */ 132 static void f0_quadratic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 133 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 134 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 135 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 136 { 137 PetscInt d; 138 f0[0] = 0; 139 for (d = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d]; 140 f0[0] += u_t[uOff[2]] - (2*t + 1 + 3*X[0]*X[0] - 2*X[0]*X[1] + X[1]*X[1]); 141 } 142 143 /* 144 CASE: cubic 145 In 2D we use exact solution: 146 147 u = t + x^3 + y^3 148 v = t + 2x^3 - 3x^2y 149 p = 3/2 x^2 + 3/2 y^2 - 1 150 T = t + 1/2 x^2 + 1/2 y^2 151 f = < t(3x^2 + 3y^2) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y) + 3x + 1, 152 t(3x^2 - 6xy) + 6x^2y^3 + 3x^4y - 6xy^4 - \nu(12x - 6y) + 3y + 1> 153 Q = x^4 + xy^3 + 2x^3y - 3x^2y^2 + xt + yt - 2\alpha + 1 154 155 so that 156 157 \nabla \cdot u = 3x^2 - 3x^2 = 0 158 159 du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p - f 160 = <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 161 162 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 163 */ 164 static PetscErrorCode cubic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 165 { 166 u[0] = time + X[0]*X[0]*X[0] + X[1]*X[1]*X[1]; 167 u[1] = time + 2.0*X[0]*X[0]*X[0] - 3.0*X[0]*X[0]*X[1]; 168 return 0; 169 } 170 static PetscErrorCode cubic_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 171 { 172 u[0] = 1.0; 173 u[1] = 1.0; 174 return 0; 175 } 176 177 static PetscErrorCode cubic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 178 { 179 p[0] = 3.0*X[0]*X[0]/2.0 + 3.0*X[1]*X[1]/2.0 - 1.0; 180 return 0; 181 } 182 183 static PetscErrorCode cubic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 184 { 185 T[0] = time + X[0]*X[0]/2.0 + X[1]*X[1]/2.0; 186 return 0; 187 } 188 static PetscErrorCode cubic_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 189 { 190 T[0] = 1.0; 191 return 0; 192 } 193 194 195 static void f0_cubic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 196 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 197 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 198 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 199 { 200 PetscInt c, d; 201 PetscInt Nc = dim; 202 const PetscReal nu = PetscRealPart(constants[0]); 203 204 for (d=0; d<dim; ++d) f0[d] = u_t[uOff[0]+d]; 205 206 for (c=0; c<Nc; ++c) { 207 for (d=0; d<dim; ++d) f0[c] += u[d]*u_x[c*dim+d]; 208 } 209 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); 210 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); 211 } 212 213 static void f0_cubic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 214 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 215 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 216 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 217 { 218 PetscInt d; 219 const PetscReal alpha = PetscRealPart(constants[1]); 220 221 for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d]; 222 f0[0] += u_t[uOff[2]] - (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); 223 } 224 225 /* 226 CASE: cubic-trigonometric 227 In 2D we use exact solution: 228 229 u = beta cos t + x^3 + y^3 230 v = beta sin t + 2x^3 - 3x^2y 231 p = 3/2 x^2 + 3/2 y^2 - 1 232 T = 20 cos t + 1/2 x^2 + 1/2 y^2 233 f = < beta cos t 3x^2 + beta sin t (3y^2 - 1) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y) + 3x, 234 beta cos t (6x^2 - 6xy) - beta sin t (3x^2) + 3x^4y + 6x^2y^3 - 6xy^4 - \nu(12x - 6y) + 3y> 235 Q = beta cos t x + beta sin t (y - 1) + x^4 + 2x^3y - 3x^2y^2 + xy^3 - 2\alpha 236 237 so that 238 239 \nabla \cdot u = 3x^2 - 3x^2 = 0 240 241 f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p 242 = <-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> 243 = <-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> 244 = <cos t (3x^2) + sin t (3y^2 - 1) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu (6x + 6y) + 3x, 245 cos t (6x^2 - 6xy) - sin t (3x^2) + 3x^4y + 6x^2y^3 - 6xy^4 - \nu (12x - 6y) + 3y> 246 247 Q = dT/dt + u \cdot \nabla T - \alpha \Delta T 248 = -sin t + <cos t + x^3 + y^3, sin t + 2x^3 - 3x^2y> . <x, y> - 2 \alpha 249 = -sin t + cos t (x) + x^4 + xy^3 + sin t (y) + 2x^3y - 3x^2y^2 - 2 \alpha 250 = cos t x + sin t (y - 1) + (x^4 + 2x^3y - 3x^2y^2 + xy^3 - 2 \alpha) 251 */ 252 static PetscErrorCode cubic_trig_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 253 { 254 u[0] = 100.*PetscCosReal(time) + X[0]*X[0]*X[0] + X[1]*X[1]*X[1]; 255 u[1] = 100.*PetscSinReal(time) + 2.0*X[0]*X[0]*X[0] - 3.0*X[0]*X[0]*X[1]; 256 return 0; 257 } 258 static PetscErrorCode cubic_trig_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 259 { 260 u[0] = -100.*PetscSinReal(time); 261 u[1] = 100.*PetscCosReal(time); 262 return 0; 263 } 264 265 static PetscErrorCode cubic_trig_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 266 { 267 p[0] = 3.0*X[0]*X[0]/2.0 + 3.0*X[1]*X[1]/2.0 - 1.0; 268 return 0; 269 } 270 271 static PetscErrorCode cubic_trig_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 272 { 273 T[0] = 100.*PetscCosReal(time) + X[0]*X[0]/2.0 + X[1]*X[1]/2.0; 274 return 0; 275 } 276 static PetscErrorCode cubic_trig_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 277 { 278 T[0] = -100.*PetscSinReal(time); 279 return 0; 280 } 281 282 static void f0_cubic_trig_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 283 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 284 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 285 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 286 { 287 const PetscReal nu = PetscRealPart(constants[0]); 288 PetscInt Nc = dim; 289 PetscInt c, d; 290 291 for (d = 0; d < dim; ++d) f0[d] = u_t[uOff[0]+d]; 292 293 for (c = 0; c < Nc; ++c) { 294 for (d = 0; d < dim; ++d) f0[c] += u[d]*u_x[c*dim+d]; 295 } 296 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]; 297 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]; 298 } 299 300 static void f0_cubic_trig_w(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 alpha = PetscRealPart(constants[1]); 306 PetscInt d; 307 308 for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d]; 309 f0[0] += u_t[uOff[2]] - (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); 310 } 311 312 static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 313 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 314 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 315 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 316 { 317 PetscInt d; 318 for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d]; 319 } 320 321 /*f1_v = \nu[grad(u) + grad(u)^T] - pI */ 322 static void f1_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 323 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 324 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 325 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 326 { 327 const PetscReal nu = PetscRealPart(constants[0]); 328 const PetscInt Nc = dim; 329 PetscInt c, d; 330 331 for (c = 0; c < Nc; ++c) { 332 for (d = 0; d < dim; ++d) { 333 f1[c*dim+d] = nu*(u_x[c*dim+d] + u_x[d*dim+c]); 334 //f1[c*dim+d] = nu*u_x[c*dim+d]; 335 } 336 f1[c*dim+c] -= u[uOff[1]]; 337 } 338 } 339 340 static void f1_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 341 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 342 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 343 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 344 { 345 const PetscReal alpha = PetscRealPart(constants[1]); 346 PetscInt d; 347 for (d = 0; d < dim; ++d) f1[d] = alpha*u_x[uOff_x[2]+d]; 348 } 349 350 /*Jacobians*/ 351 static void g1_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 352 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 353 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 354 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 355 { 356 PetscInt d; 357 for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; 358 } 359 360 static void g0_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 361 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 362 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 363 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 364 { 365 PetscInt c, d; 366 const PetscInt Nc = dim; 367 368 for (d = 0; d < dim; ++d) g0[d*dim+d] = u_tShift; 369 370 for (c = 0; c < Nc; ++c) { 371 for (d = 0; d < dim; ++d) { 372 g0[c*Nc+d] += u_x[ c*Nc+d]; 373 } 374 } 375 } 376 377 static void g1_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 378 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 379 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 380 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 381 { 382 PetscInt NcI = dim; 383 PetscInt NcJ = dim; 384 PetscInt c, d, e; 385 386 for (c = 0; c < NcI; ++c) { 387 for (d = 0; d < NcJ; ++d) { 388 for (e = 0; e < dim; ++e) { 389 if (c == d) { 390 g1[(c*NcJ+d)*dim+e] += u[e]; 391 } 392 } 393 } 394 } 395 } 396 397 398 static void g2_vp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 399 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 400 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 401 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 402 { 403 PetscInt d; 404 for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; 405 } 406 407 static void g3_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 408 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 409 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 410 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 411 { 412 const PetscReal nu = PetscRealPart(constants[0]); 413 const PetscInt Nc = dim; 414 PetscInt c, d; 415 416 for (c = 0; c < Nc; ++c) { 417 for (d = 0; d < dim; ++d) { 418 g3[((c*Nc+c)*dim+d)*dim+d] += nu; // gradU 419 g3[((c*Nc+d)*dim+d)*dim+c] += nu; // gradU transpose 420 } 421 } 422 } 423 424 static void g0_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 425 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 426 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 427 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 428 { 429 PetscInt d; 430 for (d = 0; d < dim; ++d) g0[d] = u_tShift; 431 } 432 433 static void g0_wu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 434 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 435 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 436 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 437 { 438 PetscInt d; 439 for (d = 0; d < dim; ++d) g0[d] = u_x[uOff_x[2]+d]; 440 } 441 442 static void g1_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 443 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 444 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 445 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 446 { 447 PetscInt d; 448 for (d = 0; d < dim; ++d) g1[d] = u[uOff[0]+d]; 449 } 450 451 static void g3_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 452 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 453 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 454 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 455 { 456 const PetscReal alpha = PetscRealPart(constants[1]); 457 PetscInt d; 458 459 for (d = 0; d < dim; ++d) g3[d*dim+d] = alpha; 460 } 461 462 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 463 { 464 PetscInt sol; 465 PetscErrorCode ierr; 466 467 468 PetscFunctionBeginUser; 469 options->solType = SOL_QUADRATIC; 470 471 ierr = PetscOptionsBegin(comm, "", "Low Mach flow Problem Options", "DMPLEX");CHKERRQ(ierr); 472 sol = options->solType; 473 ierr = PetscOptionsEList("-sol_type", "The solution type", "ex62.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL);CHKERRQ(ierr); 474 options->solType = (SolType) sol; 475 ierr = PetscOptionsEnd(); 476 PetscFunctionReturn(0); 477 } 478 479 static PetscErrorCode SetupParameters(AppCtx *user) 480 { 481 PetscBag bag; 482 Parameter *p; 483 PetscErrorCode ierr; 484 485 PetscFunctionBeginUser; 486 /* setup PETSc parameter bag */ 487 ierr = PetscBagGetData(user->bag, (void **) &p);CHKERRQ(ierr); 488 ierr = PetscBagSetName(user->bag, "par", "Low Mach flow parameters");CHKERRQ(ierr); 489 bag = user->bag; 490 ierr = PetscBagRegisterReal(bag, &p->nu, 1.0, "nu", "Kinematic viscosity");CHKERRQ(ierr); 491 ierr = PetscBagRegisterReal(bag, &p->alpha, 1.0, "alpha", "Thermal diffusivity");CHKERRQ(ierr); 492 ierr = PetscBagRegisterReal(bag, &p->T_in, 1.0, "T_in", "Inlet temperature");CHKERRQ(ierr); 493 PetscFunctionReturn(0); 494 } 495 496 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 497 { 498 PetscErrorCode ierr; 499 500 PetscFunctionBeginUser; 501 ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr); 502 ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 503 ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 504 PetscFunctionReturn(0); 505 } 506 507 static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 508 { 509 PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 510 PetscErrorCode (*exactFuncs_t[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 511 PetscDS prob; 512 Parameter *ctx; 513 PetscInt id; 514 PetscErrorCode ierr; 515 516 PetscFunctionBeginUser; 517 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 518 switch(user->solType){ 519 case SOL_QUADRATIC: 520 ierr = PetscDSSetResidual(prob, 0, f0_quadratic_v, f1_v);CHKERRQ(ierr); 521 ierr = PetscDSSetResidual(prob, 2, f0_quadratic_w, f1_w);CHKERRQ(ierr); 522 523 exactFuncs[0] = quadratic_u; 524 exactFuncs[1] = quadratic_p; 525 exactFuncs[2] = quadratic_T; 526 exactFuncs_t[0] = quadratic_u_t; 527 exactFuncs_t[1] = NULL; 528 exactFuncs_t[2] = quadratic_T_t; 529 break; 530 case SOL_CUBIC: 531 ierr = PetscDSSetResidual(prob, 0, f0_cubic_v, f1_v);CHKERRQ(ierr); 532 ierr = PetscDSSetResidual(prob, 2, f0_cubic_w, f1_w);CHKERRQ(ierr); 533 534 exactFuncs[0] = cubic_u; 535 exactFuncs[1] = cubic_p; 536 exactFuncs[2] = cubic_T; 537 exactFuncs_t[0] = cubic_u_t; 538 exactFuncs_t[1] = NULL; 539 exactFuncs_t[2] = cubic_T_t; 540 break; 541 case SOL_CUBIC_TRIG: 542 ierr = PetscDSSetResidual(prob, 0, f0_cubic_trig_v, f1_v);CHKERRQ(ierr); 543 ierr = PetscDSSetResidual(prob, 2, f0_cubic_trig_w, f1_w);CHKERRQ(ierr); 544 545 exactFuncs[0] = cubic_trig_u; 546 exactFuncs[1] = cubic_trig_p; 547 exactFuncs[2] = cubic_trig_T; 548 exactFuncs_t[0] = cubic_trig_u_t; 549 exactFuncs_t[1] = NULL; 550 exactFuncs_t[2] = cubic_trig_T_t; 551 break; 552 default: SETERRQ2(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType); 553 } 554 555 ierr = PetscDSSetResidual(prob, 1, f0_q, NULL);CHKERRQ(ierr); 556 557 ierr = PetscDSSetJacobian(prob, 0, 0, g0_vu, g1_vu, NULL, g3_vu);CHKERRQ(ierr); 558 ierr = PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_vp, NULL);CHKERRQ(ierr); 559 ierr = PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL, NULL);CHKERRQ(ierr); 560 ierr = PetscDSSetJacobian(prob, 2, 0, g0_wu, NULL, NULL, NULL);CHKERRQ(ierr); 561 ierr = PetscDSSetJacobian(prob, 2, 2, g0_wT, g1_wT, NULL, g3_wT);CHKERRQ(ierr); 562 /* Setup constants */ 563 { 564 Parameter *param; 565 PetscScalar constants[3]; 566 567 ierr = PetscBagGetData(user->bag, (void **) ¶m);CHKERRQ(ierr); 568 569 constants[0] = param->nu; 570 constants[1] = param->alpha; 571 constants[2] = param->T_in; 572 ierr = PetscDSSetConstants(prob, 3, constants);CHKERRQ(ierr); 573 } 574 /* Setup Boundary Conditions */ 575 ierr = PetscBagGetData(user->bag, (void **) &ctx);CHKERRQ(ierr); 576 id = 3; 577 ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall velocity", "marker", 0, 0, NULL, (void (*)(void)) exactFuncs[0], (void (*)(void)) exactFuncs_t[0], 1, &id, ctx);CHKERRQ(ierr); 578 id = 1; 579 ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall velocity", "marker", 0, 0, NULL, (void (*)(void)) exactFuncs[0], (void (*)(void)) exactFuncs_t[0], 1, &id, ctx);CHKERRQ(ierr); 580 id = 2; 581 ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall velocity", "marker", 0, 0, NULL, (void (*)(void)) exactFuncs[0], (void (*)(void)) exactFuncs_t[0], 1, &id, ctx);CHKERRQ(ierr); 582 id = 4; 583 ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall velocity", "marker", 0, 0, NULL, (void (*)(void)) exactFuncs[0], (void (*)(void)) exactFuncs_t[0], 1, &id, ctx);CHKERRQ(ierr); 584 id = 3; 585 ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall temp", "marker", 2, 0, NULL, (void (*)(void)) exactFuncs[2], (void (*)(void)) exactFuncs_t[2], 1, &id, ctx);CHKERRQ(ierr); 586 id = 1; 587 ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall temp", "marker", 2, 0, NULL, (void (*)(void)) exactFuncs[2], (void (*)(void)) exactFuncs_t[2], 1, &id, ctx);CHKERRQ(ierr); 588 id = 2; 589 ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall temp", "marker", 2, 0, NULL, (void (*)(void)) exactFuncs[2], (void (*)(void)) exactFuncs_t[2], 1, &id, ctx);CHKERRQ(ierr); 590 id = 4; 591 ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall temp", "marker", 2, 0, NULL, (void (*)(void)) exactFuncs[2], (void (*)(void)) exactFuncs_t[2], 1, &id, ctx);CHKERRQ(ierr); 592 593 /*setup exact solution.*/ 594 ierr = PetscDSSetExactSolution(prob, 0, exactFuncs[0], ctx);CHKERRQ(ierr); 595 ierr = PetscDSSetExactSolution(prob, 1, exactFuncs[1], ctx);CHKERRQ(ierr); 596 ierr = PetscDSSetExactSolution(prob, 2, exactFuncs[2], ctx);CHKERRQ(ierr); 597 ierr = PetscDSSetExactSolutionTimeDerivative(prob, 0, exactFuncs_t[0], ctx);CHKERRQ(ierr); 598 ierr = PetscDSSetExactSolutionTimeDerivative(prob, 1, exactFuncs_t[1], ctx);CHKERRQ(ierr); 599 ierr = PetscDSSetExactSolutionTimeDerivative(prob, 2, exactFuncs_t[2], ctx);CHKERRQ(ierr); 600 PetscFunctionReturn(0); 601 } 602 603 static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 604 { 605 DM cdm = dm; 606 PetscFE fe[3]; 607 Parameter *param; 608 MPI_Comm comm; 609 DMPolytopeType ct; 610 PetscInt dim, cStart; 611 PetscBool simplex; 612 PetscErrorCode ierr; 613 614 PetscFunctionBeginUser; 615 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 616 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr); 617 ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr); 618 simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 619 /* Create finite element */ 620 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 621 ierr = PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr); 622 ierr = PetscObjectSetName((PetscObject) fe[0], "velocity");CHKERRQ(ierr); 623 624 ierr = PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr); 625 ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr); 626 ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr); 627 628 ierr = PetscFECreateDefault(comm, dim, 1, simplex, "temp_", PETSC_DEFAULT, &fe[2]);CHKERRQ(ierr); 629 ierr = PetscFECopyQuadrature(fe[0], fe[2]);CHKERRQ(ierr); 630 ierr = PetscObjectSetName((PetscObject) fe[2], "temperature");CHKERRQ(ierr); 631 632 /* Set discretization and boundary conditions for each mesh */ 633 ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr); 634 ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr); 635 ierr = DMSetField(dm, 2, NULL, (PetscObject) fe[2]);CHKERRQ(ierr); 636 ierr = DMCreateDS(dm);CHKERRQ(ierr); 637 ierr = SetupProblem(dm, user);CHKERRQ(ierr); 638 ierr = PetscBagGetData(user->bag, (void **) ¶m);CHKERRQ(ierr); 639 while (cdm) { 640 ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr); 641 ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 642 } 643 ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr); 644 ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr); 645 ierr = PetscFEDestroy(&fe[2]);CHKERRQ(ierr); 646 647 { 648 PetscObject pressure; 649 MatNullSpace nullspacePres; 650 651 ierr = DMGetField(dm, 1, NULL, &pressure);CHKERRQ(ierr); 652 ierr = MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres);CHKERRQ(ierr); 653 ierr = PetscObjectCompose(pressure, "nullspace", (PetscObject) nullspacePres);CHKERRQ(ierr); 654 ierr = MatNullSpaceDestroy(&nullspacePres);CHKERRQ(ierr); 655 } 656 657 PetscFunctionReturn(0); 658 } 659 660 static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace) 661 { 662 Vec vec; 663 PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero}; 664 PetscErrorCode ierr; 665 666 PetscFunctionBeginUser; 667 if (ofield != 1) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Nullspace must be for pressure field at index 1, not %D", ofield); 668 funcs[nfield] = constant; 669 ierr = DMCreateGlobalVector(dm, &vec);CHKERRQ(ierr); 670 ierr = DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec);CHKERRQ(ierr); 671 ierr = VecNormalize(vec, NULL);CHKERRQ(ierr); 672 ierr = PetscObjectSetName((PetscObject) vec, "Pressure Null Space");CHKERRQ(ierr); 673 ierr = VecViewFromOptions(vec, NULL, "-pressure_nullspace_view");CHKERRQ(ierr); 674 ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_FALSE, 1, &vec, nullSpace);CHKERRQ(ierr); 675 ierr = VecDestroy(&vec);CHKERRQ(ierr); 676 PetscFunctionReturn(0); 677 } 678 679 static PetscErrorCode RemoveDiscretePressureNullspace_Private(TS ts, Vec u) 680 { 681 DM dm; 682 MatNullSpace nullsp; 683 PetscErrorCode ierr; 684 685 PetscFunctionBegin; 686 ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 687 ierr = CreatePressureNullSpace(dm, 1, 1, &nullsp);CHKERRQ(ierr); 688 ierr = MatNullSpaceRemove(nullsp, u);CHKERRQ(ierr); 689 ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr); 690 PetscFunctionReturn(0); 691 } 692 693 /* Make the discrete pressure discretely divergence free */ 694 static PetscErrorCode RemoveDiscretePressureNullspace(TS ts) 695 { 696 Vec u; 697 PetscErrorCode ierr; 698 699 PetscFunctionBegin; 700 ierr = TSGetSolution(ts, &u);CHKERRQ(ierr); 701 ierr = RemoveDiscretePressureNullspace_Private(ts, u);CHKERRQ(ierr); 702 PetscFunctionReturn(0); 703 } 704 705 static PetscErrorCode SetInitialConditions(TS ts, Vec u) 706 { 707 DM dm; 708 PetscReal t; 709 PetscErrorCode ierr; 710 711 PetscFunctionBegin; 712 ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 713 ierr = TSGetTime(ts, &t);CHKERRQ(ierr); 714 ierr = DMComputeExactSolution(dm, t, u, NULL);CHKERRQ(ierr); 715 ierr = RemoveDiscretePressureNullspace_Private(ts, u);CHKERRQ(ierr); 716 PetscFunctionReturn(0); 717 } 718 719 static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx) 720 { 721 PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 722 void *ctxs[3]; 723 DM dm; 724 PetscDS ds; 725 Vec v; 726 PetscReal ferrors[3]; 727 PetscInt f; 728 PetscErrorCode ierr; 729 730 PetscFunctionBeginUser; 731 ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 732 ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 733 734 for (f = 0; f < 3; ++f) {ierr = PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f]);CHKERRQ(ierr);} 735 ierr = DMComputeL2FieldDiff(dm, crtime, exactFuncs, ctxs, u, ferrors);CHKERRQ(ierr); 736 ierr = PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: [%2.3g, %2.3g, %2.3g]\n", (int) step, (double) crtime, (double) ferrors[0], (double) ferrors[1], (double) ferrors[2]);CHKERRQ(ierr); 737 738 ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr); 739 //ierr = TSGetSolution(ts, &u);CHKERRQ(ierr); 740 ierr = PetscObjectSetName((PetscObject) u, "Numerical Solution");CHKERRQ(ierr); 741 ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr); 742 ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr); 743 744 ierr = DMGetGlobalVector(dm, &v);CHKERRQ(ierr); 745 // ierr = VecSet(v, 0.0);CHKERRQ(ierr); 746 ierr = DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, v);CHKERRQ(ierr); 747 ierr = PetscObjectSetName((PetscObject) v, "Exact Solution");CHKERRQ(ierr); 748 ierr = VecViewFromOptions(v, NULL, "-exact_vec_view");CHKERRQ(ierr); 749 ierr = DMRestoreGlobalVector(dm, &v);CHKERRQ(ierr); 750 751 PetscFunctionReturn(0); 752 } 753 754 int main(int argc, char **argv) 755 { 756 DM dm; /* problem definition */ 757 TS ts; /* timestepper */ 758 Vec u; /* solution */ 759 AppCtx user; /* user-defined work context */ 760 PetscReal t; 761 PetscErrorCode ierr; 762 763 ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 764 ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 765 ierr = PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag);CHKERRQ(ierr); 766 ierr = SetupParameters(&user);CHKERRQ(ierr); 767 ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr); 768 ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 769 ierr = TSSetDM(ts, dm);CHKERRQ(ierr); 770 ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr); 771 /* Setup problem */ 772 ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr); 773 ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr); 774 775 ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 776 ierr = DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace);CHKERRQ(ierr); 777 778 ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user);CHKERRQ(ierr); 779 ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user);CHKERRQ(ierr); 780 ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user);CHKERRQ(ierr); 781 ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 782 ierr = TSSetPreStep(ts, RemoveDiscretePressureNullspace);CHKERRQ(ierr); 783 ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 784 785 ierr = TSSetComputeInitialCondition(ts, SetInitialConditions);CHKERRQ(ierr); /* Must come after SetFromOptions() */ 786 ierr = SetInitialConditions(ts, u);CHKERRQ(ierr); 787 ierr = TSGetTime(ts, &t);CHKERRQ(ierr); 788 ierr = DMSetOutputSequenceNumber(dm, 0, t);CHKERRQ(ierr); 789 ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr); 790 ierr = TSMonitorSet(ts, MonitorError, &user, NULL);CHKERRQ(ierr);CHKERRQ(ierr); 791 792 ierr = TSSolve(ts, u);CHKERRQ(ierr); 793 ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr); 794 ierr = PetscObjectSetName((PetscObject) u, "Numerical Solution");CHKERRQ(ierr); 795 796 ierr = VecDestroy(&u);CHKERRQ(ierr); 797 ierr = DMDestroy(&dm);CHKERRQ(ierr); 798 ierr = TSDestroy(&ts);CHKERRQ(ierr); 799 ierr = PetscBagDestroy(&user.bag);CHKERRQ(ierr); 800 ierr = PetscFinalize(); 801 return ierr; 802 } 803 804 /*TEST 805 806 test: 807 suffix: 2d_tri_p2_p1_p1 808 requires: triangle !single 809 args: -dm_plex_separate_marker -sol_type quadratic -dm_refine 0 \ 810 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 811 -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1 \ 812 -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 813 -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 814 -fieldsplit_0_pc_type lu \ 815 -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 816 817 # TODO Need trig t for convergence in time, also need to refine in space 818 test: 819 # Using -dm_refine 5 -convest_num_refine 2 gives L_2 convergence rate: [0.89, 0.011, 1.0] 820 suffix: 2d_tri_p2_p1_p1_tconv 821 requires: triangle !single 822 args: -dm_plex_separate_marker -sol_type cubic_trig -dm_refine 0 \ 823 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 824 -ts_max_steps 4 -ts_dt 0.1 -ts_convergence_estimate -convest_num_refine 1 \ 825 -snes_error_if_not_converged -snes_convergence_test correct_pressure \ 826 -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 827 -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 828 -fieldsplit_0_pc_type lu \ 829 -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 830 831 test: 832 # Using -dm_refine 3 -convest_num_refine 3 gives L_2 convergence rate: [3.0, 2.5, 1.9] 833 suffix: 2d_tri_p2_p1_p1_sconv 834 requires: triangle !single 835 args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \ 836 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 837 -ts_max_steps 1 -ts_dt 1e-4 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \ 838 -snes_error_if_not_converged -snes_convergence_test correct_pressure \ 839 -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_atol 1e-16 -ksp_error_if_not_converged \ 840 -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 841 -fieldsplit_0_pc_type lu \ 842 -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 843 844 test: 845 suffix: 2d_tri_p3_p2_p2 846 requires: triangle !single 847 args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \ 848 -vel_petscspace_degree 3 -pres_petscspace_degree 2 -temp_petscspace_degree 2 \ 849 -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1 \ 850 -snes_convergence_test correct_pressure \ 851 -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 852 -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 853 -fieldsplit_0_pc_type lu \ 854 -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 855 856 TEST*/ 857