1 static char help[] = "Time-dependent reactive low Mach Flow in 2d and 3d channels with finite elements.\n\ 2 We solve the reactive low Mach flow problem in a rectangular domain\n\ 3 using a parallel unstructured mesh (DMPLEX) to discretize the flow\n\ 4 and particles (DWSWARM) to discretize the chemical species.\n\n\n"; 5 6 /*F 7 This low Mach flow is time-dependent isoviscous Navier-Stokes flow. We discretize using the 8 finite element method on an unstructured mesh. The weak form equations are 9 10 \begin{align*} 11 < q, \nabla\cdot u > = 0 12 <v, du/dt> + <v, u \cdot \nabla u> + < \nabla v, \nu (\nabla u + {\nabla u}^T) > - < \nabla\cdot v, p > - < v, f > = 0 13 < w, u \cdot \nabla T > + < \nabla w, \alpha \nabla T > - < w, Q > = 0 14 \end{align*} 15 16 where $\nu$ is the kinematic viscosity and $\alpha$ is thermal diffusivity. 17 18 For visualization, use 19 20 -dm_view hdf5:$PWD/sol.h5 -sol_vec_view hdf5:$PWD/sol.h5::append -exact_vec_view hdf5:$PWD/sol.h5::append 21 22 The particles can be visualized using 23 24 -part_dm_view draw -part_dm_view_swarm_radius 0.03 25 26 F*/ 27 28 #include <petscdmplex.h> 29 #include <petscdmswarm.h> 30 #include <petscts.h> 31 #include <petscds.h> 32 #include <petscbag.h> 33 34 typedef enum {SOL_TRIG_TRIG, NUM_SOL_TYPES} SolType; 35 const char *solTypes[NUM_SOL_TYPES+1] = {"trig_trig", "unknown"}; 36 37 typedef enum {PART_LAYOUT_CELL, PART_LAYOUT_BOX, NUM_PART_LAYOUT_TYPES} PartLayoutType; 38 const char *partLayoutTypes[NUM_PART_LAYOUT_TYPES+1] = {"cell", "box", "unknown"}; 39 40 typedef struct { 41 PetscReal nu; /* Kinematic viscosity */ 42 PetscReal alpha; /* Thermal diffusivity */ 43 PetscReal T_in; /* Inlet temperature*/ 44 PetscReal omega; /* Rotation speed in MMS benchmark */ 45 } Parameter; 46 47 typedef struct { 48 /* Problem definition */ 49 PetscBag bag; /* Holds problem parameters */ 50 SolType solType; /* MMS solution type */ 51 PartLayoutType partLayout; /* Type of particle distribution */ 52 PetscInt Npc; /* The initial number of particles per cell */ 53 PetscReal partLower[3]; /* Lower left corner of particle box */ 54 PetscReal partUpper[3]; /* Upper right corner of particle box */ 55 PetscInt Npb; /* The initial number of particles per box dimension */ 56 } AppCtx; 57 58 typedef struct { 59 PetscReal ti; /* The time for ui, at the beginning of the advection solve */ 60 PetscReal tf; /* The time for uf, at the end of the advection solve */ 61 Vec ui; /* The PDE solution field at ti */ 62 Vec uf; /* The PDE solution field at tf */ 63 Vec x0; /* The initial particle positions at t = 0 */ 64 PetscErrorCode (*exact)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 65 AppCtx *ctx; /* Context for exact solution */ 66 } AdvCtx; 67 68 static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 69 { 70 PetscInt d; 71 for (d = 0; d < Nc; ++d) u[d] = 0.0; 72 return 0; 73 } 74 75 static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 76 { 77 PetscInt d; 78 for (d = 0; d < Nc; ++d) u[d] = 1.0; 79 return 0; 80 } 81 82 /* 83 CASE: trigonometric-trigonometric 84 In 2D we use exact solution: 85 86 x = r0 cos(w t + theta0) r0 = sqrt(x0^2 + y0^2) 87 y = r0 sin(w t + theta0) theta0 = arctan(y0/x0) 88 u = -w r0 sin(theta0) = -w y 89 v = w r0 cos(theta0) = w x 90 p = x + y - 1 91 T = t + x + y 92 f = <1, 1> 93 Q = 1 + w (x - y)/r 94 95 so that 96 97 \nabla \cdot u = 0 + 0 = 0 98 99 f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p 100 = <0, 0> + u_i d_i u_j - \nu 0 + <1, 1> 101 = <1, 1> + w^2 <-y, x> . <<0, 1>, <-1, 0>> 102 = <1, 1> + w^2 <-x, -y> 103 = <1, 1> - w^2 <x, y> 104 105 Q = dT/dt + u \cdot \nabla T - \alpha \Delta T 106 = 1 + <u, v> . <1, 1> - \alpha 0 107 = 1 + u + v 108 */ 109 static PetscErrorCode trig_trig_x(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *x, void *ctx) 110 { 111 const PetscReal x0 = X[0]; 112 const PetscReal y0 = X[1]; 113 const PetscReal R0 = PetscSqrtReal(x0*x0 + y0*y0); 114 const PetscReal theta0 = PetscAtan2Real(y0, x0); 115 Parameter *p = (Parameter *) ctx; 116 117 x[0] = R0*PetscCosReal(p->omega*time + theta0); 118 x[1] = R0*PetscSinReal(p->omega*time + theta0); 119 return 0; 120 } 121 static PetscErrorCode trig_trig_u(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 122 { 123 Parameter *p = (Parameter *) ctx; 124 125 u[0] = -p->omega*X[1]; 126 u[1] = p->omega*X[0]; 127 return 0; 128 } 129 static PetscErrorCode trig_trig_u_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) 130 { 131 u[0] = 0.0; 132 u[1] = 0.0; 133 return 0; 134 } 135 136 static PetscErrorCode trig_trig_p(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) 137 { 138 p[0] = X[0] + X[1] - 1.0; 139 return 0; 140 } 141 142 static PetscErrorCode trig_trig_T(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 143 { 144 T[0] = time + X[0] + X[1]; 145 return 0; 146 } 147 static PetscErrorCode trig_trig_T_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) 148 { 149 T[0] = 1.0; 150 return 0; 151 } 152 153 static void f0_trig_trig_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 154 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 155 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 156 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 157 { 158 const PetscReal omega = PetscRealPart(constants[3]); 159 PetscInt Nc = dim; 160 PetscInt c, d; 161 162 for (d = 0; d < dim; ++d) f0[d] = u_t[uOff[0]+d]; 163 164 for (c = 0; c < Nc; ++c) { 165 for (d = 0; d < dim; ++d) f0[c] += u[d]*u_x[c*dim+d]; 166 } 167 f0[0] -= 1.0 - omega*omega*X[0]; 168 f0[1] -= 1.0 - omega*omega*X[1]; 169 } 170 171 static void f0_trig_trig_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 172 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 173 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 174 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 175 { 176 const PetscReal omega = PetscRealPart(constants[3]); 177 PetscInt d; 178 179 for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d]; 180 f0[0] += u_t[uOff[2]] - (1.0 + omega*(X[0] - X[1])); 181 } 182 183 static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 184 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 185 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 186 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 187 { 188 PetscInt d; 189 for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d]; 190 } 191 192 /*f1_v = \nu[grad(u) + grad(u)^T] - pI */ 193 static void f1_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, 194 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 195 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 196 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 197 { 198 const PetscReal nu = PetscRealPart(constants[0]); 199 const PetscInt Nc = dim; 200 PetscInt c, d; 201 202 for (c = 0; c < Nc; ++c) { 203 for (d = 0; d < dim; ++d) { 204 f1[c*dim+d] = nu*(u_x[c*dim+d] + u_x[d*dim+c]); 205 } 206 f1[c*dim+c] -= u[uOff[1]]; 207 } 208 } 209 210 static void f1_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, 211 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 212 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 213 PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 214 { 215 const PetscReal alpha = PetscRealPart(constants[1]); 216 PetscInt d; 217 for (d = 0; d < dim; ++d) f1[d] = alpha*u_x[uOff_x[2]+d]; 218 } 219 220 /*Jacobians*/ 221 static void g1_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 222 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 223 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 224 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 225 { 226 PetscInt d; 227 for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; 228 } 229 230 static void g0_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 231 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 232 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 233 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 234 { 235 PetscInt c, d; 236 const PetscInt Nc = dim; 237 238 for (d = 0; d < dim; ++d) g0[d*dim+d] = u_tShift; 239 240 for (c = 0; c < Nc; ++c) { 241 for (d = 0; d < dim; ++d) { 242 g0[c*Nc+d] += u_x[c*Nc+d]; 243 } 244 } 245 } 246 247 static void g1_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 248 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 249 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 250 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 251 { 252 PetscInt NcI = dim; 253 PetscInt NcJ = dim; 254 PetscInt c, d, e; 255 256 for (c = 0; c < NcI; ++c) { 257 for (d = 0; d < NcJ; ++d) { 258 for (e = 0; e < dim; ++e) { 259 if (c == d) { 260 g1[(c*NcJ+d)*dim+e] += u[e]; 261 } 262 } 263 } 264 } 265 } 266 267 268 static void g2_vp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 269 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 270 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 271 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 272 { 273 PetscInt d; 274 for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; 275 } 276 277 static void g3_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 278 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 279 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 280 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 281 { 282 const PetscReal nu = PetscRealPart(constants[0]); 283 const PetscInt Nc = dim; 284 PetscInt c, d; 285 286 for (c = 0; c < Nc; ++c) { 287 for (d = 0; d < dim; ++d) { 288 g3[((c*Nc+c)*dim+d)*dim+d] += nu; // gradU 289 g3[((c*Nc+d)*dim+d)*dim+c] += nu; // gradU transpose 290 } 291 } 292 } 293 294 static void g0_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 295 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 296 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 297 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 298 { 299 PetscInt d; 300 for (d = 0; d < dim; ++d) g0[d] = u_tShift; 301 } 302 303 static void g0_wu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 304 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 305 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 306 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 307 { 308 PetscInt d; 309 for (d = 0; d < dim; ++d) g0[d] = u_x[uOff_x[2]+d]; 310 } 311 312 static void g1_wT(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 316 { 317 PetscInt d; 318 for (d = 0; d < dim; ++d) g1[d] = u[uOff[0]+d]; 319 } 320 321 static void g3_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 322 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 323 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 324 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 325 { 326 const PetscReal alpha = PetscRealPart(constants[1]); 327 PetscInt d; 328 329 for (d = 0; d < dim; ++d) g3[d*dim+d] = alpha; 330 } 331 332 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 333 { 334 PetscInt sol, pl, n; 335 PetscErrorCode ierr; 336 337 PetscFunctionBeginUser; 338 options->solType = SOL_TRIG_TRIG; 339 options->partLayout = PART_LAYOUT_CELL; 340 options->Npc = 1; 341 options->Npb = 1; 342 343 options->partLower[0] = options->partLower[1] = options->partLower[2] = 0.; 344 options->partUpper[0] = options->partUpper[1] = options->partUpper[2] = 1.; 345 ierr = PetscOptionsBegin(comm, "", "Low Mach flow Problem Options", "DMPLEX");CHKERRQ(ierr); 346 sol = options->solType; 347 ierr = PetscOptionsEList("-sol_type", "The solution type", "ex77.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL);CHKERRQ(ierr); 348 options->solType = (SolType) sol; 349 pl = options->partLayout; 350 ierr = PetscOptionsEList("-part_layout_type", "The particle layout type", "ex77.c", partLayoutTypes, NUM_PART_LAYOUT_TYPES, partLayoutTypes[options->partLayout], &pl, NULL);CHKERRQ(ierr); 351 options->partLayout = (PartLayoutType) pl; 352 ierr = PetscOptionsInt("-Npc", "The initial number of particles per cell", "ex77.c", options->Npc, &options->Npc, NULL);CHKERRQ(ierr); 353 n = 3; 354 ierr = PetscOptionsRealArray("-part_lower", "The lower left corner of the particle box", "ex77.c", options->partLower, &n, NULL);CHKERRQ(ierr); 355 n = 3; 356 ierr = PetscOptionsRealArray("-part_upper", "The upper right corner of the particle box", "ex77.c", options->partUpper, &n, NULL);CHKERRQ(ierr); 357 ierr = PetscOptionsInt("-Npb", "The initial number of particles per box dimension", "ex77.c", options->Npb, &options->Npb, NULL);CHKERRQ(ierr); 358 ierr = PetscOptionsEnd(); 359 PetscFunctionReturn(0); 360 } 361 362 static PetscErrorCode SetupParameters(AppCtx *user) 363 { 364 PetscBag bag; 365 Parameter *p; 366 PetscErrorCode ierr; 367 368 PetscFunctionBeginUser; 369 /* setup PETSc parameter bag */ 370 ierr = PetscBagGetData(user->bag, (void **) &p);CHKERRQ(ierr); 371 ierr = PetscBagSetName(user->bag, "par", "Low Mach flow parameters");CHKERRQ(ierr); 372 bag = user->bag; 373 ierr = PetscBagRegisterReal(bag, &p->nu, 1.0, "nu", "Kinematic viscosity");CHKERRQ(ierr); 374 ierr = PetscBagRegisterReal(bag, &p->alpha, 1.0, "alpha", "Thermal diffusivity");CHKERRQ(ierr); 375 ierr = PetscBagRegisterReal(bag, &p->T_in, 1.0, "T_in", "Inlet temperature");CHKERRQ(ierr); 376 ierr = PetscBagRegisterReal(bag, &p->omega, 1.0, "omega", "Rotation speed in MMS benchmark");CHKERRQ(ierr); 377 PetscFunctionReturn(0); 378 } 379 380 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 381 { 382 PetscErrorCode ierr; 383 384 PetscFunctionBeginUser; 385 ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr); 386 ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 387 ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 388 PetscFunctionReturn(0); 389 } 390 391 static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 392 { 393 PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 394 PetscErrorCode (*exactFuncs_t[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 395 PetscDS prob; 396 Parameter *ctx; 397 PetscInt id; 398 PetscErrorCode ierr; 399 400 PetscFunctionBeginUser; 401 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 402 switch(user->solType){ 403 case SOL_TRIG_TRIG: 404 ierr = PetscDSSetResidual(prob, 0, f0_trig_trig_v, f1_v);CHKERRQ(ierr); 405 ierr = PetscDSSetResidual(prob, 2, f0_trig_trig_w, f1_w);CHKERRQ(ierr); 406 407 exactFuncs[0] = trig_trig_u; 408 exactFuncs[1] = trig_trig_p; 409 exactFuncs[2] = trig_trig_T; 410 exactFuncs_t[0] = trig_trig_u_t; 411 exactFuncs_t[1] = NULL; 412 exactFuncs_t[2] = trig_trig_T_t; 413 break; 414 default: SETERRQ2(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType); 415 } 416 417 ierr = PetscDSSetResidual(prob, 1, f0_q, NULL);CHKERRQ(ierr); 418 419 ierr = PetscDSSetJacobian(prob, 0, 0, g0_vu, g1_vu, NULL, g3_vu);CHKERRQ(ierr); 420 ierr = PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_vp, NULL);CHKERRQ(ierr); 421 ierr = PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL, NULL);CHKERRQ(ierr); 422 ierr = PetscDSSetJacobian(prob, 2, 0, g0_wu, NULL, NULL, NULL);CHKERRQ(ierr); 423 ierr = PetscDSSetJacobian(prob, 2, 2, g0_wT, g1_wT, NULL, g3_wT);CHKERRQ(ierr); 424 /* Setup constants */ 425 { 426 Parameter *param; 427 PetscScalar constants[4]; 428 429 ierr = PetscBagGetData(user->bag, (void **) ¶m);CHKERRQ(ierr); 430 431 constants[0] = param->nu; 432 constants[1] = param->alpha; 433 constants[2] = param->T_in; 434 constants[3] = param->omega; 435 ierr = PetscDSSetConstants(prob, 4, constants);CHKERRQ(ierr); 436 } 437 /* Setup Boundary Conditions */ 438 ierr = PetscBagGetData(user->bag, (void **) &ctx);CHKERRQ(ierr); 439 id = 3; 440 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); 441 id = 1; 442 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); 443 id = 2; 444 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); 445 id = 4; 446 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); 447 id = 3; 448 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); 449 id = 1; 450 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); 451 id = 2; 452 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); 453 id = 4; 454 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); 455 456 /*setup exact solution.*/ 457 ierr = PetscDSSetExactSolution(prob, 0, exactFuncs[0], ctx);CHKERRQ(ierr); 458 ierr = PetscDSSetExactSolution(prob, 1, exactFuncs[1], ctx);CHKERRQ(ierr); 459 ierr = PetscDSSetExactSolution(prob, 2, exactFuncs[2], ctx);CHKERRQ(ierr); 460 ierr = PetscDSSetExactSolutionTimeDerivative(prob, 0, exactFuncs_t[0], ctx);CHKERRQ(ierr); 461 ierr = PetscDSSetExactSolutionTimeDerivative(prob, 1, exactFuncs_t[1], ctx);CHKERRQ(ierr); 462 ierr = PetscDSSetExactSolutionTimeDerivative(prob, 2, exactFuncs_t[2], ctx);CHKERRQ(ierr); 463 PetscFunctionReturn(0); 464 } 465 466 /* x_t = v 467 468 Note that here we use the velocity field at t_{n+1} to advect the particles from 469 t_n to t_{n+1}. If we use both of these fields, we could use Crank-Nicholson or 470 the method of characteristics. 471 */ 472 static PetscErrorCode FreeStreaming(TS ts, PetscReal t, Vec X, Vec F, void *ctx) 473 { 474 AdvCtx *adv = (AdvCtx *) ctx; 475 Vec u = adv->ui; 476 DM sdm, dm, vdm; 477 Vec vel, locvel, pvel; 478 IS vis; 479 DMInterpolationInfo ictx; 480 const PetscScalar *coords, *v; 481 PetscScalar *f; 482 PetscInt vf[1] = {0}; 483 PetscInt dim, Np; 484 PetscErrorCode ierr; 485 486 PetscFunctionBeginUser; 487 ierr = TSGetDM(ts, &sdm);CHKERRQ(ierr); 488 ierr = DMSwarmGetCellDM(sdm, &dm);CHKERRQ(ierr); 489 ierr = DMGetGlobalVector(sdm, &pvel);CHKERRQ(ierr); 490 ierr = DMSwarmGetLocalSize(sdm, &Np);CHKERRQ(ierr); 491 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 492 /* Get local velocity */ 493 ierr = DMCreateSubDM(dm, 1, vf, &vis, &vdm);CHKERRQ(ierr); 494 ierr = VecGetSubVector(u, vis, &vel);CHKERRQ(ierr); 495 ierr = DMGetLocalVector(vdm, &locvel);CHKERRQ(ierr); 496 ierr = DMPlexInsertBoundaryValues(vdm, PETSC_TRUE, locvel, adv->ti, NULL, NULL, NULL);CHKERRQ(ierr); 497 ierr = DMGlobalToLocalBegin(vdm, vel, INSERT_VALUES, locvel);CHKERRQ(ierr); 498 ierr = DMGlobalToLocalEnd(vdm, vel, INSERT_VALUES, locvel);CHKERRQ(ierr); 499 ierr = VecRestoreSubVector(u, vis, &vel);CHKERRQ(ierr); 500 ierr = ISDestroy(&vis);CHKERRQ(ierr); 501 /* Interpolate velocity */ 502 ierr = DMInterpolationCreate(PETSC_COMM_SELF, &ictx);CHKERRQ(ierr); 503 ierr = DMInterpolationSetDim(ictx, dim);CHKERRQ(ierr); 504 ierr = DMInterpolationSetDof(ictx, dim);CHKERRQ(ierr); 505 ierr = VecGetArrayRead(X, &coords);CHKERRQ(ierr); 506 ierr = DMInterpolationAddPoints(ictx, Np, (PetscReal *) coords);CHKERRQ(ierr); 507 ierr = VecRestoreArrayRead(X, &coords);CHKERRQ(ierr); 508 /* Particles that lie outside the domain should be dropped, 509 whereas particles that move to another partition should trigger a migration */ 510 ierr = DMInterpolationSetUp(ictx, vdm, PETSC_FALSE, PETSC_TRUE);CHKERRQ(ierr); 511 ierr = VecSet(pvel, 0.);CHKERRQ(ierr); 512 ierr = DMInterpolationEvaluate(ictx, vdm, locvel, pvel);CHKERRQ(ierr); 513 ierr = DMInterpolationDestroy(&ictx);CHKERRQ(ierr); 514 ierr = DMRestoreLocalVector(vdm, &locvel);CHKERRQ(ierr); 515 ierr = DMDestroy(&vdm);CHKERRQ(ierr); 516 517 ierr = VecGetArray(F, &f);CHKERRQ(ierr); 518 ierr = VecGetArrayRead(pvel, &v);CHKERRQ(ierr); 519 ierr = PetscArraycpy(f, v, Np*dim);CHKERRQ(ierr); 520 ierr = VecRestoreArrayRead(pvel, &v);CHKERRQ(ierr); 521 ierr = VecRestoreArray(F, &f);CHKERRQ(ierr); 522 ierr = DMRestoreGlobalVector(sdm, &pvel);CHKERRQ(ierr); 523 PetscFunctionReturn(0); 524 } 525 526 static PetscErrorCode SetInitialParticleConditions(TS ts, Vec u) 527 { 528 AppCtx *user; 529 void *ctx; 530 DM dm; 531 PetscScalar *coords; 532 PetscReal x[3], dx[3]; 533 PetscInt n[3]; 534 PetscInt Np, dim, d, i, j, k; 535 PetscErrorCode ierr; 536 537 PetscFunctionBegin; 538 ierr = TSGetApplicationContext(ts, &ctx);CHKERRQ(ierr); 539 user = ((AdvCtx *) ctx)->ctx; 540 ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 541 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 542 switch (user->partLayout) { 543 case PART_LAYOUT_CELL: 544 ierr = DMSwarmSetPointCoordinatesRandom(dm, user->Npc);CHKERRQ(ierr); 545 break; 546 case PART_LAYOUT_BOX: 547 Np = 1; 548 for (d = 0; d < dim; ++d) { 549 n[d] = user->Npb; 550 dx[d] = (user->partUpper[d] - user->partLower[d])/PetscMax(1, n[d] - 1); 551 Np *= n[d]; 552 } 553 ierr = VecGetArray(u, &coords);CHKERRQ(ierr); 554 switch (dim) { 555 case 2: 556 x[0] = user->partLower[0]; 557 for (i = 0; i < n[0]; ++i, x[0] += dx[0]) { 558 x[1] = user->partLower[1]; 559 for (j = 0; j < n[1]; ++j, x[1] += dx[1]) { 560 const PetscInt p = j*n[0] + i; 561 for (d = 0; d < dim; ++d) coords[p*dim + d] = x[d]; 562 } 563 } 564 break; 565 case 3: 566 x[0] = user->partLower[0]; 567 for (i = 0; i < n[0]; ++i, x[0] += dx[0]) { 568 x[1] = user->partLower[1]; 569 for (j = 0; j < n[1]; ++j, x[1] += dx[1]) { 570 x[2] = user->partLower[2]; 571 for (k = 0; k < n[2]; ++k, x[2] += dx[2]) { 572 const PetscInt p = (k*n[1] + j)*n[0] + i; 573 for (d = 0; d < dim; ++d) coords[p*dim + d] = x[d]; 574 } 575 } 576 } 577 break; 578 default: SETERRQ1(PetscObjectComm((PetscObject) ts), PETSC_ERR_SUP, "Do not support particle layout in dimension %D", dim); 579 } 580 ierr = VecRestoreArray(u, &coords);CHKERRQ(ierr); 581 break; 582 default: SETERRQ1(PetscObjectComm((PetscObject) ts), PETSC_ERR_ARG_WRONG, "Invalid particle layout type %s", partLayoutTypes[PetscMin(user->partLayout, NUM_PART_LAYOUT_TYPES)]); 583 } 584 PetscFunctionReturn(0); 585 } 586 587 static PetscErrorCode SetupDiscretization(DM dm, DM sdm, AppCtx *user) 588 { 589 DM cdm = dm; 590 PetscFE fe[3]; 591 Parameter *param; 592 PetscInt *cellid, n[3]; 593 PetscReal x[3], dx[3]; 594 PetscScalar *coords; 595 DMPolytopeType ct; 596 PetscInt dim, d, cStart, cEnd, c, Np, p, i, j, k; 597 PetscBool simplex; 598 MPI_Comm comm; 599 PetscErrorCode ierr; 600 601 PetscFunctionBeginUser; 602 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 603 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 604 ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr); 605 simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 606 /* Create finite element */ 607 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 608 ierr = PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr); 609 ierr = PetscObjectSetName((PetscObject) fe[0], "velocity");CHKERRQ(ierr); 610 611 ierr = PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr); 612 ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr); 613 ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr); 614 615 ierr = PetscFECreateDefault(comm, dim, 1, simplex, "temp_", PETSC_DEFAULT, &fe[2]);CHKERRQ(ierr); 616 ierr = PetscFECopyQuadrature(fe[0], fe[2]);CHKERRQ(ierr); 617 ierr = PetscObjectSetName((PetscObject) fe[2], "temperature");CHKERRQ(ierr); 618 619 /* Set discretization and boundary conditions for each mesh */ 620 ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr); 621 ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr); 622 ierr = DMSetField(dm, 2, NULL, (PetscObject) fe[2]);CHKERRQ(ierr); 623 ierr = DMCreateDS(dm);CHKERRQ(ierr); 624 ierr = SetupProblem(dm, user);CHKERRQ(ierr); 625 ierr = PetscBagGetData(user->bag, (void **) ¶m);CHKERRQ(ierr); 626 while (cdm) { 627 ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr); 628 ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 629 } 630 ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr); 631 ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr); 632 ierr = PetscFEDestroy(&fe[2]);CHKERRQ(ierr); 633 634 { 635 PetscObject pressure; 636 MatNullSpace nullspacePres; 637 638 ierr = DMGetField(dm, 1, NULL, &pressure);CHKERRQ(ierr); 639 ierr = MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres);CHKERRQ(ierr); 640 ierr = PetscObjectCompose(pressure, "nullspace", (PetscObject) nullspacePres);CHKERRQ(ierr); 641 ierr = MatNullSpaceDestroy(&nullspacePres);CHKERRQ(ierr); 642 } 643 644 /* Setup particle information */ 645 ierr = DMSwarmSetType(sdm, DMSWARM_PIC);CHKERRQ(ierr); 646 ierr = DMSwarmRegisterPetscDatatypeField(sdm, "mass", 1, PETSC_REAL);CHKERRQ(ierr); 647 ierr = DMSwarmFinalizeFieldRegister(sdm);CHKERRQ(ierr); 648 switch (user->partLayout) { 649 case PART_LAYOUT_CELL: 650 ierr = DMSwarmSetLocalSizes(sdm, (cEnd - cStart) * user->Npc, 0);CHKERRQ(ierr); 651 ierr = DMSetFromOptions(sdm);CHKERRQ(ierr); 652 ierr = DMSwarmGetField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr); 653 for (c = cStart; c < cEnd; ++c) { 654 for (p = 0; p < user->Npc; ++p) { 655 const PetscInt n = c*user->Npc + p; 656 657 cellid[n] = c; 658 } 659 } 660 ierr = DMSwarmRestoreField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr); 661 ierr = DMSwarmSetPointCoordinatesRandom(sdm, user->Npc);CHKERRQ(ierr); 662 break; 663 case PART_LAYOUT_BOX: 664 Np = 1; 665 for (d = 0; d < dim; ++d) { 666 n[d] = user->Npb; 667 dx[d] = (user->partUpper[d] - user->partLower[d])/PetscMax(1, n[d] - 1); 668 Np *= n[d]; 669 } 670 ierr = DMSwarmSetLocalSizes(sdm, Np, 0);CHKERRQ(ierr); 671 ierr = DMSetFromOptions(sdm);CHKERRQ(ierr); 672 ierr = DMSwarmGetField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 673 switch (dim) { 674 case 2: 675 x[0] = user->partLower[0]; 676 for (i = 0; i < n[0]; ++i, x[0] += dx[0]) { 677 x[1] = user->partLower[1]; 678 for (j = 0; j < n[1]; ++j, x[1] += dx[1]) { 679 const PetscInt p = j*n[0] + i; 680 for (d = 0; d < dim; ++d) coords[p*dim + d] = x[d]; 681 } 682 } 683 break; 684 case 3: 685 x[0] = user->partLower[0]; 686 for (i = 0; i < n[0]; ++i, x[0] += dx[0]) { 687 x[1] = user->partLower[1]; 688 for (j = 0; j < n[1]; ++j, x[1] += dx[1]) { 689 x[2] = user->partLower[2]; 690 for (k = 0; k < n[2]; ++k, x[2] += dx[2]) { 691 const PetscInt p = (k*n[1] + j)*n[0] + i; 692 for (d = 0; d < dim; ++d) coords[p*dim + d] = x[d]; 693 } 694 } 695 } 696 break; 697 default: SETERRQ1(comm, PETSC_ERR_SUP, "Do not support particle layout in dimension %D", dim); 698 } 699 ierr = DMSwarmRestoreField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 700 ierr = DMSwarmGetField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr); 701 for (p = 0; p < Np; ++p) cellid[p] = 0; 702 ierr = DMSwarmRestoreField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr); 703 ierr = DMSwarmMigrate(sdm, PETSC_TRUE);CHKERRQ(ierr); 704 break; 705 default: SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Invalid particle layout type %s", partLayoutTypes[PetscMin(user->partLayout, NUM_PART_LAYOUT_TYPES)]); 706 } 707 ierr = PetscObjectSetName((PetscObject) sdm, "Particles");CHKERRQ(ierr); 708 ierr = DMViewFromOptions(sdm, NULL, "-dm_view");CHKERRQ(ierr); 709 PetscFunctionReturn(0); 710 } 711 712 static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace) 713 { 714 Vec vec; 715 PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero}; 716 PetscErrorCode ierr; 717 718 PetscFunctionBeginUser; 719 if (ofield != 1) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Nullspace must be for pressure field at index 1, not %D", ofield); 720 funcs[nfield] = constant; 721 ierr = DMCreateGlobalVector(dm, &vec);CHKERRQ(ierr); 722 ierr = DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec);CHKERRQ(ierr); 723 ierr = VecNormalize(vec, NULL);CHKERRQ(ierr); 724 ierr = PetscObjectSetName((PetscObject) vec, "Pressure Null Space");CHKERRQ(ierr); 725 ierr = VecViewFromOptions(vec, NULL, "-pressure_nullspace_view");CHKERRQ(ierr); 726 ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_FALSE, 1, &vec, nullSpace);CHKERRQ(ierr); 727 ierr = VecDestroy(&vec);CHKERRQ(ierr); 728 PetscFunctionReturn(0); 729 } 730 731 static PetscErrorCode RemoveDiscretePressureNullspace_Private(TS ts, Vec u) 732 { 733 DM dm; 734 MatNullSpace nullsp; 735 PetscErrorCode ierr; 736 737 PetscFunctionBegin; 738 ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 739 ierr = CreatePressureNullSpace(dm, 1, 1, &nullsp);CHKERRQ(ierr); 740 ierr = MatNullSpaceRemove(nullsp, u);CHKERRQ(ierr); 741 ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr); 742 PetscFunctionReturn(0); 743 } 744 745 /* Make the discrete pressure discretely divergence free */ 746 static PetscErrorCode RemoveDiscretePressureNullspace(TS ts) 747 { 748 Vec u; 749 PetscErrorCode ierr; 750 751 PetscFunctionBegin; 752 ierr = TSGetSolution(ts, &u);CHKERRQ(ierr); 753 ierr = RemoveDiscretePressureNullspace_Private(ts, u);CHKERRQ(ierr); 754 PetscFunctionReturn(0); 755 } 756 757 static PetscErrorCode SetInitialConditions(TS ts, Vec u) 758 { 759 DM dm; 760 PetscReal t; 761 PetscErrorCode ierr; 762 763 PetscFunctionBegin; 764 ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 765 ierr = TSGetTime(ts, &t);CHKERRQ(ierr); 766 ierr = DMComputeExactSolution(dm, t, u, NULL);CHKERRQ(ierr); 767 ierr = RemoveDiscretePressureNullspace_Private(ts, u);CHKERRQ(ierr); 768 PetscFunctionReturn(0); 769 } 770 771 static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx) 772 { 773 PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 774 void *ctxs[3]; 775 DM dm; 776 PetscDS ds; 777 Vec v; 778 PetscReal ferrors[3]; 779 PetscInt tl, l, f; 780 PetscErrorCode ierr; 781 782 PetscFunctionBeginUser; 783 ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 784 ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 785 786 for (f = 0; f < 3; ++f) {ierr = PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f]);CHKERRQ(ierr);} 787 ierr = DMComputeL2FieldDiff(dm, crtime, exactFuncs, ctxs, u, ferrors);CHKERRQ(ierr); 788 ierr = PetscObjectGetTabLevel((PetscObject) ts, &tl);CHKERRQ(ierr); 789 for (l = 0; l < tl; ++l) {ierr = PetscPrintf(PETSC_COMM_WORLD, "\t");CHKERRQ(ierr);} 790 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); 791 792 ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr); 793 ierr = PetscObjectSetName((PetscObject) u, "Numerical Solution");CHKERRQ(ierr); 794 ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr); 795 ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr); 796 797 ierr = DMGetGlobalVector(dm, &v);CHKERRQ(ierr); 798 ierr = DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, v);CHKERRQ(ierr); 799 ierr = PetscObjectSetName((PetscObject) v, "Exact Solution");CHKERRQ(ierr); 800 ierr = VecViewFromOptions(v, NULL, "-exact_vec_view");CHKERRQ(ierr); 801 ierr = DMRestoreGlobalVector(dm, &v);CHKERRQ(ierr); 802 803 PetscFunctionReturn(0); 804 } 805 806 /* Note that adv->x0 will not be correct after migration */ 807 static PetscErrorCode ComputeParticleError(TS ts, Vec u, Vec e) 808 { 809 AdvCtx *adv; 810 DM sdm; 811 Parameter *param; 812 const PetscScalar *xp0, *xp; 813 PetscScalar *ep; 814 PetscReal time; 815 PetscInt dim, Np, p; 816 MPI_Comm comm; 817 PetscErrorCode ierr; 818 819 PetscFunctionBeginUser; 820 ierr = TSGetTime(ts, &time);CHKERRQ(ierr); 821 ierr = TSGetApplicationContext(ts, (void **) &adv);CHKERRQ(ierr); 822 ierr = PetscBagGetData(adv->ctx->bag, (void **) ¶m);CHKERRQ(ierr); 823 ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr); 824 ierr = TSGetDM(ts, &sdm);CHKERRQ(ierr); 825 ierr = DMGetDimension(sdm, &dim);CHKERRQ(ierr); 826 ierr = DMSwarmGetLocalSize(sdm, &Np);CHKERRQ(ierr); 827 ierr = VecGetArrayRead(adv->x0, &xp0);CHKERRQ(ierr); 828 ierr = VecGetArrayRead(u, &xp);CHKERRQ(ierr); 829 ierr = VecGetArrayWrite(e, &ep);CHKERRQ(ierr); 830 for (p = 0; p < Np; ++p) { 831 PetscScalar x[3]; 832 PetscReal x0[3]; 833 PetscInt d; 834 835 for (d = 0; d < dim; ++d) x0[d] = PetscRealPart(xp0[p*dim+d]); 836 ierr = adv->exact(dim, time, x0, 1, x, param);CHKERRQ(ierr); 837 for (d = 0; d < dim; ++d) ep[p*dim+d] += x[d] - xp[p*dim+d]; 838 } 839 ierr = VecRestoreArrayRead(adv->x0, &xp0);CHKERRQ(ierr); 840 ierr = VecRestoreArrayRead(u, &xp);CHKERRQ(ierr); 841 ierr = VecRestoreArrayWrite(e, &ep);CHKERRQ(ierr); 842 PetscFunctionReturn(0); 843 } 844 845 static PetscErrorCode MonitorParticleError(TS ts, PetscInt step, PetscReal time, Vec u, void *ctx) 846 { 847 AdvCtx *adv = (AdvCtx *) ctx; 848 DM sdm; 849 Parameter *param; 850 const PetscScalar *xp0, *xp; 851 PetscReal error = 0.0; 852 PetscInt dim, tl, l, Np, p; 853 MPI_Comm comm; 854 PetscErrorCode ierr; 855 856 PetscFunctionBeginUser; 857 ierr = PetscBagGetData(adv->ctx->bag, (void **) ¶m);CHKERRQ(ierr); 858 ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr); 859 ierr = TSGetDM(ts, &sdm);CHKERRQ(ierr); 860 ierr = DMGetDimension(sdm, &dim);CHKERRQ(ierr); 861 ierr = DMSwarmGetLocalSize(sdm, &Np);CHKERRQ(ierr); 862 ierr = VecGetArrayRead(adv->x0, &xp0);CHKERRQ(ierr); 863 ierr = VecGetArrayRead(u, &xp);CHKERRQ(ierr); 864 for (p = 0; p < Np; ++p) { 865 PetscScalar x[3]; 866 PetscReal x0[3]; 867 PetscReal perror = 0.0; 868 PetscInt d; 869 870 for (d = 0; d < dim; ++d) x0[d] = PetscRealPart(xp0[p*dim+d]); 871 ierr = adv->exact(dim, time, x0, 1, x, param);CHKERRQ(ierr); 872 for (d = 0; d < dim; ++d) perror += PetscSqr(PetscRealPart(x[d] - xp[p*dim+d])); 873 error += perror; 874 } 875 ierr = VecRestoreArrayRead(adv->x0, &xp0);CHKERRQ(ierr); 876 ierr = VecRestoreArrayRead(u, &xp);CHKERRQ(ierr); 877 ierr = PetscObjectGetTabLevel((PetscObject) ts, &tl);CHKERRQ(ierr); 878 for (l = 0; l < tl; ++l) {ierr = PetscPrintf(PETSC_COMM_WORLD, "\t");CHKERRQ(ierr);} 879 ierr = PetscPrintf(comm, "Timestep: %04d time = %-8.4g \t L_2 Particle Error: [%2.3g]\n", (int) step, (double) time, (double) error);CHKERRQ(ierr); 880 PetscFunctionReturn(0); 881 } 882 883 static PetscErrorCode AdvectParticles(TS ts) 884 { 885 TS sts; 886 DM sdm; 887 Vec coordinates; 888 AdvCtx *adv; 889 PetscReal time; 890 PetscBool lreset, reset; 891 PetscInt dim, n, N, newn, newN; 892 PetscErrorCode ierr; 893 894 PetscFunctionBeginUser; 895 ierr = PetscObjectQuery((PetscObject) ts, "_SwarmTS", (PetscObject *) &sts);CHKERRQ(ierr); 896 ierr = TSGetDM(sts, &sdm);CHKERRQ(ierr); 897 ierr = TSGetRHSFunction(sts, NULL, NULL, (void **) &adv);CHKERRQ(ierr); 898 ierr = DMGetDimension(sdm, &dim);CHKERRQ(ierr); 899 ierr = DMSwarmGetSize(sdm, &N);CHKERRQ(ierr); 900 ierr = DMSwarmGetLocalSize(sdm, &n);CHKERRQ(ierr); 901 ierr = DMSwarmCreateGlobalVectorFromField(sdm, DMSwarmPICField_coor, &coordinates);CHKERRQ(ierr); 902 ierr = TSGetTime(ts, &time);CHKERRQ(ierr); 903 ierr = TSSetMaxTime(sts, time);CHKERRQ(ierr); 904 adv->tf = time; 905 ierr = TSSolve(sts, coordinates);CHKERRQ(ierr); 906 ierr = DMSwarmDestroyGlobalVectorFromField(sdm, DMSwarmPICField_coor, &coordinates);CHKERRQ(ierr); 907 ierr = VecCopy(adv->uf, adv->ui);CHKERRQ(ierr); 908 adv->ti = adv->tf; 909 910 ierr = DMSwarmMigrate(sdm, PETSC_TRUE);CHKERRQ(ierr); 911 ierr = DMSwarmGetSize(sdm, &newN);CHKERRQ(ierr); 912 ierr = DMSwarmGetLocalSize(sdm, &newn);CHKERRQ(ierr); 913 lreset = (n != newn || N != newN) ? PETSC_TRUE : PETSC_FALSE; 914 ierr = MPI_Allreduce(&lreset, &reset, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject) sts));CHKERRMPI(ierr); 915 if (reset) { 916 ierr = TSReset(sts);CHKERRQ(ierr); 917 ierr = DMSwarmVectorDefineField(sdm, DMSwarmPICField_coor);CHKERRQ(ierr); 918 } 919 ierr = DMViewFromOptions(sdm, NULL, "-dm_view");CHKERRQ(ierr); 920 PetscFunctionReturn(0); 921 } 922 923 int main(int argc, char **argv) 924 { 925 DM dm, sdm; 926 TS ts, sts; 927 Vec u, xtmp; 928 AppCtx user; 929 AdvCtx adv; 930 PetscReal t; 931 PetscInt dim; 932 PetscErrorCode ierr; 933 934 ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 935 ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 936 ierr = PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag);CHKERRQ(ierr); 937 ierr = SetupParameters(&user);CHKERRQ(ierr); 938 ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr); 939 ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 940 ierr = TSSetDM(ts, dm);CHKERRQ(ierr); 941 ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr); 942 /* Discretize chemical species */ 943 ierr = DMCreate(PETSC_COMM_WORLD, &sdm);CHKERRQ(ierr); 944 ierr = PetscObjectSetOptionsPrefix((PetscObject) sdm, "part_");CHKERRQ(ierr); 945 ierr = DMSetType(sdm, DMSWARM);CHKERRQ(ierr); 946 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 947 ierr = DMSetDimension(sdm, dim);CHKERRQ(ierr); 948 ierr = DMSwarmSetCellDM(sdm, dm);CHKERRQ(ierr); 949 /* Setup problem */ 950 ierr = SetupDiscretization(dm, sdm, &user);CHKERRQ(ierr); 951 ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr); 952 953 ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 954 ierr = DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace);CHKERRQ(ierr); 955 956 ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user);CHKERRQ(ierr); 957 ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user);CHKERRQ(ierr); 958 ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user);CHKERRQ(ierr); 959 ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 960 ierr = TSSetPreStep(ts, RemoveDiscretePressureNullspace);CHKERRQ(ierr); 961 ierr = TSMonitorSet(ts, MonitorError, &user, NULL);CHKERRQ(ierr);CHKERRQ(ierr); 962 ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 963 964 ierr = TSSetComputeInitialCondition(ts, SetInitialConditions);CHKERRQ(ierr); /* Must come after SetFromOptions() */ 965 ierr = SetInitialConditions(ts, u);CHKERRQ(ierr); 966 ierr = TSGetTime(ts, &t);CHKERRQ(ierr); 967 ierr = DMSetOutputSequenceNumber(dm, 0, t);CHKERRQ(ierr); 968 ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr); 969 970 /* Setup particle position integrator */ 971 ierr = TSCreate(PETSC_COMM_WORLD, &sts);CHKERRQ(ierr); 972 ierr = PetscObjectSetOptionsPrefix((PetscObject) sts, "part_");CHKERRQ(ierr); 973 ierr = PetscObjectIncrementTabLevel((PetscObject) sts, (PetscObject) ts, 1);CHKERRQ(ierr); 974 ierr = TSSetDM(sts, sdm);CHKERRQ(ierr); 975 ierr = TSSetProblemType(sts, TS_NONLINEAR);CHKERRQ(ierr); 976 ierr = TSSetExactFinalTime(sts, TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 977 ierr = TSMonitorSet(sts, MonitorParticleError, &adv, NULL);CHKERRQ(ierr);CHKERRQ(ierr); 978 ierr = TSSetFromOptions(sts);CHKERRQ(ierr); 979 ierr = TSSetApplicationContext(sts, &adv);CHKERRQ(ierr); 980 ierr = TSSetComputeExactError(sts, ComputeParticleError);CHKERRQ(ierr); 981 ierr = TSSetComputeInitialCondition(sts, SetInitialParticleConditions);CHKERRQ(ierr); 982 adv.ti = t; 983 adv.uf = u; 984 ierr = VecDuplicate(adv.uf, &adv.ui); 985 ierr = VecCopy(u, adv.ui);CHKERRQ(ierr); 986 ierr = TSSetRHSFunction(sts, NULL, FreeStreaming, &adv);CHKERRQ(ierr); 987 ierr = TSSetPostStep(ts, AdvectParticles);CHKERRQ(ierr); 988 ierr = PetscObjectCompose((PetscObject) ts, "_SwarmTS", (PetscObject) sts);CHKERRQ(ierr); 989 ierr = DMSwarmVectorDefineField(sdm, DMSwarmPICField_coor);CHKERRQ(ierr); 990 ierr = DMCreateGlobalVector(sdm, &adv.x0);CHKERRQ(ierr); 991 ierr = DMSwarmCreateGlobalVectorFromField(sdm, DMSwarmPICField_coor, &xtmp);CHKERRQ(ierr); 992 ierr = VecCopy(xtmp, adv.x0);CHKERRQ(ierr); 993 ierr = DMSwarmDestroyGlobalVectorFromField(sdm, DMSwarmPICField_coor, &xtmp);CHKERRQ(ierr); 994 switch(user.solType){ 995 case SOL_TRIG_TRIG: adv.exact = trig_trig_x;break; 996 default: SETERRQ2(PetscObjectComm((PetscObject) sdm), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user.solType, NUM_SOL_TYPES)], user.solType); 997 } 998 adv.ctx = &user; 999 1000 ierr = TSSolve(ts, u);CHKERRQ(ierr); 1001 ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr); 1002 ierr = PetscObjectSetName((PetscObject) u, "Numerical Solution");CHKERRQ(ierr); 1003 1004 ierr = VecDestroy(&u);CHKERRQ(ierr); 1005 ierr = VecDestroy(&adv.x0);CHKERRQ(ierr); 1006 ierr = VecDestroy(&adv.ui);CHKERRQ(ierr); 1007 ierr = DMDestroy(&dm);CHKERRQ(ierr); 1008 ierr = DMDestroy(&sdm);CHKERRQ(ierr); 1009 ierr = TSDestroy(&ts);CHKERRQ(ierr); 1010 ierr = TSDestroy(&sts);CHKERRQ(ierr); 1011 ierr = PetscBagDestroy(&user.bag);CHKERRQ(ierr); 1012 ierr = PetscFinalize(); 1013 return ierr; 1014 } 1015 1016 /*TEST 1017 1018 # Swarm does not work with complex 1019 test: 1020 suffix: 2d_tri_p2_p1_p1_tconvp 1021 requires: triangle !single !complex 1022 args: -dm_plex_separate_marker -sol_type trig_trig -dm_refine 2 \ 1023 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 1024 -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1 -ts_monitor_cancel \ 1025 -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 1026 -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 1027 -fieldsplit_0_pc_type lu \ 1028 -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi \ 1029 -omega 0.5 -part_layout_type box -part_lower 0.25,0.25 -part_upper 0.75,0.75 -Npb 5 \ 1030 -part_ts_max_steps 2 -part_ts_dt 0.05 -part_ts_convergence_estimate -convest_num_refine 1 -part_ts_monitor_cancel 1031 test: 1032 suffix: 2d_tri_p2_p1_p1_exit 1033 requires: triangle !single !complex 1034 args: -dm_plex_separate_marker -sol_type trig_trig -dm_refine 1 \ 1035 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 1036 -dmts_check .001 -ts_max_steps 10 -ts_dt 0.1 \ 1037 -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 1038 -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 1039 -fieldsplit_0_pc_type lu \ 1040 -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi \ 1041 -omega 0.5 -part_layout_type box -part_lower 0.25,0.25 -part_upper 0.75,0.75 -Npb 5 \ 1042 -part_ts_max_steps 20 -part_ts_dt 0.05 1043 1044 TEST*/ 1045