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 332 PetscFunctionBeginUser; 333 options->solType = SOL_QUADRATIC; 334 options->showError = PETSC_FALSE; 335 PetscOptionsBegin(comm, "", "Stokes Problem Options", "DMPLEX"); 336 sol = options->solType; 337 PetscCall(PetscOptionsEList("-sol_type", "The solution type", "ex62.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL)); 338 options->solType = (SolType) sol; 339 PetscCall(PetscOptionsBool("-show_error", "Output the error for verification", "ex62.c", options->showError, &options->showError, NULL)); 340 PetscOptionsEnd(); 341 PetscFunctionReturn(0); 342 } 343 344 static PetscErrorCode SetupParameters(AppCtx *user) 345 { 346 PetscBag bag; 347 Parameter *p; 348 349 PetscFunctionBeginUser; 350 /* setup PETSc parameter bag */ 351 PetscCall(PetscBagGetData(user->bag, (void **) &p)); 352 PetscCall(PetscBagSetName(user->bag, "par", "Poiseuille flow parameters")); 353 bag = user->bag; 354 PetscCall(PetscBagRegisterReal(bag, &p->nu, 1.0, "nu", "Kinematic viscosity")); 355 PetscCall(PetscBagRegisterReal(bag, &p->alpha, 1.0, "alpha", "Thermal diffusivity")); 356 PetscCall(PetscBagRegisterReal(bag, &p->theta, 0.0, "theta", "Angle of pipe wall to x-axis")); 357 PetscFunctionReturn(0); 358 } 359 360 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 361 { 362 PetscFunctionBeginUser; 363 PetscCall(DMCreate(comm, dm)); 364 PetscCall(DMSetType(*dm, DMPLEX)); 365 PetscCall(DMSetFromOptions(*dm)); 366 { 367 Parameter *param; 368 Vec coordinates; 369 PetscScalar *coords; 370 PetscReal theta; 371 PetscInt cdim, N, bs, i; 372 373 PetscCall(DMGetCoordinateDim(*dm, &cdim)); 374 PetscCall(DMGetCoordinates(*dm, &coordinates)); 375 PetscCall(VecGetLocalSize(coordinates, &N)); 376 PetscCall(VecGetBlockSize(coordinates, &bs)); 377 PetscCheck(bs == cdim,comm, PETSC_ERR_ARG_WRONG, "Invalid coordinate blocksize %" PetscInt_FMT " != embedding dimension %" PetscInt_FMT, bs, cdim); 378 PetscCall(VecGetArray(coordinates, &coords)); 379 PetscCall(PetscBagGetData(user->bag, (void **) ¶m)); 380 theta = param->theta; 381 for (i = 0; i < N; i += cdim) { 382 PetscScalar x = coords[i+0]; 383 PetscScalar y = coords[i+1]; 384 385 coords[i+0] = PetscCosReal(theta)*x - PetscSinReal(theta)*y; 386 coords[i+1] = PetscSinReal(theta)*x + PetscCosReal(theta)*y; 387 } 388 PetscCall(VecRestoreArray(coordinates, &coords)); 389 PetscCall(DMSetCoordinates(*dm, coordinates)); 390 } 391 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 392 PetscFunctionReturn(0); 393 } 394 395 static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 396 { 397 PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 398 PetscDS prob; 399 DMLabel label; 400 Parameter *ctx; 401 PetscInt id; 402 403 PetscFunctionBeginUser; 404 PetscCall(DMGetLabel(dm, "marker", &label)); 405 PetscCall(DMGetDS(dm, &prob)); 406 switch(user->solType) { 407 case SOL_QUADRATIC: 408 PetscCall(PetscDSSetResidual(prob, 0, f0_quadratic_v, f1_v)); 409 PetscCall(PetscDSSetResidual(prob, 2, f0_quadratic_w, f1_w)); 410 411 exactFuncs[0] = quadratic_u; 412 exactFuncs[1] = linear_p; 413 exactFuncs[2] = linear_T; 414 break; 415 case SOL_CUBIC: 416 PetscCall(PetscDSSetResidual(prob, 0, f0_cubic_v, f1_v)); 417 PetscCall(PetscDSSetResidual(prob, 2, f0_cubic_w, f1_w)); 418 419 exactFuncs[0] = cubic_u; 420 exactFuncs[1] = quadratic_p; 421 exactFuncs[2] = quadratic_T; 422 break; 423 default: SETERRQ(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%d)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType); 424 } 425 426 PetscCall(PetscDSSetResidual(prob, 1, f0_q, NULL)); 427 428 PetscCall(PetscDSSetJacobian(prob, 0, 0, g0_vu, g1_vu, NULL, g3_vu)); 429 PetscCall(PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_vp, NULL)); 430 PetscCall(PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL, NULL)); 431 PetscCall(PetscDSSetJacobian(prob, 2, 0, g0_wu, NULL, NULL, NULL)); 432 PetscCall(PetscDSSetJacobian(prob, 2, 2, NULL, g1_wT, NULL, g3_wT)); 433 /* Setup constants */ 434 { 435 Parameter *param; 436 PetscScalar constants[3]; 437 438 PetscCall(PetscBagGetData(user->bag, (void **) ¶m)); 439 440 constants[0] = param->nu; 441 constants[1] = param->alpha; 442 constants[2] = param->theta; 443 PetscCall(PetscDSSetConstants(prob, 3, constants)); 444 } 445 /* Setup Boundary Conditions */ 446 PetscCall(PetscBagGetData(user->bag, (void **) &ctx)); 447 id = 3; 448 PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL)); 449 id = 1; 450 PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL)); 451 id = 2; 452 PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL)); 453 id = 4; 454 PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL)); 455 id = 3; 456 PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL)); 457 id = 1; 458 PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL)); 459 id = 2; 460 PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL)); 461 id = 4; 462 PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL)); 463 464 /*setup exact solution.*/ 465 PetscCall(PetscDSSetExactSolution(prob, 0, exactFuncs[0], ctx)); 466 PetscCall(PetscDSSetExactSolution(prob, 1, exactFuncs[1], ctx)); 467 PetscCall(PetscDSSetExactSolution(prob, 2, exactFuncs[2], ctx)); 468 PetscFunctionReturn(0); 469 } 470 471 static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 472 { 473 DM cdm = dm; 474 PetscFE fe[3]; 475 Parameter *param; 476 MPI_Comm comm; 477 PetscInt dim; 478 PetscBool simplex; 479 480 PetscFunctionBeginUser; 481 PetscCall(DMGetDimension(dm, &dim)); 482 PetscCall(DMPlexIsSimplex(dm, &simplex)); 483 /* Create finite element */ 484 PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 485 PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0])); 486 PetscCall(PetscObjectSetName((PetscObject) fe[0], "velocity")); 487 488 PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1])); 489 PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 490 PetscCall(PetscObjectSetName((PetscObject) fe[1], "pressure")); 491 492 PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "temp_", PETSC_DEFAULT, &fe[2])); 493 PetscCall(PetscFECopyQuadrature(fe[0], fe[2])); 494 PetscCall(PetscObjectSetName((PetscObject) fe[2], "temperature")); 495 496 /* Set discretization and boundary conditions for each mesh */ 497 PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe[0])); 498 PetscCall(DMSetField(dm, 1, NULL, (PetscObject) fe[1])); 499 PetscCall(DMSetField(dm, 2, NULL, (PetscObject) fe[2])); 500 PetscCall(DMCreateDS(dm)); 501 PetscCall(SetupProblem(dm, user)); 502 PetscCall(PetscBagGetData(user->bag, (void **) ¶m)); 503 while (cdm) { 504 PetscCall(DMCopyDisc(dm, cdm)); 505 PetscCall(DMPlexCreateBasisRotation(cdm, param->theta, 0.0, 0.0)); 506 PetscCall(DMGetCoarseDM(cdm, &cdm)); 507 } 508 PetscCall(PetscFEDestroy(&fe[0])); 509 PetscCall(PetscFEDestroy(&fe[1])); 510 PetscCall(PetscFEDestroy(&fe[2])); 511 PetscFunctionReturn(0); 512 } 513 514 static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace) 515 { 516 Vec vec; 517 PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero}; 518 519 PetscFunctionBeginUser; 520 PetscCheck(ofield == 1,PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Nullspace must be for pressure field at index 1, not %" PetscInt_FMT, ofield); 521 funcs[nfield] = constant; 522 PetscCall(DMCreateGlobalVector(dm, &vec)); 523 PetscCall(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec)); 524 PetscCall(VecNormalize(vec, NULL)); 525 PetscCall(PetscObjectSetName((PetscObject) vec, "Pressure Null Space")); 526 PetscCall(VecViewFromOptions(vec, NULL, "-pressure_nullspace_view")); 527 PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_FALSE, 1, &vec, nullSpace)); 528 PetscCall(VecDestroy(&vec)); 529 PetscFunctionReturn(0); 530 } 531 532 int main(int argc, char **argv) 533 { 534 SNES snes; /* nonlinear solver */ 535 DM dm; /* problem definition */ 536 Vec u, r; /* solution, residual vectors */ 537 AppCtx user; /* user-defined work context */ 538 539 PetscFunctionBeginUser; 540 PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 541 PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 542 PetscCall(PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag)); 543 PetscCall(SetupParameters(&user)); 544 PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 545 PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 546 PetscCall(SNESSetDM(snes, dm)); 547 PetscCall(DMSetApplicationContext(dm, &user)); 548 /* Setup problem */ 549 PetscCall(SetupDiscretization(dm, &user)); 550 PetscCall(DMPlexCreateClosureIndex(dm, NULL)); 551 552 PetscCall(DMCreateGlobalVector(dm, &u)); 553 PetscCall(PetscObjectSetName((PetscObject) u, "Solution")); 554 PetscCall(VecDuplicate(u, &r)); 555 556 PetscCall(DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace)); 557 PetscCall(DMPlexSetSNESLocalFEM(dm,&user,&user,&user)); 558 559 PetscCall(SNESSetFromOptions(snes)); 560 { 561 PetscDS ds; 562 PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 563 void *ctxs[3]; 564 565 PetscCall(DMGetDS(dm, &ds)); 566 PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], &ctxs[0])); 567 PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], &ctxs[1])); 568 PetscCall(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], &ctxs[2])); 569 PetscCall(DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, u)); 570 PetscCall(PetscObjectSetName((PetscObject) u, "Exact Solution")); 571 PetscCall(VecViewFromOptions(u, NULL, "-exact_vec_view")); 572 } 573 PetscCall(DMSNESCheckFromOptions(snes, u)); 574 PetscCall(VecSet(u, 0.0)); 575 PetscCall(SNESSolve(snes, NULL, u)); 576 577 if (user.showError) { 578 PetscDS ds; 579 Vec r; 580 PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 581 void *ctxs[3]; 582 583 PetscCall(DMGetDS(dm, &ds)); 584 PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], &ctxs[0])); 585 PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], &ctxs[1])); 586 PetscCall(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], &ctxs[2])); 587 PetscCall(DMGetGlobalVector(dm, &r)); 588 PetscCall(DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, r)); 589 PetscCall(VecAXPY(r, -1.0, u)); 590 PetscCall(PetscObjectSetName((PetscObject) r, "Solution Error")); 591 PetscCall(VecViewFromOptions(r, NULL, "-error_vec_view")); 592 PetscCall(DMRestoreGlobalVector(dm, &r)); 593 } 594 PetscCall(PetscObjectSetName((PetscObject) u, "Numerical Solution")); 595 PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 596 597 PetscCall(VecDestroy(&u)); 598 PetscCall(VecDestroy(&r)); 599 PetscCall(DMDestroy(&dm)); 600 PetscCall(SNESDestroy(&snes)); 601 PetscCall(PetscBagDestroy(&user.bag)); 602 PetscCall(PetscFinalize()); 603 return 0; 604 } 605 606 /*TEST 607 608 test: 609 suffix: 2d_tri_p2_p1_p1 610 requires: triangle !single 611 args: -dm_plex_separate_marker -sol_type quadratic -dm_refine 0 \ 612 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 613 -dmsnes_check .001 -snes_error_if_not_converged \ 614 -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 615 -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 616 -fieldsplit_0_pc_type lu \ 617 -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 618 619 test: 620 # Using -dm_refine 2 -convest_num_refine 3 gives L_2 convergence rate: [2.9, 2.3, 1.9] 621 suffix: 2d_tri_p2_p1_p1_conv 622 requires: triangle !single 623 args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \ 624 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 625 -snes_error_if_not_converged -snes_convergence_test correct_pressure -snes_convergence_estimate -convest_num_refine 1 \ 626 -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 627 -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 628 -fieldsplit_0_pc_type lu \ 629 -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 630 631 test: 632 suffix: 2d_tri_p3_p2_p2 633 requires: triangle !single 634 args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \ 635 -vel_petscspace_degree 3 -pres_petscspace_degree 2 -temp_petscspace_degree 2 \ 636 -dmsnes_check .001 -snes_error_if_not_converged \ 637 -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 638 -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 639 -fieldsplit_0_pc_type lu \ 640 -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 641 642 TEST*/ 643