<html><head><meta http-equiv="Content-Type" content="text/html; charset=us-ascii"></head><body style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;" class=""><div class=""><br class=""></div> You can run the hyper case with -fp_trap in a debugger and it may stop exactly in hypre where the inf/nan appears. Or if you know how to catch floating point exceptions in your debugger run in the debugger with that option to find where in hypre the inf/nan appears.<div class=""><br class=""></div><div class=""> Barry</div><div class=""><br class=""><div><br class=""><blockquote type="cite" class=""><div class="">On Jul 1, 2020, at 8:22 AM, Andrea Iob <<a href="mailto:andrea_iob@hotmail.com" class="">andrea_iob@hotmail.com</a>> wrote:</div><br class="Apple-interchange-newline"><div class="">
<meta http-equiv="Content-Type" content="text/html; charset=iso-8859-1" class="">
<style type="text/css" style="display:none;" class=""> P {margin-top:0;margin-bottom:0;} </style>
<div dir="ltr" class="">
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
Thanks for the answers.</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
<br class="">
</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
The matrix contains both the advection and the diffusion terms of the Navier-Stokes equations. The case is subsonic (M=0.5), Reynolds number is in the order of 5000 so I expect the hyperbolic features to be the dominant ones. I understand this is not the type
of problems AMG is designed for, but I was hoping to still be able to solve the system. This is the first time I use BoomerAMG and I'm not familar with all the options, however all the tests I've done gave me INFs/NANs.</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
<br class="">
</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
You suggestion to use GAMG with custom smoothers seems to do the trick: after setting the GAMG smoother of each level to FGMRES/ILU, the solver converges in 93 iterations (overall it seems faster than ILU). Great! I've set the smoothers like this:<br class="">
</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
<br class="">
</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
<span style="font-family: Consolas, Courier, monospace;" class=""> for (int i = 0; i < nLevels; ++i) {</span><span class=""><br class="">
</span>
<div class=""><span style="font-family: Consolas, Courier, monospace;" class=""> KSP levelSmoother;</span><br class="">
</div>
<div class=""><span style="font-family: Consolas, Courier, monospace;" class=""> PCMGGetSmoother(preconditioner, i, &levelSmoother);</span><br class="">
</div>
<div class=""><br class="">
</div>
<div class=""><span style="font-family: Consolas, Courier, monospace;" class=""> PC levelSmootherPreconditioner;</span><br class="">
</div>
<div class=""><span style="font-family: Consolas, Courier, monospace;" class=""> KSPGetPC(levelSmoother, &levelSmootherPreconditioner);</span><br class="">
</div>
<div class=""><span style="font-family: Consolas, Courier, monospace;" class=""> PCSetType(levelSmootherPreconditioner, PCILU);</span><br class="">
</div>
<div class=""><span style="font-family: Consolas, Courier, monospace;" class=""> </span>
<br class="">
</div>
<div class=""><span style="font-family: Consolas, Courier, monospace;" class=""> KSPSetType(levelSmoother, KSPFGMRES);</span><br class="">
</div>
<span style="font-family: Consolas, Courier, monospace;" class=""> }</span></div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
<span class=""><br class="">
</span></div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
<span class="">Thank you very much!</span></div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
<span class="">Best regards,</span></div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
<span class="">Andrea</span><br class="">
</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
<br class="">
</div>
<div class="">
<hr tabindex="-1" style="display:inline-block; width:98%" class="">
<div id="divRplyFwdMsg" dir="ltr" class=""><font style="font-size:11pt" face="Calibri, sans-serif" class=""><b class="">From:</b> Mark Adams <<a href="mailto:mfadams@lbl.gov" class="">mfadams@lbl.gov</a>><br class="">
<b class="">Sent:</b> Wednesday, July 1, 2020 2:15 PM<br class="">
<b class="">To:</b> Matthew Knepley <<a href="mailto:knepley@gmail.com" class="">knepley@gmail.com</a>><br class="">
<b class="">Cc:</b> Andrea Iob <<a href="mailto:andrea_iob@hotmail.com" class="">andrea_iob@hotmail.com</a>>; <a href="mailto:petsc-users@mcs.anl.gov" class="">petsc-users@mcs.anl.gov</a> <<a href="mailto:petsc-users@mcs.anl.gov" class="">petsc-users@mcs.anl.gov</a>><br class="">
<b class="">Subject:</b> Re: [petsc-users] DIVERGED_NANORINF when using HYPRE/BoomerAMG</font>
<div class=""> </div>
</div>
<div class="">
<div dir="ltr" class="">As Matt said AMG does not work for everything out-of-the-box.<br class="">
<div class=""><br class="">
</div>
<div class="">Hypre is pretty robust and I am surprised that every option gives INFs/NANs. But if it solves your stretched Lapacian then it's probably hooked up correctly.</div>
<div class=""><br class="">
</div>
You might try a non-adjoint Navier-Stokes problem, if you can.
<div class=""><br class="">
</div>
<div class="">You could try 'gamg' with the fgmres/ilu solver that works, as your smoother (use -ksp_view to check that you get everything set correctly). Use fgmres as the outer Krylov method and use -pc_gamg_agg_nsmooths 0.</div>
<div class=""><br class="">
</div>
<div class="">Mark<br class="">
<div class=""><span style="font-family: Calibri, Helvetica, sans-serif; font-size: 16px;" class=""><br class="">
</span></div>
</div>
</div>
<br class="">
<div class="x_gmail_quote">
<div dir="ltr" class="x_gmail_attr">On Tue, Jun 30, 2020 at 9:16 AM Matthew Knepley <<a href="mailto:knepley@gmail.com" class="">knepley@gmail.com</a>> wrote:<br class="">
</div>
<blockquote class="x_gmail_quote" style="margin:0px 0px 0px 0.8ex; border-left:1px solid rgb(204,204,204); padding-left:1ex">
<div dir="ltr" class="">
<div dir="ltr" class="">On Tue, Jun 30, 2020 at 7:59 AM Andrea Iob <<a href="mailto:andrea_iob@hotmail.com" target="_blank" class="">andrea_iob@hotmail.com</a>> wrote:<br class="">
</div>
<div class="x_gmail_quote">
<blockquote class="x_gmail_quote" style="margin:0px 0px 0px 0.8ex; border-left:1px solid rgb(204,204,204); padding-left:1ex">
<div dir="ltr" class="">
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
Hi,</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
<br class="">
</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
I'm trying to solve a linear system using HYPRE/BoomerAMG as preconditioner. The system comes from a two-dimensional adjoint Navier-Stokes problem.</div>
</div>
</blockquote>
<div class=""><br class="">
</div>
<div class="">Do you expect the system you are solving to be elliptic? BoomerAMG is designed for elliptic systems, and can easily fail if applied</div>
<div class="">to a more general system, say one with zeros on the diagonal.</div>
<div class=""><br class="">
</div>
<div class=""> Thanks,</div>
<div class=""><br class="">
</div>
<div class=""> Matt</div>
<div class=""> </div>
<blockquote class="x_gmail_quote" style="margin:0px 0px 0px 0.8ex; border-left:1px solid rgb(204,204,204); padding-left:1ex">
<div dir="ltr" class="">
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
The mesh is structured (6400 cells) and there are 24 DoFs for each cell (the matrix has a total of 153600 rows). The condition number of the matrix should be in the order of 1e9. Using ILU preconditioner and FGMRES, the system is correctly solved (it takes
800 iterations to reach a residual of 1e-12). However, replacing ILU with HYPRE/BoomerAMG, the solver stops right after the first iteration (DIVERGED_NANORINF). I've tried different BoomerAMG options and different solvers (including Richardson to use BoomerAMG
without a Krylov method), but I always get the same DIVERGED_NANORINF error.<br class="">
</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
<br class="">
</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
Using the same HYPRE/BoomerAMG + FGMERS setup on another problem (Laplacian on a highly stretched grid), the solver reaches convergence very quickly without any problems.
<br class="">
</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
<br class="">
</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
At the bottom of this mail, you'll find the log of a run that stops with DIVERGED_NANORINF. Do you spot something in my setup that may explain why the solver is not converging?
<br class="">
</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
<br class="">
</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
Thanks. Best regards,</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
Andrea</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
<br class="">
</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
<font size="2" class=""><span style="font-size:11pt" class="">-----------------------------------------</span></font><br class="">
</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt;" class="">
<font size="2" class=""><span style="font-size:11pt" class="">Assembly started...<br class="">
Reading matrix file...<br class="">
- Number of rows (N) = 153600<br class="">
- Number of non-zero entries (NNZ) = 18247680<br class="">
Assembly completed.<br class="">
Time elapsed 83.3231s<br class="">
Set up started...<br class="">
Set up completed.<br class="">
Time elapsed 6.62441s<br class="">
Solution started...<br class="">
0 KSP unpreconditioned resid norm 7.736650641501e-01 true resid norm 7.736650641501e-01 ||r(i)||/||b|| 1.000000000000e+00<br class="">
0 KSP Residual norm 7.736650641501e-01 % max 1.000000000000e+00 min 1.000000000000e+00 max/min 1.000000000000e+00<br class="">
[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------<br class="">
[0]PETSC ERROR:<br class="">
[0]PETSC ERROR: KSPSolve has not converged due to Nan or Inf inner product<br class="">
[0]PETSC ERROR: See <a href="http://www.mcs.anl.gov/petsc/documentation/faq.html" rel="noopener noreferrer" target="_blank" class="">
http://www.mcs.anl.gov/petsc/documentation/faq.html</a> for trouble shooting.<br class="">
[0]PETSC ERROR: Petsc Release Version 3.10.3, unknown<br class="">
[0]PETSC ERROR: bitpit_system_solver on a named xxx by xxx Tue Jun 30 09:42:09 2020<br class="">
[0]PETSC ERROR: Configure options PETSC_ARCH=arch-linux2-c-opt-gcc7-hypre --with-blaslapack-dir=/opt/lapack/3.8.0-gcc.7.4.0/lib64 --with-debugging=0 COPTFLAGS=-O3 CXXOPTFLAGS=-O3 FOPTFLAGS=-O3 --with-valgrind=1 --with-valgrind-dir=/opt/valgrind/3.14.0/ --prefix=/opt/petsc/3.10.3_gcc-7.4.0
--with-mpi=1 --download-hypre<br class="">
[0]PETSC ERROR: #1 KSPGMRESClassicalGramSchmidtOrthogonalization() line 67 in /root/InstallSources/petsc/src/ksp/ksp/impls/gmres/borthog2.c<br class="">
[0]PETSC ERROR: #2 KSPFGMRESCycle() line 175 in /root/InstallSources/petsc/src/ksp/ksp/impls/gmres/fgmres/fgmres.c<br class="">
[0]PETSC ERROR: #3 KSPSolve_FGMRES() line 291 in /root/InstallSources/petsc/src/ksp/ksp/impls/gmres/fgmres/fgmres.c<br class="">
[0]PETSC ERROR: #4 KSPSolve() line 780 in /root/InstallSources/petsc/src/ksp/ksp/interface/itfunc.c<br class="">
KSP Object: 1 MPI processes<br class="">
type: fgmres<br class="">
restart=45, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement<br class="">
happy breakdown tolerance 1e-30<br class="">
maximum iterations=10000, nonzero initial guess<br class="">
tolerances: relative=1e-12, absolute=1e-50, divergence=10000.<br class="">
right preconditioning<br class="">
using UNPRECONDITIONED norm type for convergence test<br class="">
PC Object: 1 MPI processes<br class="">
type: hypre<br class="">
HYPRE BoomerAMG preconditioning<br class="">
Cycle type V<br class="">
Maximum number of levels 25<br class="">
Maximum number of iterations PER hypre call 1<br class="">
Convergence tolerance PER hypre call 0.<br class="">
Threshold for strong coupling 0.7<br class="">
Interpolation truncation factor 0.3<br class="">
Interpolation: max elements per row 2<br class="">
Number of levels of aggressive coarsening 4<br class="">
Number of paths for aggressive coarsening 5<br class="">
Maximum row sums 0.9<br class="">
Sweeps down 1<br class="">
Sweeps up 1<br class="">
Sweeps on coarse 1<br class="">
Relax down sequential-Gauss-Seidel<br class="">
Relax up sequential-Gauss-Seidel<br class="">
Relax on coarse Gaussian-elimination<br class="">
Relax weight (all) 1.<br class="">
Outer relax weight (all) 1.<br class="">
Using CF-relaxation<br class="">
Smooth type Euclid<br class="">
Smooth num levels 25<br class="">
Euclid ILU(k) levels 0<br class="">
Euclid ILU(k) drop tolerance 0.<br class="">
Euclid ILU use Block-Jacobi? 1<br class="">
Measure type local<br class="">
Coarsen type HMIS<br class="">
Interpolation type ext+i<br class="">
linear system matrix = precond matrix:<br class="">
Mat Object: Mat_0x1e88040_0 1 MPI processes<br class="">
type: seqaij<br class="">
rows=153600, cols=153600, bs=24<br class="">
total: nonzeros=18247680, allocated nonzeros=18247680<br class="">
total number of mallocs used during MatSetValues calls =0<br class="">
using I-node routines: found 32000 nodes, limit used is 5<br class="">
Solution completed.<br class="">
Time elapsed 6.93336s<br class="">
Solution information...<br class="">
Error (L2) : 32.8359<br class="">
Error (Linf) : 1.93278<br class="">
************************************************************************************************************************<br class="">
*** WIDEN YOUR WINDOW TO 120 CHARACTERS. Use 'enscript -r -fCourier9' to print this document ***<br class="">
************************************************************************************************************************<br class="">
<br class="">
---------------------------------------------- PETSc Performance Summary: ----------------------------------------------<br class="">
<br class="">
bitpit_system_solver on a named xxx with 1 processor, by xxx Tue Jun 30 09:43:46 2020<br class="">
Using Petsc Release Version 3.10.3, unknown<br class="">
<br class="">
Max Max/Min Avg Total<br class="">
Time (sec): 9.718e+01 1.000 9.718e+01<br class="">
Objects: 5.300e+01 1.000 5.300e+01<br class="">
Flop: 1.107e+08 1.000 1.107e+08 1.107e+08<br class="">
Flop/sec: 1.139e+06 1.000 1.139e+06 1.139e+06<br class="">
MPI Messages: 0.000e+00 0.000 0.000e+00 0.000e+00<br class="">
MPI Message Lengths: 0.000e+00 0.000 0.000e+00 0.000e+00<br class="">
MPI Reductions: 0.000e+00 0.000<br class="">
<br class="">
Flop counting convention: 1 flop = 1 real number operation of type (multiply/divide/add/subtract)<br class="">
e.g., VecAXPY() for real vectors of length N --> 2N flop<br class="">
and VecAXPY() for complex vectors of length N --> 8N flop<br class="">
<br class="">
Summary of Stages: ----- Time ------ ----- Flop ------ --- Messages --- -- Message Lengths -- -- Reductions --<br class="">
Avg %Total Avg %Total Count %Total Avg %Total Count %Total<br class="">
0: Main Stage: 9.7178e+01 100.0% 1.1071e+08 100.0% 0.000e+00 0.0% 0.000e+00 0.0% 0.000e+00 0.0%<br class="">
<br class="">
------------------------------------------------------------------------------------------------------------------------<br class="">
See the 'Profiling' chapter of the users' manual for details on interpreting output.<br class="">
Phase summary info:<br class="">
Count: number of times phase was executed<br class="">
Time and Flop: Max - maximum over all processors<br class="">
Ratio - ratio of maximum to minimum over all processors<br class="">
Mess: number of messages sent<br class="">
AvgLen: average message length (bytes)<br class="">
Reduct: number of global reductions<br class="">
Global: entire computation<br class="">
Stage: stages of a computation. Set stages with PetscLogStagePush() and PetscLogStagePop().<br class="">
%T - percent time in this phase %F - percent flop in this phase<br class="">
%M - percent messages in this phase %L - percent message lengths in this phase<br class="">
%R - percent reductions in this phase<br class="">
Total Mflop/s: 10e-6 * (sum of flop over all processors)/(max time over all processors)<br class="">
------------------------------------------------------------------------------------------------------------------------<br class="">
Event Count Time (sec) Flop --- Global --- --- Stage ---- Total<br class="">
Max Ratio Max Ratio Max Ratio Mess AvgLen Reduct %T %F %M %L %R %T %F %M %L %R Mflop/s<br class="">
------------------------------------------------------------------------------------------------------------------------<br class="">
<br class="">
--- Event Stage 0: Main Stage<br class="">
<br class="">
BuildTwoSidedF 2 1.0 1.8471e-04 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0<br class="">
MatMult 3 1.0 6.8688e-02 1.0 1.09e+08 1.0 0.0e+00 0.0e+00 0.0e+00 0 98 0 0 0 0 98 0 0 0 1587<br class="">
MatConvert 2 1.0 3.8534e-01 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0<br class="">
MatAssemblyBegin 3 1.0 2.8112e-04 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0<br class="">
MatAssemblyEnd 3 1.0 7.8563e-02 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0<br class="">
MatGetRowIJ 2 1.0 3.6100e-06 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0<br class="">
MatView 2 1.0 5.4618e+01 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 56 0 0 0 0 56 0 0 0 0 0<br class="">
VecView 2 1.0 7.6981e-01 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 1 0 0 0 0 1 0 0 0 0 0<br class="">
VecMDot 1 1.0 1.7635e-04 1.0 3.07e+05 1.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 1742<br class="">
VecNorm 3 1.0 5.4832e-04 1.0 9.22e+05 1.0 0.0e+00 0.0e+00 0.0e+00 0 1 0 0 0 0 1 0 0 0 1681<br class="">
VecScale 1 1.0 1.5216e-04 1.0 1.54e+05 1.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 1009<br class="">
VecCopy 1 1.0 1.5406e-04 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0<br class="">
VecSet 21 1.0 7.1106e-03 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0<br class="">
VecAYPX 1 1.0 1.9940e-04 1.0 1.54e+05 1.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 770<br class="">
VecWAXPY 1 1.0 7.8802e-04 1.0 1.54e+05 1.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 195<br class="">
KSPSetUp 2 1.0 6.8995e-03 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0<br class="">
PCSetUp 2 1.0 1.3190e+01 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 14 0 0 0 0 14 0 0 0 0 0<br class="">
PCApply 1 1.0 2.2728e-01 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 0.0e+00 0 0 0 0 0 0 0 0 0 0 0<br class="">
------------------------------------------------------------------------------------------------------------------------<br class="">
<br class="">
Memory usage is given in bytes:<br class="">
<br class="">
Object Type Creations Destructions Memory Descendants' Mem.<br class="">
Reports information only for process 0.<br class="">
<br class="">
--- Event Stage 0: Main Stage<br class="">
<br class="">
Matrix 3 1 2872 0.<br class="">
Vector 32 20 17234800 0.<br class="">
Index Set 4 4 3200 0.<br class="">
IS L to G Mapping 2 0 0 0.<br class="">
Vec Scatter 2 0 0 0.<br class="">
Viewer 6 4 3392 0.<br class="">
Krylov Solver 2 1 61676 0.<br class="">
Preconditioner 2 1 1432 0.<br class="">
========================================================================================================================<br class="">
Average time to get PetscTime(): 3.36e-08<br class="">
#PETSc Option Table entries:<br class="">
-ksp_converged_reason<br class="">
-ksp_error_if_not_converged<br class="">
-ksp_monitor_singular_value<br class="">
-ksp_monitor_true_residual<br class="">
-log_view<br class="">
-pc_hypre_boomeramg_agg_nl 4<br class="">
-pc_hypre_boomeramg_agg_num_paths 5<br class="">
-pc_hypre_boomeramg_coarsen_type HMIS<br class="">
-pc_hypre_boomeramg_eu_bj true<br class="">
-pc_hypre_boomeramg_interp_type ext+i<br class="">
-pc_hypre_boomeramg_P_max 2<br class="">
-pc_hypre_boomeramg_relax_type_all sequential-Gauss-Seidel<br class="">
-pc_hypre_boomeramg_smooth_type Euclid<br class="">
-pc_hypre_boomeramg_strong_threshold 0.7<br class="">
-pc_hypre_boomeramg_truncfactor 0.3<br class="">
#End of PETSc Option Table entries<br class="">
Compiled without FORTRAN kernels<br class="">
Compiled with full precision matrices (default)<br class="">
sizeof(short) 2 sizeof(int) 4 sizeof(long) 8 sizeof(void*) 8 sizeof(PetscScalar) 8 sizeof(PetscInt) 4<br class="">
Configure options: PETSC_ARCH=arch-linux2-c-opt-gcc7-hypre --with-blaslapack-dir=/opt/lapack/3.8.0-gcc.7.4.0/lib64 --with-debugging=0 COPTFLAGS=-O3 CXXOPTFLAGS=-O3 FOPTFLAGS=-O3 --with-valgrind=1 --with-valgrind-dir=/opt/valgrind/3.14.0/ --prefix=/opt/petsc/3.10.3_gcc-7.4.0
--with-mpi=1 --download-hypre<br class="">
-----------------------------------------<br class="">
Libraries compiled on 2020-06-26 12:35:46 on xxx<br class="">
Machine characteristics: Linux-3.10.0-1127.8.2.el7.x86_64-x86_64-with-centos-7.6.1810-Core<br class="">
Using PETSc directory: /opt/petsc/3.10.3_gcc-7.4.0<br class="">
Using PETSc arch:<br class="">
-----------------------------------------<br class="">
<br class="">
Using C compiler: mpicc -fPIC -Wall -Wwrite-strings -Wno-strict-aliasing -Wno-unknown-pragmas -fstack-protector -fvisibility=hidden -O3<br class="">
Using Fortran compiler: mpif90 -fPIC -Wall -ffree-line-length-0 -Wno-unused-dummy-argument -O3<br class="">
-----------------------------------------<br class="">
<br class="">
Using include paths: -I/opt/petsc/3.10.3_gcc-7.4.0/include -I/opt/valgrind/3.14.0/include<br class="">
-----------------------------------------<br class="">
<br class="">
Using C linker: mpicc<br class="">
Using Fortran linker: mpif90<br class="">
Using libraries: -Wl,-rpath,/opt/petsc/3.10.3_gcc-7.4.0/lib -L/opt/petsc/3.10.3_gcc-7.4.0/lib -lpetsc -Wl,-rpath,/opt/petsc/3.10.3_gcc-7.4.0/lib -L/opt/petsc/3.10.3_gcc-7.4.0/lib -Wl,-rpath,/opt/lapack/3.8.0-gcc.7.4.0/lib64 -L/opt/lapack/3.8.0-gcc.7.4.0/lib64
-Wl,-rpath,/opt/openmpi4.0.0-gcc7.4.0_slurm/lib -L/opt/openmpi4.0.0-gcc7.4.0_slurm/lib -Wl,-rpath,/opt/gnu/7.4/lib/gcc/x86_64-pc-linux-gnu/7.4.0 -L/opt/gnu/7.4/lib/gcc/x86_64-pc-linux-gnu/7.4.0 -Wl,-rpath,/opt/gnu/7.4/lib/gcc -L/opt/gnu/7.4/lib/gcc -Wl,-rpath,/opt/gnu/7.4/lib64
-L/opt/gnu/7.4/lib64 -Wl,-rpath,/opt/gnu/7.4/lib -L/opt/gnu/7.4/lib -lHYPRE -llapack -lblas -lm -lX11 -lstdc++ -ldl -lmpi_usempif08 -lmpi_usempi_ignore_tkr -lmpi_mpifh -lmpi -lgfortran -lm -lgfortran -lm -lgcc_s -lquadmath -lpthread -lstdc++ -ldl<br class="">
-----------------------------------------</span></font><br class="">
</div>
</div>
</blockquote>
</div>
<br clear="all" class="">
<div class=""><br class="">
</div>
-- <br class="">
<div dir="ltr" class="">
<div dir="ltr" class="">
<div class="">
<div dir="ltr" class="">
<div class="">
<div dir="ltr" class="">
<div class="">What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br class="">
-- Norbert Wiener</div>
<div class=""><br class="">
</div>
<div class=""><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank" class="">https://www.cse.buffalo.edu/~knepley/</a><br class="">
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</blockquote>
</div>
</div>
</div>
</div>
</div></blockquote></div><br class=""></div></body></html>