Lines Matching +full:- +full:j

6 The 2D driven cavity problem is solved in a velocity-vorticity formulation.\n\
8 -lidvelocity &ltlid&gt, where &ltlid&gt = dimensionless velocity of lid\n\
9 -grashof &ltgr&gt, where &ltgr&gt = dimensionless temperature gradent\n\
10 -prandtl &ltpr&gt, where &ltpr&gt = dimensionless thermal/momentum diffusity ratio\n\
11 -contours : draw contour plots of solution\n\n";
18 /*F-----------------------------------------------------------------------
25 - \triangle U - \nabla_y \Omega & = & 0 \\
26 - \triangle V + \nabla_x\Omega & = & 0 \\
27 - \triangle \Omega + \nabla \cdot ([U*\Omega,V*\Omega]) - GR* \nabla_x T & = & 0 \\
28 - \triangle T + PR* \nabla \cdot ([U*T,V*T]) & = & 0
33 No-slip, rigid-wall Dirichlet conditions are used for $ [U,V]$.
35 vorticity: $ \Omega = - \nabla_y U + \nabla_x V$, where along each
41 A finite difference approximation with the usual 5-point stencil
48 * applied matrix-free via the option -snes_mf
50 without a preconditioner due to ill-conditioning).
52 ------------------------------------------------------------------------F*/
58 petscsys.h - base PETSc routines petscvec.h - vectors
59 petscmat.h - matrices
60 petscis.h - index sets petscksp.h - Krylov subspace methods
61 petscviewer.h - viewers petscpc.h - preconditioners
62 petscksp.h - linear solvers
74 User-defined routines and data structures
84 PetscBool draw_contours; /* flag - 1 indicates drawing contours */
92 AppCtx user; /* user-defined work context */ in main()
122 PetscCall(PetscOptionsGetReal(NULL, NULL, "-lidvelocity", &user.lidvelocity, NULL)); in main()
123 PetscCall(PetscOptionsGetReal(NULL, NULL, "-prandtl", &user.prandtl, NULL)); in main()
124 PetscCall(PetscOptionsGetReal(NULL, NULL, "-grashof", &user.grashof, NULL)); in main()
125 PetscCall(PetscOptionsHasName(NULL, NULL, "-contours", &user.draw_contours)); in main()
132 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - in main()
135 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ in main()
137 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - in main()
140 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ in main()
146 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - in main()
148 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ in main()
162 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - in main()
165 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ in main()
173 /* ------------------------------------------------------------------- */
176 FormInitialGuess - Forms initial approximation.
179 user - user-defined application context
180 X - vector
183 X - vector
187 PetscInt i, j, mx, xs, ys, xm, ym; in FormInitialGuess() local
192 grashof = user->grashof; in FormInitialGuess()
195 dx = 1.0 / (mx - 1); in FormInitialGuess()
198 Get local grid boundaries (for 2-dimensional DMDA): in FormInitialGuess()
199 xs, ys - starting grid indices (no ghost points) in FormInitialGuess()
200 xm, ym - widths of local grid (no ghost points) in FormInitialGuess()
206 - For default PETSc vectors, VecGetArray() returns a pointer to in FormInitialGuess()
208 - You MUST call VecRestoreArray() when you no longer need access to in FormInitialGuess()
217 for (j = ys; j < ys + ym; j++) { in FormInitialGuess()
219 x[j][i].u = 0.0; in FormInitialGuess()
220 x[j][i].v = 0.0; in FormInitialGuess()
221 x[j][i].omega = 0.0; in FormInitialGuess()
222 x[j][i].temp = (grashof > 0) * i * dx; in FormInitialGuess()
236 PetscInt xints, xinte, yints, yinte, i, j; in FormFunctionLocal() local
242 grashof = user->grashof; in FormFunctionLocal()
243 prandtl = user->prandtl; in FormFunctionLocal()
244 lid = user->lidvelocity; in FormFunctionLocal()
253 dhx = (PetscReal)(info->mx - 1); in FormFunctionLocal()
254 dhy = (PetscReal)(info->my - 1); in FormFunctionLocal()
260 xints = info->xs; in FormFunctionLocal()
261 xinte = info->xs + info->xm; in FormFunctionLocal()
262 yints = info->ys; in FormFunctionLocal()
263 yinte = info->ys + info->ym; in FormFunctionLocal()
267 j = 0; in FormFunctionLocal()
270 for (i = info->xs; i < info->xs + info->xm; i++) { in FormFunctionLocal()
271 f[j][i].u = x[j][i].u; in FormFunctionLocal()
272 f[j][i].v = x[j][i].v; in FormFunctionLocal()
273 f[j][i].omega = x[j][i].omega + (x[j + 1][i].u - x[j][i].u) * dhy; in FormFunctionLocal()
274 f[j][i].temp = x[j][i].temp - x[j + 1][i].temp; in FormFunctionLocal()
279 if (yinte == info->my) { in FormFunctionLocal()
280 j = info->my - 1; in FormFunctionLocal()
281 yinte = yinte - 1; in FormFunctionLocal()
283 for (i = info->xs; i < info->xs + info->xm; i++) { in FormFunctionLocal()
284 f[j][i].u = x[j][i].u - lid; in FormFunctionLocal()
285 f[j][i].v = x[j][i].v; in FormFunctionLocal()
286 f[j][i].omega = x[j][i].omega + (x[j][i].u - x[j - 1][i].u) * dhy; in FormFunctionLocal()
287 f[j][i].temp = x[j][i].temp - x[j - 1][i].temp; in FormFunctionLocal()
296 for (j = info->ys; j < info->ys + info->ym; j++) { in FormFunctionLocal()
297 f[j][i].u = x[j][i].u; in FormFunctionLocal()
298 f[j][i].v = x[j][i].v; in FormFunctionLocal()
299 f[j][i].omega = x[j][i].omega - (x[j][i + 1].v - x[j][i].v) * dhx; in FormFunctionLocal()
300 f[j][i].temp = x[j][i].temp; in FormFunctionLocal()
305 if (xinte == info->mx) { in FormFunctionLocal()
306 i = info->mx - 1; in FormFunctionLocal()
307 xinte = xinte - 1; in FormFunctionLocal()
309 for (j = info->ys; j < info->ys + info->ym; j++) { in FormFunctionLocal()
310 f[j][i].u = x[j][i].u; in FormFunctionLocal()
311 f[j][i].v = x[j][i].v; in FormFunctionLocal()
312 f[j][i].omega = x[j][i].omega - (x[j][i].v - x[j][i - 1].v) * dhx; in FormFunctionLocal()
313 f[j][i].temp = x[j][i].temp - (PetscReal)(grashof > 0); in FormFunctionLocal()
318 for (j = yints; j < yinte; j++) { in FormFunctionLocal()
323 vx = x[j][i].u; in FormFunctionLocal()
326 vxm = .5 * (vx - avx); in FormFunctionLocal()
327 vy = x[j][i].v; in FormFunctionLocal()
330 vym = .5 * (vy - avy); in FormFunctionLocal()
333 u = x[j][i].u; in FormFunctionLocal()
334 uxx = (2.0 * u - x[j][i - 1].u - x[j][i + 1].u) * hydhx; in FormFunctionLocal()
335 uyy = (2.0 * u - x[j - 1][i].u - x[j + 1][i].u) * hxdhy; in FormFunctionLocal()
336 f[j][i].u = uxx + uyy - .5 * (x[j + 1][i].omega - x[j - 1][i].omega) * hx; in FormFunctionLocal()
339 u = x[j][i].v; in FormFunctionLocal()
340 uxx = (2.0 * u - x[j][i - 1].v - x[j][i + 1].v) * hydhx; in FormFunctionLocal()
341 uyy = (2.0 * u - x[j - 1][i].v - x[j + 1][i].v) * hxdhy; in FormFunctionLocal()
342 f[j][i].v = uxx + uyy + .5 * (x[j][i + 1].omega - x[j][i - 1].omega) * hy; in FormFunctionLocal()
345 u = x[j][i].omega; in FormFunctionLocal()
346 uxx = (2.0 * u - x[j][i - 1].omega - x[j][i + 1].omega) * hydhx; in FormFunctionLocal()
347 uyy = (2.0 * u - x[j - 1][i].omega - x[j + 1][i].omega) * hxdhy; in FormFunctionLocal()
348j][i].omega = uxx + uyy + (vxp * (u - x[j][i - 1].omega) + vxm * (x[j][i + 1].omega - u)) * hy + (… in FormFunctionLocal()
351 u = x[j][i].temp; in FormFunctionLocal()
352 uxx = (2.0 * u - x[j][i - 1].temp - x[j][i + 1].temp) * hydhx; in FormFunctionLocal()
353 uyy = (2.0 * u - x[j - 1][i].temp - x[j + 1][i].temp) * hxdhy; in FormFunctionLocal()
354j][i].temp = uxx + uyy + prandtl * ((vxp * (u - x[j][i - 1].temp) + vxm * (x[j][i + 1].temp - u)) … in FormFunctionLocal()
359 Flop count (multiply-adds are counted as 2 operations) in FormFunctionLocal()
361 PetscCall(PetscLogFlops(84.0 * info->ym * info->xm)); in FormFunctionLocal()
366 Performs sweeps of point block nonlinear Gauss-Seidel on all the local grid points
374 PetscInt xints, xinte, yints, yinte, i, j, k, l; in NonlinearGS() local
393 grashof = user->grashof; in NonlinearGS()
394 prandtl = user->prandtl; in NonlinearGS()
395 lid = user->lidvelocity; in NonlinearGS()
403 Scatter ghost points to local vector, using the 2-step process in NonlinearGS()
416 dhx = (PetscReal)(info.mx - 1); in NonlinearGS()
417 dhy = (PetscReal)(info.my - 1); in NonlinearGS()
431 j = 0; in NonlinearGS()
435 bjiu = b[j][i].u; in NonlinearGS()
436 bjiv = b[j][i].v; in NonlinearGS()
441 x[j][i].u = 0.0 + bjiu; in NonlinearGS()
442 x[j][i].v = 0.0 + bjiv; in NonlinearGS()
448 j = info.my - 1; in NonlinearGS()
452 bjiu = b[j][i].u; in NonlinearGS()
453 bjiv = b[j][i].v; in NonlinearGS()
458 x[j][i].u = lid + bjiu; in NonlinearGS()
459 x[j][i].v = bjiv; in NonlinearGS()
467 for (j = info.ys; j < info.ys + info.ym; j++) { in NonlinearGS()
469 bjiu = b[j][i].u; in NonlinearGS()
470 bjiv = b[j][i].v; in NonlinearGS()
475 x[j][i].u = 0.0 + bjiu; in NonlinearGS()
476 x[j][i].v = 0.0 + bjiv; in NonlinearGS()
482 i = info.mx - 1; in NonlinearGS()
484 for (j = info.ys; j < info.ys + info.ym; j++) { in NonlinearGS()
486 bjiu = b[j][i].u; in NonlinearGS()
487 bjiv = b[j][i].v; in NonlinearGS()
492 x[j][i].u = 0.0 + bjiu; in NonlinearGS()
493 x[j][i].v = 0.0 + bjiv; in NonlinearGS()
498 for (j = info.ys; j < info.ys + info.ym; j++) { in NonlinearGS()
509 bjiu = b[j][i].u; in NonlinearGS()
510 bjiv = b[j][i].v; in NonlinearGS()
511 bjiomega = b[j][i].omega; in NonlinearGS()
512 bjitemp = b[j][i].temp; in NonlinearGS()
520 if (i != 0 && i != info.mx - 1 && j != 0 && j != info.my - 1) { in NonlinearGS()
522 u = x[j][i].u; in NonlinearGS()
523 uxx = (2.0 * u - x[j][i - 1].u - x[j][i + 1].u) * hydhx; in NonlinearGS()
524 uyy = (2.0 * u - x[j - 1][i].u - x[j + 1][i].u) * hxdhy; in NonlinearGS()
525 fu = uxx + uyy - .5 * (x[j + 1][i].omega - x[j - 1][i].omega) * hx - bjiu; in NonlinearGS()
528 u = x[j][i].v; in NonlinearGS()
529 uxx = (2.0 * u - x[j][i - 1].v - x[j][i + 1].v) * hydhx; in NonlinearGS()
530 uyy = (2.0 * u - x[j - 1][i].v - x[j + 1][i].v) * hxdhy; in NonlinearGS()
531 fv = uxx + uyy + .5 * (x[j][i + 1].omega - x[j][i - 1].omega) * hy - bjiv; in NonlinearGS()
536 vx = x[j][i].u; in NonlinearGS()
539 vxm = .5 * (vx - avx); in NonlinearGS()
540 vy = x[j][i].v; in NonlinearGS()
543 vym = .5 * (vy - avy); in NonlinearGS()
545 u = x[j][i].omega; in NonlinearGS()
546 uxx = (2.0 * u - x[j][i - 1].omega - x[j][i + 1].omega) * hydhx; in NonlinearGS()
547 uyy = (2.0 * u - x[j - 1][i].omega - x[j + 1][i].omega) * hxdhy; in NonlinearGS()
548- x[j][i - 1].omega) + vxm * (x[j][i + 1].omega - u)) * hy + (vyp * (u - x[j - 1][i].omega) + vym … in NonlinearGS()
550 dfodo = 2.0 * (hydhx + hxdhy) + ((vxp - vxm) * hy + (vyp - vym) * hx); in NonlinearGS()
551 if (PetscRealPart(vx) > 0.0) dfodu = (u - x[j][i - 1].omega) * hy; in NonlinearGS()
552 else dfodu = (x[j][i + 1].omega - u) * hy; in NonlinearGS()
554 if (PetscRealPart(vy) > 0.0) dfodv = (u - x[j - 1][i].omega) * hx; in NonlinearGS()
555 else dfodv = (x[j + 1][i].omega - u) * hx; in NonlinearGS()
558 u = x[j][i].temp; in NonlinearGS()
559 uxx = (2.0 * u - x[j][i - 1].temp - x[j][i + 1].temp) * hydhx; in NonlinearGS()
560 uyy = (2.0 * u - x[j - 1][i].temp - x[j + 1][i].temp) * hxdhy; in NonlinearGS()
561 … ((vxp * (u - x[j][i - 1].temp) + vxm * (x[j][i + 1].temp - u)) * hy + (vyp * (u - x[j - 1][i].tem… in NonlinearGS()
562 dftdt = 2.0 * (hydhx + hxdhy) + prandtl * ((vxp - vxm) * hy + (vyp - vym) * hx); in NonlinearGS()
563 if (PetscRealPart(vx) > 0.0) dftdu = prandtl * (u - x[j][i - 1].temp) * hy; in NonlinearGS()
564 else dftdu = prandtl * (x[j][i + 1].temp - u) * hy; in NonlinearGS()
566 if (PetscRealPart(vy) > 0.0) dftdv = prandtl * (u - x[j - 1][i].temp) * hx; in NonlinearGS()
567 else dftdv = prandtl * (x[j + 1][i].temp - u) * hx; in NonlinearGS()
574 by simple back-substitution in NonlinearGS()
578 yo = (fomega - (dfodu * yu + dfodv * yv)) / dfodo; in NonlinearGS()
579 yt = (ftemp - (dftdu * yu + dftdv * yv)) / dftdt; in NonlinearGS()
581 x[j][i].u = x[j][i].u - yu; in NonlinearGS()
582 x[j][i].v = x[j][i].v - yv; in NonlinearGS()
583 x[j][i].temp = x[j][i].temp - yt; in NonlinearGS()
584 x[j][i].omega = x[j][i].omega - yo; in NonlinearGS()
587 fomega = x[j][i].omega - (x[j][i + 1].v - x[j][i].v) * dhx - bjiomega; in NonlinearGS()
588 ftemp = x[j][i].temp - bjitemp; in NonlinearGS()
591 x[j][i].omega = x[j][i].omega - fomega; in NonlinearGS()
592 x[j][i].temp = x[j][i].temp - ftemp; in NonlinearGS()
594 if (i == info.mx - 1) { in NonlinearGS()
595 fomega = x[j][i].omega - (x[j][i].v - x[j][i - 1].v) * dhx - bjiomega; in NonlinearGS()
596 ftemp = x[j][i].temp - (PetscReal)(grashof > 0) - bjitemp; in NonlinearGS()
599 x[j][i].omega = x[j][i].omega - fomega; in NonlinearGS()
600 x[j][i].temp = x[j][i].temp - ftemp; in NonlinearGS()
602 if (j == 0) { in NonlinearGS()
603 fomega = x[j][i].omega + (x[j + 1][i].u - x[j][i].u) * dhy - bjiomega; in NonlinearGS()
604 ftemp = x[j][i].temp - x[j + 1][i].temp - bjitemp; in NonlinearGS()
607 x[j][i].omega = x[j][i].omega - fomega; in NonlinearGS()
608 x[j][i].temp = x[j][i].temp - ftemp; in NonlinearGS()
610 if (j == info.my - 1) { in NonlinearGS()
611 fomega = x[j][i].omega + (x[j][i].u - x[j - 1][i].u) * dhy - bjiomega; in NonlinearGS()
612 ftemp = x[j][i].temp - x[j - 1][i].temp - bjitemp; in NonlinearGS()
615 x[j][i].omega = x[j][i].omega - fomega; in NonlinearGS()
616 x[j][i].temp = x[j][i].temp - ftemp; in NonlinearGS()
623 …pxnorm = PetscRealPart(x[j][i].u * x[j][i].u + x[j][i].v * x[j][i].v + x[j][i].omega * x[j][i].ome… in NonlinearGS()
645 args: -da_refine 3 -snes_monitor_short -pc_type mg -ksp_type fgmres -pc_mg_type full
651 …args: -snes_monitor_short -ksp_monitor_short -pc_type fieldsplit -pc_fieldsplit_type symmetric_mul…
658-snes_monitor_short -pc_type redundant -dm_mat_type mpiaij -redundant_pc_factor_mat_solver_type pa…
664-snes_monitor_short -pc_type redundant -dm_mat_type mpiaij -redundant_pc_factor_mat_solver_type pa…
669 …args: -snes_monitor_short -ksp_monitor_short -pc_type fieldsplit -pc_fieldsplit_type multiplicativ…
675 …args: -snes_monitor_short -pc_type mg -dm_mat_type baij -mg_coarse_pc_type bjacobi -da_refine 3 -k…
681 …args: -snes_converged_reason -pc_type mg -dm_mat_type baij -mg_coarse_pc_type bjacobi -da_refine 3…
687 args: -snes_monitor_short -ksp_pc_side right
692 args: -snes_monitor_ksp draw::draw_lg -ksp_pc_side right
698 …args: -da_refine 3 -snes_monitor_short -pc_type mg -ksp_type fgmres -pc_mg_type full -snes_type ne…
704 …args: -da_refine 3 -snes_monitor_short -pc_type mg -ksp_type fgmres -pc_mg_type full -snes_type ne…
710 args: -da_refine 3 -snes_converged_reason -pc_type mg -mat_fd_type ds
716 args: -da_refine 3 -snes_converged_reason -pc_type mg -mat_fd_type ds -mat_fd_coloring_bcols
724 …gs: -da_refine 3 -snes_monitor_short -pc_type redundant -dm_mat_type mpiaij -redundant_ksp_type pr…
730 …gs: -da_refine 3 -snes_monitor_short -pc_type redundant -dm_mat_type mpiaij -redundant_ksp_type pr…
735 …args: -snes_monitor_short -ksp_monitor_short -pc_type fieldsplit -snes_view -ksp_type fgmres -da_r…
741 …args: -snes_monitor_short -ksp_monitor_short -pc_type fieldsplit -snes_view -da_refine 1 -ksp_type…
746-snes_monitor_short -ksp_monitor_short -pc_type fieldsplit -pc_fieldsplit_block_size 2 -pc_fieldsp…
752 …args: -snes_monitor_short -ksp_monitor_short -pc_type fieldsplit -pc_fieldsplit_type multiplicativ…
758 …args: -da_refine 3 -da_overlap 2 -snes_monitor_short -snes_type aspin -grashof 4e4 -lidvelocity 10…
764 args: -ksp_type bcgsl -ksp_monitor_short -da_refine 2 -ksp_bcgsl_ell 3 -snes_view
770 …args: -da_refine 3 -snes_monitor_short -pc_type mg -ksp_type fgmres -pc_mg_type full -mat_fd_color…
777-da_refine 4 -ksp_type fgmres -pc_type bjacobi -pc_bjacobi_blocks 2 -sub_ksp_type gmres -sub_ksp_m…
782 …args: -da_refine 2 -pc_type lu -ksp_type cgne -ksp_monitor_short -ksp_converged_reason -ksp_view -
783 filter: grep -v HERMITIAN
788 args: -da_refine 1 -ksp_monitor_short -ksp_type cgs
793-ksp_type fgmres -pc_type composite -pc_composite_type MULTIPLICATIVE -pc_composite_pcs fieldsplit…
798-ksp_type fgmres -pc_type composite -pc_composite_type MULTIPLICATIVE -pc_composite_pcs fieldsplit…
804-ksp_type fgmres -pc_type composite -pc_composite_type MULTIPLICATIVE -pc_composite_pcs fieldsplit…
810-da_refine 3 -grashof 4e4 -lidvelocity 100 -snes_monitor_short -snes_type composite -snes_composit…
816 …args: -dm_vec_type cuda -dm_mat_type aijcusparse -pc_type none -ksp_type fgmres -snes_monitor_shor…
821 …args: -dm_vec_type hip -dm_mat_type aijhipsparse -pc_type none -ksp_type fgmres -snes_monitor_shor…
825-pc_type fieldsplit -snes_view draw -fieldsplit_x_velocity_pc_type mg -fieldsplit_x_velocity_pc_mg…
830 args: -snes_monitor_solution draw::draw_ports -da_refine 1
836-da_refine 4 -snes_monitor_short -snes_type fas -fas_levels_snes_type ngs -fas_levels_snes_ngs_swe…
841-da_refine 4 -snes_monitor_short -snes_type fas -snes_fas_type full -snes_fas_full_downsweep -fas_…
846 args: -da_refine 3 -snes_converged_reason -pc_type mg -mat_fd_type ds
852 args: -da_refine 3 -snes_converged_reason -pc_type mg -mat_fd_type ds -dm_mat_type baij
858 args: -da_refine 3 -snes_converged_reason -pc_type mg -mat_fd_type ds -mat_fd_coloring_bcols 1
864 args: -da_refine 3 -snes_monitor_short -pc_type mg
869 args: -da_refine 3 -snes_monitor_short -pc_type mg -dm_mat_type baij
875 args: -da_refine 3 -snes_monitor_short -pc_type mg -mat_fd_coloring_bcols 1
881-ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_block_size 4 -pc_fieldsplit_type additive -pc_…
886-ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_block_size 4 -pc_fieldsplit_type additive -pc_…
891-ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_block_size 4 -pc_fieldsplit_type SCHUR -pc_fie…
899-pc_type fieldsplit -pc_fieldsplit_block_size 4 -pc_fieldsplit_type SCHUR -pc_fieldsplit_0_fields …
905-pc_type fieldsplit -pc_fieldsplit_block_size 4 -pc_fieldsplit_type SCHUR -pc_fieldsplit_0_fields …
911 …args: -da_refine 3 -snes_monitor_short -snes_fd_color -snes_fd_color_use_mat -mat_coloring_type gr…
919 args: -da_refine 3 -snes_monitor_short -pc_type hypre -ksp_norm_type unpreconditioned
925 args: -ksp_type ibcgs -ksp_monitor_short -da_refine 2 -snes_view
931 args: -pc_type kaczmarz -ksp_monitor_short -snes_monitor_short -snes_view
937 args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type klu
943 …args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type klu -pc_factor_mat_orderi…
949 … args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type klu -mat_klu_use_btf 0
956 args: -da_refine 3 -snes_monitor_short -pc_type ml
960-da_refine 4 -snes_monitor_short -snes_type ngmres -npc_fas_levels_snes_type ngs -npc_fas_levels_s…
965-da_refine 3 -snes_monitor_short -snes_type ngmres -npc_snes_type fas -npc_fas_levels_snes_type ng…
971-snes_grid_sequence 2 -lidvelocity 200 -grashof 1e4 -snes_monitor_short -snes_view -snes_converged…
977 …args: -da_refine 4 -da_overlap 2 -snes_monitor_short -snes_type ngmres -snes_max_it 10 -npc_snes_t…
982 args: -snes_type ngs -snes_view -snes_monitor -snes_rtol 1e-4
987 args: -snes_type ngs -snes_ngs_secant -snes_view -snes_monitor -snes_rtol 1e-4
994 args: -pc_type parms -ksp_monitor_short -snes_view
999 args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type superlu
1004 …args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type superlu -dm_mat_type sell…
1010 args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type superlu_dist
1017 args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type superlu_dist
1024 filter: grep -v iam | grep -v openMP
1025 …: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type superlu_dist -mat_superlu_dis…
1031 …args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type superlu_dist -pc_precisio…
1038 …args: -da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type mumps -pc_precision {{sin…
1045 filter: grep -v iam | grep -v openMP
1046-da_grid_x 20 -da_grid_y 20 -pc_type lu -pc_factor_mat_solver_type superlu_dist -mat_superlu_dist_…
1051 …args: -da_grid_x 20 -da_grid_y 20 -{snes,ksp}_monitor_short -pc_type lu -pc_factor_mat_solver_type…
1056 …args: -da_grid_x 20 -da_grid_y 20 -{snes,ksp}_monitor_short -pc_type lu -pc_factor_mat_solver_type…
1061 args: -da_refine 1 -ksp_monitor_short -ksp_type tcqmr
1066 args: -da_refine 1 -ksp_monitor_short -ksp_type tfqmr
1072 …args: -da_refine 2 -pc_type lu -pc_factor_mat_solver_type umfpack -snes_view -snes_monitor_short -
1078 args: -da_refine 5 -snes_monitor -ksp_monitor -snes_view
1084 args: -da_refine 5 -snes_monitor -ksp_monitor -snes_view -pc_type mg
1091 … args: -da_refine 5 -snes_monitor -snes_converged_reason -pc_type hypre -dm_mat_type {{aij baij}}
1097 …args: -da_refine 1 -snes_monitor -snes_converged_reason -pc_type hypre -dm_mat_type {{seqaij mpiai…
1103 args: -da_refine 5 -snes_monitor -ksp_monitor -snes_view -pc_type ml
1109 args: -da_refine 5 -log_view
1110 filter: head -n 2
1111 filter_output: head -n 2
1117 args: -da_refine 5 -log_view -pc_type mg
1118 filter: head -n 2
1119 filter_output: head -n 2
1125 args: -da_refine 5 -log_view
1126 filter: head -n 2
1127 filter_output: head -n 2
1133 args: -da_refine 5 -log_view -pc_type mg
1134 filter: head -n 2
1135 filter_output: head -n 2
1141 …args: -snes_monitor -dm_mat_type seqaijcusparse -dm_vec_type seqcuda -pc_type gamg -ksp_monitor -m…
1147 …args: -snes_monitor -dm_mat_type mpiaijcusparse -dm_vec_type mpicuda -pc_type gamg -ksp_monitor -m…
1153 …s: -dm_mat_type aijcusparse -dm_vec_type cuda -da_refine 3 -pc_type mg -mg_levels_ksp_type chebysh…
1160 …args: -snes_monitor -dm_mat_type mpiaijhipsparse -dm_vec_type hip -pc_type gamg -ksp_monitor -mg_l…
1166 …args: -snes_monitor -dm_mat_type mpiaijhipsparse -dm_vec_type mpihip -pc_type gamg -ksp_monitor -m…
1172 …s: -dm_mat_type aijhipsparse -dm_vec_type hip -da_refine 3 -pc_type mg -mg_levels_ksp_type chebysh…
1179-dm_mat_type aijviennacl -dm_vec_type viennacl -da_refine 3 -pc_type mg -mg_levels_ksp_type cheby…
1186 args: -dm_mat_type baij -snes_monitor -ksp_monitor -snes_view
1192 args: -dm_mat_type baij -snes_monitor -ksp_monitor -snes_view
1198 args: -pc_type lu -pc_factor_mat_solver_type mkl_cpardiso -ksp_monitor
1203 args: -log_view -log_view_memory -da_refine 4
1204 filter: grep MatFDColorSetUp | wc -w | xargs -I % sh -c "expr % \> 21"
1209-pc_type fieldsplit -da_refine 3 -all_ksp_monitor -fieldsplit_y_velocity_pc_type lu -fieldsplit_te…
1213 args: -mat_type aij -pc_type asm -pc_asm_sub_mat_type dense -snes_view
1219 args: -da_refine 2 -ksp_monitor -snes_monitor -snes_view -pc_type hypre -pc_hypre_type euclid
1225 …args: -da_refine 2 -ksp_monitor -snes_monitor -snes_view -pc_type hypre -pc_hypre_type euclid -pc_…
1231 …args: -da_refine 2 -ksp_monitor -snes_monitor -snes_view -pc_type hypre -pc_hypre_type euclid -pc_…
1237 args: -da_refine 100 -petsc_ci_portable_error_output -error_output_stdout
1238 …filter: grep -E -v "(memory block|leaked context|not freed before MPI_Finalize|Could be the progra…
1242 …args: -snes_monitor -ksp_converged_reason -ksp_type hpddm -pc_type jacobi -dm_mat_type aijcusparse…
1245 …filter: sed -e "s/Linear solve converged due to CONVERGED_RTOL iterations 15/Linear solve converge…
1246 …args: -ksp_hpddm_type {{gmres gcrodr}separate output} -ksp_hpddm_precision {{single double}shared …
1249 …filter: sed -e "s/Linear solve converged due to CONVERGED_RTOL iterations 15/Linear solve converge…
1250 args: -ksp_hpddm_type gcrodr -ksp_pc_side right
1251 output_file: output/ex19_hpddm_cuda_ksp_hpddm_type-gcrodr.out