1 static char help[] = "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 a steady-state 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, 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 <petscds.h> 25 #include <petscbag.h> 26 27 typedef enum {SOL_QUADRATIC, SOL_CUBIC, NUM_SOL_TYPES} SolType; 28 const char *solTypes[NUM_SOL_TYPES+1] = {"quadratic", "cubic", "unknown"}; 29 30 typedef struct { 31 PetscReal nu; /* Kinematic viscosity */ 32 PetscReal theta; /* Angle of pipe wall to x-axis */ 33 PetscReal alpha; /* Thermal diffusivity */ 34 PetscReal T_in; /* Inlet temperature*/ 35 } Parameter; 36 37 typedef struct { 38 PetscBool showError; 39 PetscBag bag; 40 SolType solType; 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 = x^2 + y^2 62 v = 2x^2 - 2xy 63 p = x + y - 1 64 T = x + y 65 f = <2x^3 + 4x^2y - 2xy^2 -4\nu + 1, 4xy^2 + 2x^2y - 2y^3 -4\nu + 1> 66 Q = 3x^2 + y^2 - 2xy 67 68 so that 69 70 (1) \nabla \cdot u = 2x - 2x = 0 71 72 (2) u \cdot \nabla u - \nu \Delta u + \nabla p - f 73 = <2x^3 + 4x^2y -2xy^2, 4xy^2 + 2x^2y - 2y^3> -\nu <4, 4> + <1, 1> - <2x^3 + 4x^2y - 2xy^2 -4\nu + 1, 4xy^2 + 2x^2y - 2y^3 - 4\nu + 1> = 0 74 75 (3) u \cdot \nabla T - \alpha \Delta T - Q = 3x^2 + y^2 - 2xy - \alpha*0 - 3x^2 - y^2 + 2xy = 0 76 */ 77 78 static PetscErrorCode quadratic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 79 { 80 u[0] = X[0]*X[0] + X[1]*X[1]; 81 u[1] = 2.0*X[0]*X[0] - 2.0*X[0]*X[1]; 82 return 0; 83 } 84 85 static PetscErrorCode linear_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 86 { 87 p[0] = X[0] + X[1] - 1.0; 88 return 0; 89 } 90 91 static PetscErrorCode linear_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 92 { 93 T[0] = X[0] + X[1]; 94 return 0; 95 } 96 97 static void f0_quadratic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 98 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 99 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 100 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 101 { 102 PetscInt c, d; 103 PetscInt Nc = dim; 104 const PetscReal nu = PetscRealPart(constants[0]); 105 106 for (c=0; c<Nc; ++c) { 107 for (d=0; d<dim; ++d) f0[c] += u[d]*u_x[c*dim+d]; 108 } 109 f0[0] -= (2*X[0]*X[0]*X[0] + 4*X[0]*X[0]*X[1] - 2*X[0]*X[1]*X[1] - 4.0*nu + 1); 110 f0[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 + 1); 111 } 112 113 static void f0_quadratic_w(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 PetscInt d; 119 for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d]; 120 f0[0] -= (3*X[0]*X[0] + X[1]*X[1] - 2*X[0]*X[1]); 121 } 122 123 /* 124 CASE: cubic 125 In 2D we use exact solution: 126 127 u = x^3 + y^3 128 v = 2x^3 - 3x^2y 129 p = 3/2 x^2 + 3/2 y^2 - 1 130 T = 1/2 x^2 + 1/2 y^2 131 f = <3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y), 6x^2y^3 + 3x^4y - 6xy^4 - \nu(12x - 6y) + 3y> 132 Q = x^4 + xy^3 + 2x^3y - 3x^2y^2 - 2 133 134 so that 135 136 \nabla \cdot u = 3x^2 - 3x^2 = 0 137 138 u \cdot \nabla u - \nu \Delta u + \nabla p - f 139 = <3x^5 + 6x^3y^2 - 6x^2y^3, 6x^2y^3 + 3x^4y - 6xy^4> - \nu<6x + 6y, 12x - 6y> + <3x, 3y> - <3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y), 6x^2y^3 + 3x^4y - 6xy^4 - \nu(12x - 6y) + 3y> = 0 140 141 u \cdot \nabla T - \alpha\Delta T - Q = (x^3 + y^3) x + (2x^3 - 3x^2y) y - 2*\alpha - (x^4 + xy^3 + 2x^3y - 3x^2y^2 - 2) = 0 142 */ 143 144 static PetscErrorCode cubic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 145 { 146 u[0] = X[0]*X[0]*X[0] + X[1]*X[1]*X[1]; 147 u[1] = 2.0*X[0]*X[0]*X[0] - 3.0*X[0]*X[0]*X[1]; 148 return 0; 149 } 150 151 static PetscErrorCode quadratic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 152 { 153 p[0] = 3.0*X[0]*X[0]/2.0 + 3.0*X[1]*X[1]/2.0 - 1.0; 154 return 0; 155 } 156 157 static PetscErrorCode quadratic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 158 { 159 T[0] = X[0]*X[0]/2.0 + X[1]*X[1]/2.0; 160 return 0; 161 } 162 163 static void f0_cubic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 164 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 165 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 166 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 167 { 168 PetscInt c, d; 169 PetscInt Nc = dim; 170 const PetscReal nu = PetscRealPart(constants[0]); 171 172 for (c=0; c<Nc; ++c) { 173 for (d=0; d<dim; ++d) f0[c] += u[d]*u_x[c*dim+d]; 174 } 175 f0[0] -= (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]); 176 f0[1] -= (6*X[0]*X[0]*X[1]*X[1]*X[1] + 3*X[0]*X[0]*X[0]*X[0]*X[1] - 6*X[0]*X[1]*X[1]*X[1]*X[1] - (12*X[0] - 6*X[1])*nu + 3*X[1]); 177 } 178 179 static void f0_cubic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 180 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 181 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 182 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 183 { 184 const PetscReal alpha = PetscRealPart(constants[1]); 185 PetscInt d; 186 187 for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d]; 188 f0[0] -= (X[0]*X[0]*X[0]*X[0] + X[0]*X[1]*X[1]*X[1] + 2.0*X[0]*X[0]*X[0]*X[1] - 3.0*X[0]*X[0]*X[1]*X[1] - 2.0*alpha); 189 } 190 191 static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 192 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 193 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 194 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 195 { 196 PetscInt d; 197 for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d]; 198 } 199 200 static void f1_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 201 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 202 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 203 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 204 { 205 const PetscReal nu = PetscRealPart(constants[0]); 206 const PetscInt Nc = dim; 207 PetscInt c, d; 208 209 for (c = 0; c < Nc; ++c) { 210 for (d = 0; d < dim; ++d) { 211 f1[c*dim+d] = nu*(u_x[c*dim+d] + u_x[d*dim+c]); 212 //f1[c*dim+d] = nu*u_x[c*dim+d]; 213 } 214 f1[c*dim+c] -= u[uOff[1]]; 215 } 216 } 217 218 static void f1_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 219 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 220 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 221 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 222 { 223 const PetscReal alpha = PetscRealPart(constants[1]); 224 PetscInt d; 225 for (d = 0; d < dim; ++d) f1[d] = alpha*u_x[uOff_x[2]+d]; 226 } 227 228 /*Jacobians*/ 229 static void g1_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 230 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 231 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 232 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 233 { 234 PetscInt d; 235 for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; 236 } 237 238 static void g0_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 239 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 240 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 241 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 242 { 243 const PetscInt Nc = dim; 244 PetscInt c, d; 245 246 for (c = 0; c < Nc; ++c) { 247 for (d = 0; d < dim; ++d) { 248 g0[c*Nc+d] = u_x[ c*Nc+d]; 249 } 250 } 251 } 252 253 static void g1_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 254 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 255 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 256 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 257 { 258 PetscInt NcI = dim; 259 PetscInt NcJ = dim; 260 PetscInt c, d, e; 261 262 for (c = 0; c < NcI; ++c) { 263 for (d = 0; d < NcJ; ++d) { 264 for (e = 0; e < dim; ++e) { 265 if (c == d) { 266 g1[(c*NcJ+d)*dim+e] = u[e]; 267 } 268 } 269 } 270 } 271 } 272 273 static void g2_vp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 274 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 275 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 276 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 277 { 278 PetscInt d; 279 for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; 280 } 281 282 static void g3_vu(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 286 { 287 const PetscReal nu = PetscRealPart(constants[0]); 288 const PetscInt Nc = dim; 289 PetscInt c, d; 290 291 for (c = 0; c < Nc; ++c) { 292 for (d = 0; d < dim; ++d) { 293 g3[((c*Nc+c)*dim+d)*dim+d] += nu; // gradU 294 g3[((c*Nc+d)*dim+d)*dim+c] += nu; // gradU transpose 295 } 296 } 297 } 298 299 static void g0_wu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 300 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 301 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 302 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 303 { 304 PetscInt d; 305 for (d = 0; d < dim; ++d) g0[d] = u_x[uOff_x[2]+d]; 306 } 307 308 static void g1_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 309 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 310 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 311 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 312 { 313 PetscInt d; 314 for (d = 0; d < dim; ++d) g1[d] = u[uOff[0]+d]; 315 } 316 317 static void g3_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 318 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 319 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 320 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 321 { 322 const PetscReal alpha = PetscRealPart(constants[1]); 323 PetscInt d; 324 325 for (d = 0; d < dim; ++d) g3[d*dim+d] = alpha; 326 } 327 328 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 329 { 330 PetscInt sol; 331 PetscErrorCode ierr; 332 333 PetscFunctionBeginUser; 334 options->solType = SOL_QUADRATIC; 335 options->showError = PETSC_FALSE; 336 337 ierr = PetscOptionsBegin(comm, "", "Stokes Problem Options", "DMPLEX");CHKERRQ(ierr); 338 sol = options->solType; 339 CHKERRQ(PetscOptionsEList("-sol_type", "The solution type", "ex62.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL)); 340 options->solType = (SolType) sol; 341 CHKERRQ(PetscOptionsBool("-show_error", "Output the error for verification", "ex62.c", options->showError, &options->showError, NULL)); 342 ierr = PetscOptionsEnd();CHKERRQ(ierr); 343 PetscFunctionReturn(0); 344 } 345 346 static PetscErrorCode SetupParameters(AppCtx *user) 347 { 348 PetscBag bag; 349 Parameter *p; 350 351 PetscFunctionBeginUser; 352 /* setup PETSc parameter bag */ 353 CHKERRQ(PetscBagGetData(user->bag, (void **) &p)); 354 CHKERRQ(PetscBagSetName(user->bag, "par", "Poiseuille flow parameters")); 355 bag = user->bag; 356 CHKERRQ(PetscBagRegisterReal(bag, &p->nu, 1.0, "nu", "Kinematic viscosity")); 357 CHKERRQ(PetscBagRegisterReal(bag, &p->alpha, 1.0, "alpha", "Thermal diffusivity")); 358 CHKERRQ(PetscBagRegisterReal(bag, &p->theta, 0.0, "theta", "Angle of pipe wall to x-axis")); 359 PetscFunctionReturn(0); 360 } 361 362 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 363 { 364 PetscFunctionBeginUser; 365 CHKERRQ(DMCreate(comm, dm)); 366 CHKERRQ(DMSetType(*dm, DMPLEX)); 367 CHKERRQ(DMSetFromOptions(*dm)); 368 { 369 Parameter *param; 370 Vec coordinates; 371 PetscScalar *coords; 372 PetscReal theta; 373 PetscInt cdim, N, bs, i; 374 375 CHKERRQ(DMGetCoordinateDim(*dm, &cdim)); 376 CHKERRQ(DMGetCoordinates(*dm, &coordinates)); 377 CHKERRQ(VecGetLocalSize(coordinates, &N)); 378 CHKERRQ(VecGetBlockSize(coordinates, &bs)); 379 PetscCheckFalse(bs != cdim,comm, PETSC_ERR_ARG_WRONG, "Invalid coordinate blocksize %D != embedding dimension %D", bs, cdim); 380 CHKERRQ(VecGetArray(coordinates, &coords)); 381 CHKERRQ(PetscBagGetData(user->bag, (void **) ¶m)); 382 theta = param->theta; 383 for (i = 0; i < N; i += cdim) { 384 PetscScalar x = coords[i+0]; 385 PetscScalar y = coords[i+1]; 386 387 coords[i+0] = PetscCosReal(theta)*x - PetscSinReal(theta)*y; 388 coords[i+1] = PetscSinReal(theta)*x + PetscCosReal(theta)*y; 389 } 390 CHKERRQ(VecRestoreArray(coordinates, &coords)); 391 CHKERRQ(DMSetCoordinates(*dm, coordinates)); 392 } 393 CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 394 PetscFunctionReturn(0); 395 } 396 397 static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 398 { 399 PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 400 PetscDS prob; 401 DMLabel label; 402 Parameter *ctx; 403 PetscInt id; 404 405 PetscFunctionBeginUser; 406 CHKERRQ(DMGetLabel(dm, "marker", &label)); 407 CHKERRQ(DMGetDS(dm, &prob)); 408 switch(user->solType) { 409 case SOL_QUADRATIC: 410 CHKERRQ(PetscDSSetResidual(prob, 0, f0_quadratic_v, f1_v)); 411 CHKERRQ(PetscDSSetResidual(prob, 2, f0_quadratic_w, f1_w)); 412 413 exactFuncs[0] = quadratic_u; 414 exactFuncs[1] = linear_p; 415 exactFuncs[2] = linear_T; 416 break; 417 case SOL_CUBIC: 418 CHKERRQ(PetscDSSetResidual(prob, 0, f0_cubic_v, f1_v)); 419 CHKERRQ(PetscDSSetResidual(prob, 2, f0_cubic_w, f1_w)); 420 421 exactFuncs[0] = cubic_u; 422 exactFuncs[1] = quadratic_p; 423 exactFuncs[2] = quadratic_T; 424 break; 425 default: SETERRQ(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType); 426 } 427 428 CHKERRQ(PetscDSSetResidual(prob, 1, f0_q, NULL)); 429 430 CHKERRQ(PetscDSSetJacobian(prob, 0, 0, g0_vu, g1_vu, NULL, g3_vu)); 431 CHKERRQ(PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_vp, NULL)); 432 CHKERRQ(PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL, NULL)); 433 CHKERRQ(PetscDSSetJacobian(prob, 2, 0, g0_wu, NULL, NULL, NULL)); 434 CHKERRQ(PetscDSSetJacobian(prob, 2, 2, NULL, g1_wT, NULL, g3_wT)); 435 /* Setup constants */ 436 { 437 Parameter *param; 438 PetscScalar constants[3]; 439 440 CHKERRQ(PetscBagGetData(user->bag, (void **) ¶m)); 441 442 constants[0] = param->nu; 443 constants[1] = param->alpha; 444 constants[2] = param->theta; 445 CHKERRQ(PetscDSSetConstants(prob, 3, constants)); 446 } 447 /* Setup Boundary Conditions */ 448 CHKERRQ(PetscBagGetData(user->bag, (void **) &ctx)); 449 id = 3; 450 CHKERRQ(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL)); 451 id = 1; 452 CHKERRQ(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL)); 453 id = 2; 454 CHKERRQ(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL)); 455 id = 4; 456 CHKERRQ(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL)); 457 id = 3; 458 CHKERRQ(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL)); 459 id = 1; 460 CHKERRQ(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL)); 461 id = 2; 462 CHKERRQ(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL)); 463 id = 4; 464 CHKERRQ(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL)); 465 466 /*setup exact solution.*/ 467 CHKERRQ(PetscDSSetExactSolution(prob, 0, exactFuncs[0], ctx)); 468 CHKERRQ(PetscDSSetExactSolution(prob, 1, exactFuncs[1], ctx)); 469 CHKERRQ(PetscDSSetExactSolution(prob, 2, exactFuncs[2], ctx)); 470 PetscFunctionReturn(0); 471 } 472 473 static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 474 { 475 DM cdm = dm; 476 PetscFE fe[3]; 477 Parameter *param; 478 MPI_Comm comm; 479 PetscInt dim; 480 PetscBool simplex; 481 482 PetscFunctionBeginUser; 483 CHKERRQ(DMGetDimension(dm, &dim)); 484 CHKERRQ(DMPlexIsSimplex(dm, &simplex)); 485 /* Create finite element */ 486 CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm)); 487 CHKERRQ(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0])); 488 CHKERRQ(PetscObjectSetName((PetscObject) fe[0], "velocity")); 489 490 CHKERRQ(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1])); 491 CHKERRQ(PetscFECopyQuadrature(fe[0], fe[1])); 492 CHKERRQ(PetscObjectSetName((PetscObject) fe[1], "pressure")); 493 494 CHKERRQ(PetscFECreateDefault(comm, dim, 1, simplex, "temp_", PETSC_DEFAULT, &fe[2])); 495 CHKERRQ(PetscFECopyQuadrature(fe[0], fe[2])); 496 CHKERRQ(PetscObjectSetName((PetscObject) fe[2], "temperature")); 497 498 /* Set discretization and boundary conditions for each mesh */ 499 CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe[0])); 500 CHKERRQ(DMSetField(dm, 1, NULL, (PetscObject) fe[1])); 501 CHKERRQ(DMSetField(dm, 2, NULL, (PetscObject) fe[2])); 502 CHKERRQ(DMCreateDS(dm)); 503 CHKERRQ(SetupProblem(dm, user)); 504 CHKERRQ(PetscBagGetData(user->bag, (void **) ¶m)); 505 while (cdm) { 506 CHKERRQ(DMCopyDisc(dm, cdm)); 507 CHKERRQ(DMPlexCreateBasisRotation(cdm, param->theta, 0.0, 0.0)); 508 CHKERRQ(DMGetCoarseDM(cdm, &cdm)); 509 } 510 CHKERRQ(PetscFEDestroy(&fe[0])); 511 CHKERRQ(PetscFEDestroy(&fe[1])); 512 CHKERRQ(PetscFEDestroy(&fe[2])); 513 PetscFunctionReturn(0); 514 } 515 516 static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace) 517 { 518 Vec vec; 519 PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero}; 520 521 PetscFunctionBeginUser; 522 PetscCheckFalse(ofield != 1,PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Nullspace must be for pressure field at index 1, not %D", ofield); 523 funcs[nfield] = constant; 524 CHKERRQ(DMCreateGlobalVector(dm, &vec)); 525 CHKERRQ(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec)); 526 CHKERRQ(VecNormalize(vec, NULL)); 527 CHKERRQ(PetscObjectSetName((PetscObject) vec, "Pressure Null Space")); 528 CHKERRQ(VecViewFromOptions(vec, NULL, "-pressure_nullspace_view")); 529 CHKERRQ(MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_FALSE, 1, &vec, nullSpace)); 530 CHKERRQ(VecDestroy(&vec)); 531 PetscFunctionReturn(0); 532 } 533 534 int main(int argc, char **argv) 535 { 536 SNES snes; /* nonlinear solver */ 537 DM dm; /* problem definition */ 538 Vec u, r; /* solution, residual vectors */ 539 AppCtx user; /* user-defined work context */ 540 PetscErrorCode ierr; 541 542 ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 543 CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user)); 544 CHKERRQ(PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag)); 545 CHKERRQ(SetupParameters(&user)); 546 CHKERRQ(SNESCreate(PETSC_COMM_WORLD, &snes)); 547 CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 548 CHKERRQ(SNESSetDM(snes, dm)); 549 CHKERRQ(DMSetApplicationContext(dm, &user)); 550 /* Setup problem */ 551 CHKERRQ(SetupDiscretization(dm, &user)); 552 CHKERRQ(DMPlexCreateClosureIndex(dm, NULL)); 553 554 CHKERRQ(DMCreateGlobalVector(dm, &u)); 555 CHKERRQ(PetscObjectSetName((PetscObject) u, "Solution")); 556 CHKERRQ(VecDuplicate(u, &r)); 557 558 CHKERRQ(DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace)); 559 CHKERRQ(DMPlexSetSNESLocalFEM(dm,&user,&user,&user)); 560 561 CHKERRQ(SNESSetFromOptions(snes)); 562 { 563 PetscDS ds; 564 PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 565 void *ctxs[3]; 566 567 CHKERRQ(DMGetDS(dm, &ds)); 568 CHKERRQ(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], &ctxs[0])); 569 CHKERRQ(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], &ctxs[1])); 570 CHKERRQ(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], &ctxs[2])); 571 CHKERRQ(DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, u)); 572 CHKERRQ(PetscObjectSetName((PetscObject) u, "Exact Solution")); 573 CHKERRQ(VecViewFromOptions(u, NULL, "-exact_vec_view")); 574 } 575 CHKERRQ(DMSNESCheckFromOptions(snes, u)); 576 CHKERRQ(VecSet(u, 0.0)); 577 CHKERRQ(SNESSolve(snes, NULL, u)); 578 579 if (user.showError) { 580 PetscDS ds; 581 Vec r; 582 PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 583 void *ctxs[3]; 584 585 CHKERRQ(DMGetDS(dm, &ds)); 586 CHKERRQ(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], &ctxs[0])); 587 CHKERRQ(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], &ctxs[1])); 588 CHKERRQ(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], &ctxs[2])); 589 CHKERRQ(DMGetGlobalVector(dm, &r)); 590 CHKERRQ(DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, r)); 591 CHKERRQ(VecAXPY(r, -1.0, u)); 592 CHKERRQ(PetscObjectSetName((PetscObject) r, "Solution Error")); 593 CHKERRQ(VecViewFromOptions(r, NULL, "-error_vec_view")); 594 CHKERRQ(DMRestoreGlobalVector(dm, &r)); 595 } 596 CHKERRQ(PetscObjectSetName((PetscObject) u, "Numerical Solution")); 597 CHKERRQ(VecViewFromOptions(u, NULL, "-sol_vec_view")); 598 599 CHKERRQ(VecDestroy(&u)); 600 CHKERRQ(VecDestroy(&r)); 601 CHKERRQ(DMDestroy(&dm)); 602 CHKERRQ(SNESDestroy(&snes)); 603 CHKERRQ(PetscBagDestroy(&user.bag)); 604 ierr = PetscFinalize(); 605 return ierr; 606 } 607 608 /*TEST 609 610 test: 611 suffix: 2d_tri_p2_p1_p1 612 requires: triangle !single 613 args: -dm_plex_separate_marker -sol_type quadratic -dm_refine 0 \ 614 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 615 -dmsnes_check .001 -snes_error_if_not_converged \ 616 -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 617 -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 618 -fieldsplit_0_pc_type lu \ 619 -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 620 621 test: 622 # Using -dm_refine 2 -convest_num_refine 3 gives L_2 convergence rate: [2.9, 2.3, 1.9] 623 suffix: 2d_tri_p2_p1_p1_conv 624 requires: triangle !single 625 args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \ 626 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 627 -snes_error_if_not_converged -snes_convergence_test correct_pressure -snes_convergence_estimate -convest_num_refine 1 \ 628 -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 629 -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 630 -fieldsplit_0_pc_type lu \ 631 -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 632 633 test: 634 suffix: 2d_tri_p3_p2_p2 635 requires: triangle !single 636 args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \ 637 -vel_petscspace_degree 3 -pres_petscspace_degree 2 -temp_petscspace_degree 2 \ 638 -dmsnes_check .001 -snes_error_if_not_converged \ 639 -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 640 -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 641 -fieldsplit_0_pc_type lu \ 642 -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 643 644 TEST*/ 645