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 static void g2_vp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 268 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 269 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 270 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 271 { 272 PetscInt d; 273 for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; 274 } 275 276 static void g3_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 277 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 278 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 279 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 280 { 281 const PetscReal nu = PetscRealPart(constants[0]); 282 const PetscInt Nc = dim; 283 PetscInt c, d; 284 285 for (c = 0; c < Nc; ++c) { 286 for (d = 0; d < dim; ++d) { 287 g3[((c*Nc+c)*dim+d)*dim+d] += nu; // gradU 288 g3[((c*Nc+d)*dim+d)*dim+c] += nu; // gradU transpose 289 } 290 } 291 } 292 293 static void g0_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 294 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 295 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 296 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 297 { 298 PetscInt d; 299 for (d = 0; d < dim; ++d) g0[d] = u_tShift; 300 } 301 302 static void g0_wu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 303 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 304 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 305 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 306 { 307 PetscInt d; 308 for (d = 0; d < dim; ++d) g0[d] = u_x[uOff_x[2]+d]; 309 } 310 311 static void g1_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 312 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 313 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 314 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 315 { 316 PetscInt d; 317 for (d = 0; d < dim; ++d) g1[d] = u[uOff[0]+d]; 318 } 319 320 static void g3_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, 321 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 322 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 323 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 324 { 325 const PetscReal alpha = PetscRealPart(constants[1]); 326 PetscInt d; 327 328 for (d = 0; d < dim; ++d) g3[d*dim+d] = alpha; 329 } 330 331 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 332 { 333 PetscInt sol, pl, n; 334 PetscErrorCode ierr; 335 336 PetscFunctionBeginUser; 337 options->solType = SOL_TRIG_TRIG; 338 options->partLayout = PART_LAYOUT_CELL; 339 options->Npc = 1; 340 options->Npb = 1; 341 342 options->partLower[0] = options->partLower[1] = options->partLower[2] = 0.; 343 options->partUpper[0] = options->partUpper[1] = options->partUpper[2] = 1.; 344 ierr = PetscOptionsBegin(comm, "", "Low Mach flow Problem Options", "DMPLEX");PetscCall(ierr); 345 sol = options->solType; 346 PetscCall(PetscOptionsEList("-sol_type", "The solution type", "ex77.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL)); 347 options->solType = (SolType) sol; 348 pl = options->partLayout; 349 PetscCall(PetscOptionsEList("-part_layout_type", "The particle layout type", "ex77.c", partLayoutTypes, NUM_PART_LAYOUT_TYPES, partLayoutTypes[options->partLayout], &pl, NULL)); 350 options->partLayout = (PartLayoutType) pl; 351 PetscCall(PetscOptionsInt("-Npc", "The initial number of particles per cell", "ex77.c", options->Npc, &options->Npc, NULL)); 352 n = 3; 353 PetscCall(PetscOptionsRealArray("-part_lower", "The lower left corner of the particle box", "ex77.c", options->partLower, &n, NULL)); 354 n = 3; 355 PetscCall(PetscOptionsRealArray("-part_upper", "The upper right corner of the particle box", "ex77.c", options->partUpper, &n, NULL)); 356 PetscCall(PetscOptionsInt("-Npb", "The initial number of particles per box dimension", "ex77.c", options->Npb, &options->Npb, NULL)); 357 ierr = PetscOptionsEnd();PetscCall(ierr); 358 PetscFunctionReturn(0); 359 } 360 361 static PetscErrorCode SetupParameters(AppCtx *user) 362 { 363 PetscBag bag; 364 Parameter *p; 365 366 PetscFunctionBeginUser; 367 /* setup PETSc parameter bag */ 368 PetscCall(PetscBagGetData(user->bag, (void **) &p)); 369 PetscCall(PetscBagSetName(user->bag, "par", "Low Mach flow parameters")); 370 bag = user->bag; 371 PetscCall(PetscBagRegisterReal(bag, &p->nu, 1.0, "nu", "Kinematic viscosity")); 372 PetscCall(PetscBagRegisterReal(bag, &p->alpha, 1.0, "alpha", "Thermal diffusivity")); 373 PetscCall(PetscBagRegisterReal(bag, &p->T_in, 1.0, "T_in", "Inlet temperature")); 374 PetscCall(PetscBagRegisterReal(bag, &p->omega, 1.0, "omega", "Rotation speed in MMS benchmark")); 375 PetscFunctionReturn(0); 376 } 377 378 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 379 { 380 PetscFunctionBeginUser; 381 PetscCall(DMCreate(comm, dm)); 382 PetscCall(DMSetType(*dm, DMPLEX)); 383 PetscCall(DMSetFromOptions(*dm)); 384 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 385 PetscFunctionReturn(0); 386 } 387 388 static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 389 { 390 PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 391 PetscErrorCode (*exactFuncs_t[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 392 PetscDS prob; 393 DMLabel label; 394 Parameter *ctx; 395 PetscInt id; 396 397 PetscFunctionBeginUser; 398 PetscCall(DMGetLabel(dm, "marker", &label)); 399 PetscCall(DMGetDS(dm, &prob)); 400 switch(user->solType) { 401 case SOL_TRIG_TRIG: 402 PetscCall(PetscDSSetResidual(prob, 0, f0_trig_trig_v, f1_v)); 403 PetscCall(PetscDSSetResidual(prob, 2, f0_trig_trig_w, f1_w)); 404 405 exactFuncs[0] = trig_trig_u; 406 exactFuncs[1] = trig_trig_p; 407 exactFuncs[2] = trig_trig_T; 408 exactFuncs_t[0] = trig_trig_u_t; 409 exactFuncs_t[1] = NULL; 410 exactFuncs_t[2] = trig_trig_T_t; 411 break; 412 default: SETERRQ(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType); 413 } 414 415 PetscCall(PetscDSSetResidual(prob, 1, f0_q, NULL)); 416 417 PetscCall(PetscDSSetJacobian(prob, 0, 0, g0_vu, g1_vu, NULL, g3_vu)); 418 PetscCall(PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_vp, NULL)); 419 PetscCall(PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL, NULL)); 420 PetscCall(PetscDSSetJacobian(prob, 2, 0, g0_wu, NULL, NULL, NULL)); 421 PetscCall(PetscDSSetJacobian(prob, 2, 2, g0_wT, g1_wT, NULL, g3_wT)); 422 /* Setup constants */ 423 { 424 Parameter *param; 425 PetscScalar constants[4]; 426 427 PetscCall(PetscBagGetData(user->bag, (void **) ¶m)); 428 429 constants[0] = param->nu; 430 constants[1] = param->alpha; 431 constants[2] = param->T_in; 432 constants[3] = param->omega; 433 PetscCall(PetscDSSetConstants(prob, 4, constants)); 434 } 435 /* Setup Boundary Conditions */ 436 PetscCall(PetscBagGetData(user->bag, (void **) &ctx)); 437 id = 3; 438 PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], (void (*)(void)) exactFuncs_t[0], ctx, NULL)); 439 id = 1; 440 PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], (void (*)(void)) exactFuncs_t[0], ctx, NULL)); 441 id = 2; 442 PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], (void (*)(void)) exactFuncs_t[0], ctx, NULL)); 443 id = 4; 444 PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], (void (*)(void)) exactFuncs_t[0], ctx, NULL)); 445 id = 3; 446 PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], (void (*)(void)) exactFuncs_t[2], ctx, NULL)); 447 id = 1; 448 PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], (void (*)(void)) exactFuncs_t[2], ctx, NULL)); 449 id = 2; 450 PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], (void (*)(void)) exactFuncs_t[2], ctx, NULL)); 451 id = 4; 452 PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], (void (*)(void)) exactFuncs_t[2], ctx, NULL)); 453 454 /*setup exact solution.*/ 455 PetscCall(PetscDSSetExactSolution(prob, 0, exactFuncs[0], ctx)); 456 PetscCall(PetscDSSetExactSolution(prob, 1, exactFuncs[1], ctx)); 457 PetscCall(PetscDSSetExactSolution(prob, 2, exactFuncs[2], ctx)); 458 PetscCall(PetscDSSetExactSolutionTimeDerivative(prob, 0, exactFuncs_t[0], ctx)); 459 PetscCall(PetscDSSetExactSolutionTimeDerivative(prob, 1, exactFuncs_t[1], ctx)); 460 PetscCall(PetscDSSetExactSolutionTimeDerivative(prob, 2, exactFuncs_t[2], ctx)); 461 PetscFunctionReturn(0); 462 } 463 464 /* x_t = v 465 466 Note that here we use the velocity field at t_{n+1} to advect the particles from 467 t_n to t_{n+1}. If we use both of these fields, we could use Crank-Nicholson or 468 the method of characteristics. 469 */ 470 static PetscErrorCode FreeStreaming(TS ts, PetscReal t, Vec X, Vec F, void *ctx) 471 { 472 AdvCtx *adv = (AdvCtx *) ctx; 473 Vec u = adv->ui; 474 DM sdm, dm, vdm; 475 Vec vel, locvel, pvel; 476 IS vis; 477 DMInterpolationInfo ictx; 478 const PetscScalar *coords, *v; 479 PetscScalar *f; 480 PetscInt vf[1] = {0}; 481 PetscInt dim, Np; 482 483 PetscFunctionBeginUser; 484 PetscCall(TSGetDM(ts, &sdm)); 485 PetscCall(DMSwarmGetCellDM(sdm, &dm)); 486 PetscCall(DMGetGlobalVector(sdm, &pvel)); 487 PetscCall(DMSwarmGetLocalSize(sdm, &Np)); 488 PetscCall(DMGetDimension(dm, &dim)); 489 /* Get local velocity */ 490 PetscCall(DMCreateSubDM(dm, 1, vf, &vis, &vdm)); 491 PetscCall(VecGetSubVector(u, vis, &vel)); 492 PetscCall(DMGetLocalVector(vdm, &locvel)); 493 PetscCall(DMPlexInsertBoundaryValues(vdm, PETSC_TRUE, locvel, adv->ti, NULL, NULL, NULL)); 494 PetscCall(DMGlobalToLocalBegin(vdm, vel, INSERT_VALUES, locvel)); 495 PetscCall(DMGlobalToLocalEnd(vdm, vel, INSERT_VALUES, locvel)); 496 PetscCall(VecRestoreSubVector(u, vis, &vel)); 497 PetscCall(ISDestroy(&vis)); 498 /* Interpolate velocity */ 499 PetscCall(DMInterpolationCreate(PETSC_COMM_SELF, &ictx)); 500 PetscCall(DMInterpolationSetDim(ictx, dim)); 501 PetscCall(DMInterpolationSetDof(ictx, dim)); 502 PetscCall(VecGetArrayRead(X, &coords)); 503 PetscCall(DMInterpolationAddPoints(ictx, Np, (PetscReal *) coords)); 504 PetscCall(VecRestoreArrayRead(X, &coords)); 505 /* Particles that lie outside the domain should be dropped, 506 whereas particles that move to another partition should trigger a migration */ 507 PetscCall(DMInterpolationSetUp(ictx, vdm, PETSC_FALSE, PETSC_TRUE)); 508 PetscCall(VecSet(pvel, 0.)); 509 PetscCall(DMInterpolationEvaluate(ictx, vdm, locvel, pvel)); 510 PetscCall(DMInterpolationDestroy(&ictx)); 511 PetscCall(DMRestoreLocalVector(vdm, &locvel)); 512 PetscCall(DMDestroy(&vdm)); 513 514 PetscCall(VecGetArray(F, &f)); 515 PetscCall(VecGetArrayRead(pvel, &v)); 516 PetscCall(PetscArraycpy(f, v, Np*dim)); 517 PetscCall(VecRestoreArrayRead(pvel, &v)); 518 PetscCall(VecRestoreArray(F, &f)); 519 PetscCall(DMRestoreGlobalVector(sdm, &pvel)); 520 PetscFunctionReturn(0); 521 } 522 523 static PetscErrorCode SetInitialParticleConditions(TS ts, Vec u) 524 { 525 AppCtx *user; 526 void *ctx; 527 DM dm; 528 PetscScalar *coords; 529 PetscReal x[3], dx[3]; 530 PetscInt n[3]; 531 PetscInt dim, d, i, j, k; 532 533 PetscFunctionBegin; 534 PetscCall(TSGetApplicationContext(ts, &ctx)); 535 user = ((AdvCtx *) ctx)->ctx; 536 PetscCall(TSGetDM(ts, &dm)); 537 PetscCall(DMGetDimension(dm, &dim)); 538 switch (user->partLayout) { 539 case PART_LAYOUT_CELL: 540 PetscCall(DMSwarmSetPointCoordinatesRandom(dm, user->Npc)); 541 break; 542 case PART_LAYOUT_BOX: 543 for (d = 0; d < dim; ++d) { 544 n[d] = user->Npb; 545 dx[d] = (user->partUpper[d] - user->partLower[d])/PetscMax(1, n[d] - 1); 546 } 547 PetscCall(VecGetArray(u, &coords)); 548 switch (dim) { 549 case 2: 550 x[0] = user->partLower[0]; 551 for (i = 0; i < n[0]; ++i, x[0] += dx[0]) { 552 x[1] = user->partLower[1]; 553 for (j = 0; j < n[1]; ++j, x[1] += dx[1]) { 554 const PetscInt p = j*n[0] + i; 555 for (d = 0; d < dim; ++d) coords[p*dim + d] = x[d]; 556 } 557 } 558 break; 559 case 3: 560 x[0] = user->partLower[0]; 561 for (i = 0; i < n[0]; ++i, x[0] += dx[0]) { 562 x[1] = user->partLower[1]; 563 for (j = 0; j < n[1]; ++j, x[1] += dx[1]) { 564 x[2] = user->partLower[2]; 565 for (k = 0; k < n[2]; ++k, x[2] += dx[2]) { 566 const PetscInt p = (k*n[1] + j)*n[0] + i; 567 for (d = 0; d < dim; ++d) coords[p*dim + d] = x[d]; 568 } 569 } 570 } 571 break; 572 default: SETERRQ(PetscObjectComm((PetscObject) ts), PETSC_ERR_SUP, "Do not support particle layout in dimension %D", dim); 573 } 574 PetscCall(VecRestoreArray(u, &coords)); 575 break; 576 default: SETERRQ(PetscObjectComm((PetscObject) ts), PETSC_ERR_ARG_WRONG, "Invalid particle layout type %s", partLayoutTypes[PetscMin(user->partLayout, NUM_PART_LAYOUT_TYPES)]); 577 } 578 PetscFunctionReturn(0); 579 } 580 581 static PetscErrorCode SetupDiscretization(DM dm, DM sdm, AppCtx *user) 582 { 583 DM cdm = dm; 584 PetscFE fe[3]; 585 Parameter *param; 586 PetscInt *cellid, n[3]; 587 PetscReal x[3], dx[3]; 588 PetscScalar *coords; 589 DMPolytopeType ct; 590 PetscInt dim, d, cStart, cEnd, c, Np, p, i, j, k; 591 PetscBool simplex; 592 MPI_Comm comm; 593 594 PetscFunctionBeginUser; 595 PetscCall(DMGetDimension(dm, &dim)); 596 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 597 PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 598 simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 599 /* Create finite element */ 600 PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 601 PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0])); 602 PetscCall(PetscObjectSetName((PetscObject) fe[0], "velocity")); 603 604 PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1])); 605 PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 606 PetscCall(PetscObjectSetName((PetscObject) fe[1], "pressure")); 607 608 PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "temp_", PETSC_DEFAULT, &fe[2])); 609 PetscCall(PetscFECopyQuadrature(fe[0], fe[2])); 610 PetscCall(PetscObjectSetName((PetscObject) fe[2], "temperature")); 611 612 /* Set discretization and boundary conditions for each mesh */ 613 PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe[0])); 614 PetscCall(DMSetField(dm, 1, NULL, (PetscObject) fe[1])); 615 PetscCall(DMSetField(dm, 2, NULL, (PetscObject) fe[2])); 616 PetscCall(DMCreateDS(dm)); 617 PetscCall(SetupProblem(dm, user)); 618 PetscCall(PetscBagGetData(user->bag, (void **) ¶m)); 619 while (cdm) { 620 PetscCall(DMCopyDisc(dm, cdm)); 621 PetscCall(DMGetCoarseDM(cdm, &cdm)); 622 } 623 PetscCall(PetscFEDestroy(&fe[0])); 624 PetscCall(PetscFEDestroy(&fe[1])); 625 PetscCall(PetscFEDestroy(&fe[2])); 626 627 { 628 PetscObject pressure; 629 MatNullSpace nullspacePres; 630 631 PetscCall(DMGetField(dm, 1, NULL, &pressure)); 632 PetscCall(MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres)); 633 PetscCall(PetscObjectCompose(pressure, "nullspace", (PetscObject) nullspacePres)); 634 PetscCall(MatNullSpaceDestroy(&nullspacePres)); 635 } 636 637 /* Setup particle information */ 638 PetscCall(DMSwarmSetType(sdm, DMSWARM_PIC)); 639 PetscCall(DMSwarmRegisterPetscDatatypeField(sdm, "mass", 1, PETSC_REAL)); 640 PetscCall(DMSwarmFinalizeFieldRegister(sdm)); 641 switch (user->partLayout) { 642 case PART_LAYOUT_CELL: 643 PetscCall(DMSwarmSetLocalSizes(sdm, (cEnd - cStart) * user->Npc, 0)); 644 PetscCall(DMSetFromOptions(sdm)); 645 PetscCall(DMSwarmGetField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid)); 646 for (c = cStart; c < cEnd; ++c) { 647 for (p = 0; p < user->Npc; ++p) { 648 const PetscInt n = c*user->Npc + p; 649 650 cellid[n] = c; 651 } 652 } 653 PetscCall(DMSwarmRestoreField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid)); 654 PetscCall(DMSwarmSetPointCoordinatesRandom(sdm, user->Npc)); 655 break; 656 case PART_LAYOUT_BOX: 657 Np = 1; 658 for (d = 0; d < dim; ++d) { 659 n[d] = user->Npb; 660 dx[d] = (user->partUpper[d] - user->partLower[d])/PetscMax(1, n[d] - 1); 661 Np *= n[d]; 662 } 663 PetscCall(DMSwarmSetLocalSizes(sdm, Np, 0)); 664 PetscCall(DMSetFromOptions(sdm)); 665 PetscCall(DMSwarmGetField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords)); 666 switch (dim) { 667 case 2: 668 x[0] = user->partLower[0]; 669 for (i = 0; i < n[0]; ++i, x[0] += dx[0]) { 670 x[1] = user->partLower[1]; 671 for (j = 0; j < n[1]; ++j, x[1] += dx[1]) { 672 const PetscInt p = j*n[0] + i; 673 for (d = 0; d < dim; ++d) coords[p*dim + d] = x[d]; 674 } 675 } 676 break; 677 case 3: 678 x[0] = user->partLower[0]; 679 for (i = 0; i < n[0]; ++i, x[0] += dx[0]) { 680 x[1] = user->partLower[1]; 681 for (j = 0; j < n[1]; ++j, x[1] += dx[1]) { 682 x[2] = user->partLower[2]; 683 for (k = 0; k < n[2]; ++k, x[2] += dx[2]) { 684 const PetscInt p = (k*n[1] + j)*n[0] + i; 685 for (d = 0; d < dim; ++d) coords[p*dim + d] = x[d]; 686 } 687 } 688 } 689 break; 690 default: SETERRQ(comm, PETSC_ERR_SUP, "Do not support particle layout in dimension %D", dim); 691 } 692 PetscCall(DMSwarmRestoreField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords)); 693 PetscCall(DMSwarmGetField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid)); 694 for (p = 0; p < Np; ++p) cellid[p] = 0; 695 PetscCall(DMSwarmRestoreField(sdm, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid)); 696 PetscCall(DMSwarmMigrate(sdm, PETSC_TRUE)); 697 break; 698 default: SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid particle layout type %s", partLayoutTypes[PetscMin(user->partLayout, NUM_PART_LAYOUT_TYPES)]); 699 } 700 PetscCall(PetscObjectSetName((PetscObject) sdm, "Particles")); 701 PetscCall(DMViewFromOptions(sdm, NULL, "-dm_view")); 702 PetscFunctionReturn(0); 703 } 704 705 static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace) 706 { 707 Vec vec; 708 PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero}; 709 710 PetscFunctionBeginUser; 711 PetscCheck(ofield == 1,PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Nullspace must be for pressure field at index 1, not %D", ofield); 712 funcs[nfield] = constant; 713 PetscCall(DMCreateGlobalVector(dm, &vec)); 714 PetscCall(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec)); 715 PetscCall(VecNormalize(vec, NULL)); 716 PetscCall(PetscObjectSetName((PetscObject) vec, "Pressure Null Space")); 717 PetscCall(VecViewFromOptions(vec, NULL, "-pressure_nullspace_view")); 718 PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_FALSE, 1, &vec, nullSpace)); 719 PetscCall(VecDestroy(&vec)); 720 PetscFunctionReturn(0); 721 } 722 723 static PetscErrorCode RemoveDiscretePressureNullspace_Private(TS ts, Vec u) 724 { 725 DM dm; 726 MatNullSpace nullsp; 727 728 PetscFunctionBegin; 729 PetscCall(TSGetDM(ts, &dm)); 730 PetscCall(CreatePressureNullSpace(dm, 1, 1, &nullsp)); 731 PetscCall(MatNullSpaceRemove(nullsp, u)); 732 PetscCall(MatNullSpaceDestroy(&nullsp)); 733 PetscFunctionReturn(0); 734 } 735 736 /* Make the discrete pressure discretely divergence free */ 737 static PetscErrorCode RemoveDiscretePressureNullspace(TS ts) 738 { 739 Vec u; 740 741 PetscFunctionBegin; 742 PetscCall(TSGetSolution(ts, &u)); 743 PetscCall(RemoveDiscretePressureNullspace_Private(ts, u)); 744 PetscFunctionReturn(0); 745 } 746 747 static PetscErrorCode SetInitialConditions(TS ts, Vec u) 748 { 749 DM dm; 750 PetscReal t; 751 752 PetscFunctionBegin; 753 PetscCall(TSGetDM(ts, &dm)); 754 PetscCall(TSGetTime(ts, &t)); 755 PetscCall(DMComputeExactSolution(dm, t, u, NULL)); 756 PetscCall(RemoveDiscretePressureNullspace_Private(ts, u)); 757 PetscFunctionReturn(0); 758 } 759 760 static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx) 761 { 762 PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 763 void *ctxs[3]; 764 DM dm; 765 PetscDS ds; 766 Vec v; 767 PetscReal ferrors[3]; 768 PetscInt tl, l, f; 769 770 PetscFunctionBeginUser; 771 PetscCall(TSGetDM(ts, &dm)); 772 PetscCall(DMGetDS(dm, &ds)); 773 774 for (f = 0; f < 3; ++f) PetscCall(PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f])); 775 PetscCall(DMComputeL2FieldDiff(dm, crtime, exactFuncs, ctxs, u, ferrors)); 776 PetscCall(PetscObjectGetTabLevel((PetscObject) ts, &tl)); 777 for (l = 0; l < tl; ++l) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\t")); 778 PetscCall(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])); 779 780 PetscCall(DMGetGlobalVector(dm, &u)); 781 PetscCall(PetscObjectSetName((PetscObject) u, "Numerical Solution")); 782 PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 783 PetscCall(DMRestoreGlobalVector(dm, &u)); 784 785 PetscCall(DMGetGlobalVector(dm, &v)); 786 PetscCall(DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, v)); 787 PetscCall(PetscObjectSetName((PetscObject) v, "Exact Solution")); 788 PetscCall(VecViewFromOptions(v, NULL, "-exact_vec_view")); 789 PetscCall(DMRestoreGlobalVector(dm, &v)); 790 791 PetscFunctionReturn(0); 792 } 793 794 /* Note that adv->x0 will not be correct after migration */ 795 static PetscErrorCode ComputeParticleError(TS ts, Vec u, Vec e) 796 { 797 AdvCtx *adv; 798 DM sdm; 799 Parameter *param; 800 const PetscScalar *xp0, *xp; 801 PetscScalar *ep; 802 PetscReal time; 803 PetscInt dim, Np, p; 804 MPI_Comm comm; 805 806 PetscFunctionBeginUser; 807 PetscCall(TSGetTime(ts, &time)); 808 PetscCall(TSGetApplicationContext(ts, &adv)); 809 PetscCall(PetscBagGetData(adv->ctx->bag, (void **) ¶m)); 810 PetscCall(PetscObjectGetComm((PetscObject) ts, &comm)); 811 PetscCall(TSGetDM(ts, &sdm)); 812 PetscCall(DMGetDimension(sdm, &dim)); 813 PetscCall(DMSwarmGetLocalSize(sdm, &Np)); 814 PetscCall(VecGetArrayRead(adv->x0, &xp0)); 815 PetscCall(VecGetArrayRead(u, &xp)); 816 PetscCall(VecGetArrayWrite(e, &ep)); 817 for (p = 0; p < Np; ++p) { 818 PetscScalar x[3]; 819 PetscReal x0[3]; 820 PetscInt d; 821 822 for (d = 0; d < dim; ++d) x0[d] = PetscRealPart(xp0[p*dim+d]); 823 PetscCall(adv->exact(dim, time, x0, 1, x, param)); 824 for (d = 0; d < dim; ++d) ep[p*dim+d] += x[d] - xp[p*dim+d]; 825 } 826 PetscCall(VecRestoreArrayRead(adv->x0, &xp0)); 827 PetscCall(VecRestoreArrayRead(u, &xp)); 828 PetscCall(VecRestoreArrayWrite(e, &ep)); 829 PetscFunctionReturn(0); 830 } 831 832 static PetscErrorCode MonitorParticleError(TS ts, PetscInt step, PetscReal time, Vec u, void *ctx) 833 { 834 AdvCtx *adv = (AdvCtx *) ctx; 835 DM sdm; 836 Parameter *param; 837 const PetscScalar *xp0, *xp; 838 PetscReal error = 0.0; 839 PetscInt dim, tl, l, Np, p; 840 MPI_Comm comm; 841 842 PetscFunctionBeginUser; 843 PetscCall(PetscBagGetData(adv->ctx->bag, (void **) ¶m)); 844 PetscCall(PetscObjectGetComm((PetscObject) ts, &comm)); 845 PetscCall(TSGetDM(ts, &sdm)); 846 PetscCall(DMGetDimension(sdm, &dim)); 847 PetscCall(DMSwarmGetLocalSize(sdm, &Np)); 848 PetscCall(VecGetArrayRead(adv->x0, &xp0)); 849 PetscCall(VecGetArrayRead(u, &xp)); 850 for (p = 0; p < Np; ++p) { 851 PetscScalar x[3]; 852 PetscReal x0[3]; 853 PetscReal perror = 0.0; 854 PetscInt d; 855 856 for (d = 0; d < dim; ++d) x0[d] = PetscRealPart(xp0[p*dim+d]); 857 PetscCall(adv->exact(dim, time, x0, 1, x, param)); 858 for (d = 0; d < dim; ++d) perror += PetscSqr(PetscRealPart(x[d] - xp[p*dim+d])); 859 error += perror; 860 } 861 PetscCall(VecRestoreArrayRead(adv->x0, &xp0)); 862 PetscCall(VecRestoreArrayRead(u, &xp)); 863 PetscCall(PetscObjectGetTabLevel((PetscObject) ts, &tl)); 864 for (l = 0; l < tl; ++l) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\t")); 865 PetscCall(PetscPrintf(comm, "Timestep: %04d time = %-8.4g \t L_2 Particle Error: [%2.3g]\n", (int) step, (double) time, (double) error)); 866 PetscFunctionReturn(0); 867 } 868 869 static PetscErrorCode AdvectParticles(TS ts) 870 { 871 TS sts; 872 DM sdm; 873 Vec coordinates; 874 AdvCtx *adv; 875 PetscReal time; 876 PetscBool lreset, reset; 877 PetscInt dim, n, N, newn, newN; 878 879 PetscFunctionBeginUser; 880 PetscCall(PetscObjectQuery((PetscObject) ts, "_SwarmTS", (PetscObject *) &sts)); 881 PetscCall(TSGetDM(sts, &sdm)); 882 PetscCall(TSGetRHSFunction(sts, NULL, NULL, (void **) &adv)); 883 PetscCall(DMGetDimension(sdm, &dim)); 884 PetscCall(DMSwarmGetSize(sdm, &N)); 885 PetscCall(DMSwarmGetLocalSize(sdm, &n)); 886 PetscCall(DMSwarmCreateGlobalVectorFromField(sdm, DMSwarmPICField_coor, &coordinates)); 887 PetscCall(TSGetTime(ts, &time)); 888 PetscCall(TSSetMaxTime(sts, time)); 889 adv->tf = time; 890 PetscCall(TSSolve(sts, coordinates)); 891 PetscCall(DMSwarmDestroyGlobalVectorFromField(sdm, DMSwarmPICField_coor, &coordinates)); 892 PetscCall(VecCopy(adv->uf, adv->ui)); 893 adv->ti = adv->tf; 894 895 PetscCall(DMSwarmMigrate(sdm, PETSC_TRUE)); 896 PetscCall(DMSwarmGetSize(sdm, &newN)); 897 PetscCall(DMSwarmGetLocalSize(sdm, &newn)); 898 lreset = (n != newn || N != newN) ? PETSC_TRUE : PETSC_FALSE; 899 PetscCallMPI(MPI_Allreduce(&lreset, &reset, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject) sts))); 900 if (reset) { 901 PetscCall(TSReset(sts)); 902 PetscCall(DMSwarmVectorDefineField(sdm, DMSwarmPICField_coor)); 903 } 904 PetscCall(DMViewFromOptions(sdm, NULL, "-dm_view")); 905 PetscFunctionReturn(0); 906 } 907 908 int main(int argc, char **argv) 909 { 910 DM dm, sdm; 911 TS ts, sts; 912 Vec u, xtmp; 913 AppCtx user; 914 AdvCtx adv; 915 PetscReal t; 916 PetscInt dim; 917 918 PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 919 PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 920 PetscCall(PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag)); 921 PetscCall(SetupParameters(&user)); 922 PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 923 PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 924 PetscCall(TSSetDM(ts, dm)); 925 PetscCall(DMSetApplicationContext(dm, &user)); 926 /* Discretize chemical species */ 927 PetscCall(DMCreate(PETSC_COMM_WORLD, &sdm)); 928 PetscCall(PetscObjectSetOptionsPrefix((PetscObject) sdm, "part_")); 929 PetscCall(DMSetType(sdm, DMSWARM)); 930 PetscCall(DMGetDimension(dm, &dim)); 931 PetscCall(DMSetDimension(sdm, dim)); 932 PetscCall(DMSwarmSetCellDM(sdm, dm)); 933 /* Setup problem */ 934 PetscCall(SetupDiscretization(dm, sdm, &user)); 935 PetscCall(DMPlexCreateClosureIndex(dm, NULL)); 936 937 PetscCall(DMCreateGlobalVector(dm, &u)); 938 PetscCall(DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace)); 939 940 PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user)); 941 PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user)); 942 PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user)); 943 PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 944 PetscCall(TSSetPreStep(ts, RemoveDiscretePressureNullspace)); 945 PetscCall(TSMonitorSet(ts, MonitorError, &user, NULL)); 946 PetscCall(TSSetFromOptions(ts)); 947 948 PetscCall(TSSetComputeInitialCondition(ts, SetInitialConditions)); /* Must come after SetFromOptions() */ 949 PetscCall(SetInitialConditions(ts, u)); 950 PetscCall(TSGetTime(ts, &t)); 951 PetscCall(DMSetOutputSequenceNumber(dm, 0, t)); 952 PetscCall(DMTSCheckFromOptions(ts, u)); 953 954 /* Setup particle position integrator */ 955 PetscCall(TSCreate(PETSC_COMM_WORLD, &sts)); 956 PetscCall(PetscObjectSetOptionsPrefix((PetscObject) sts, "part_")); 957 PetscCall(PetscObjectIncrementTabLevel((PetscObject) sts, (PetscObject) ts, 1)); 958 PetscCall(TSSetDM(sts, sdm)); 959 PetscCall(TSSetProblemType(sts, TS_NONLINEAR)); 960 PetscCall(TSSetExactFinalTime(sts, TS_EXACTFINALTIME_MATCHSTEP)); 961 PetscCall(TSMonitorSet(sts, MonitorParticleError, &adv, NULL)); 962 PetscCall(TSSetFromOptions(sts)); 963 PetscCall(TSSetApplicationContext(sts, &adv)); 964 PetscCall(TSSetComputeExactError(sts, ComputeParticleError)); 965 PetscCall(TSSetComputeInitialCondition(sts, SetInitialParticleConditions)); 966 adv.ti = t; 967 adv.uf = u; 968 PetscCall(VecDuplicate(adv.uf, &adv.ui)); 969 PetscCall(VecCopy(u, adv.ui)); 970 PetscCall(TSSetRHSFunction(sts, NULL, FreeStreaming, &adv)); 971 PetscCall(TSSetPostStep(ts, AdvectParticles)); 972 PetscCall(PetscObjectCompose((PetscObject) ts, "_SwarmTS", (PetscObject) sts)); 973 PetscCall(DMSwarmVectorDefineField(sdm, DMSwarmPICField_coor)); 974 PetscCall(DMCreateGlobalVector(sdm, &adv.x0)); 975 PetscCall(DMSwarmCreateGlobalVectorFromField(sdm, DMSwarmPICField_coor, &xtmp)); 976 PetscCall(VecCopy(xtmp, adv.x0)); 977 PetscCall(DMSwarmDestroyGlobalVectorFromField(sdm, DMSwarmPICField_coor, &xtmp)); 978 switch(user.solType) { 979 case SOL_TRIG_TRIG: adv.exact = trig_trig_x;break; 980 default: SETERRQ(PetscObjectComm((PetscObject) sdm), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user.solType, NUM_SOL_TYPES)], user.solType); 981 } 982 adv.ctx = &user; 983 984 PetscCall(TSSolve(ts, u)); 985 PetscCall(DMTSCheckFromOptions(ts, u)); 986 PetscCall(PetscObjectSetName((PetscObject) u, "Numerical Solution")); 987 988 PetscCall(VecDestroy(&u)); 989 PetscCall(VecDestroy(&adv.x0)); 990 PetscCall(VecDestroy(&adv.ui)); 991 PetscCall(DMDestroy(&dm)); 992 PetscCall(DMDestroy(&sdm)); 993 PetscCall(TSDestroy(&ts)); 994 PetscCall(TSDestroy(&sts)); 995 PetscCall(PetscBagDestroy(&user.bag)); 996 PetscCall(PetscFinalize()); 997 return 0; 998 } 999 1000 /*TEST 1001 1002 # Swarm does not work with complex 1003 test: 1004 suffix: 2d_tri_p2_p1_p1_tconvp 1005 requires: triangle !single !complex 1006 args: -dm_plex_separate_marker -sol_type trig_trig -dm_refine 2 \ 1007 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 1008 -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1 -ts_monitor_cancel \ 1009 -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 1010 -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 1011 -fieldsplit_0_pc_type lu \ 1012 -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi \ 1013 -omega 0.5 -part_layout_type box -part_lower 0.25,0.25 -part_upper 0.75,0.75 -Npb 5 \ 1014 -part_ts_max_steps 2 -part_ts_dt 0.05 -part_ts_convergence_estimate -convest_num_refine 1 -part_ts_monitor_cancel 1015 test: 1016 suffix: 2d_tri_p2_p1_p1_exit 1017 requires: triangle !single !complex 1018 args: -dm_plex_separate_marker -sol_type trig_trig -dm_refine 1 \ 1019 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \ 1020 -dmts_check .001 -ts_max_steps 10 -ts_dt 0.1 \ 1021 -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 1022 -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 1023 -fieldsplit_0_pc_type lu \ 1024 -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi \ 1025 -omega 0.5 -part_layout_type box -part_lower 0.25,0.25 -part_upper 0.75,0.75 -Npb 5 \ 1026 -part_ts_max_steps 20 -part_ts_dt 0.05 1027 1028 TEST*/ 1029