<div dir="ltr"><div>Hello Barry, Matthew, thanks for the replies !</div><div><br></div><div>Yes, it is our custom code, and it also happens when setting -pc_type bjacobi. Before testing an iterative solver, we were using MUMPS (-ksp_type preonly -ksp_pc_type lu -pc_factor_mat_solver_type mumps) without issues.</div><div><br></div><div>Running the ex19 (as "mpirun -n 4 ex19 -da_refine 5") did not produce any problem.<br><br></div><div>To reproduce the situation on my computer, I was able to reproduce the error for a small case and -pc_type bjacobi. For that particular case, when running in the cluster the error appears at the very last iteration:</div><div><br></div><div>=====<br></div><div> 27 KSP Residual norm 8.230378644666e-06 <br>[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------<br>[0]PETSC ERROR: Invalid argument<br>[0]PETSC ERROR: Scalar value must be same on all processes, argument # 3</div><div>====</div><div><br></div><div>whereas running on my computer the error is not launched and convergence is reached instead:</div><div><br></div><div>====<br>Linear interp_ solve converged due to CONVERGED_RTOL iterations 27</div><div>====</div><div><br></div><div>I will run valgrind to seek for possible memory corruptions.</div><div><br></div><div>thank you<br></div><div>Alfredo<br></div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Mon, Aug 24, 2020 at 9:00 PM Barry Smith <<a href="mailto:bsmith@petsc.dev" target="_blank">bsmith@petsc.dev</a>> wrote:<br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><div dir="auto"><div><br></div>   Oh yes, it could happen with Nan. <div><br></div><div>   KSPGMRESClassicalGramSchmidtOrthogonalization() calls  KSPCheckDot(ksp,lhh[j]); so should detect any NAN that appear and set ksp->convergedreason  but the call to MAXPY() is still made before returning and hence producing the error message.</div><div><br></div><div>   We should circuit the orthogonalization as soon as it sees a Nan/Inf and return immediately for GMRES to cleanup and produce a very useful error message. </div><div><br></div><div>  Alfredo,</div><div><br></div><div>    It is also possible that the hypre preconditioners are producing a Nan because your matrix is too difficult for them to handle, but it would be odd to happen after many iterations.</div><div><br></div><div>   As I suggested before run with -pc_type bjacobi to see if you get the same problem.</div><div><br></div><div>  Barry</div><div><br></div><div><div><br><blockquote type="cite"><div>On Aug 24, 2020, at 6:38 PM, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>> wrote:</div><br><div><div dir="ltr"><div dir="ltr">On Mon, Aug 24, 2020 at 6:27 PM Barry Smith <<a href="mailto:bsmith@petsc.dev" target="_blank">bsmith@petsc.dev</a>> wrote:<br></div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><br>
   Alfredo,<br>
<br>
      This should never happen. The input to the VecMAXPY in gmres is computed via VMDot which produces the same result on all processes.<br>
<br>
       If you run with -pc_type bjacobi does it also happen?<br>
<br>
       Is this your custom code or does it happen in PETSc examples also? Like src/snes/tutorials/ex19 -da_refine 5 <br>
<br>
      Could be memory corruption, can you run under valgrind?<br></blockquote><div><br></div><div>Couldn't it happen if something generates a NaN? That also should not happen, but I was allowing that pilut might do it.</div><div><br></div><div>  Thanks,</div><div><br></div><div>    Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex">
    Barry<br>
<br>
<br>
> On Aug 24, 2020, at 4:05 PM, Alfredo Jaramillo <<a href="mailto:ajaramillopalma@gmail.com" target="_blank">ajaramillopalma@gmail.com</a>> wrote:<br>
> <br>
> Dear PETSc developers,<br>
> <br>
> I'm trying to solve a linear problem with GMRES preconditioned with pilut from HYPRE. For this I'm using the options:<br>
> <br>
> -ksp_type gmres -pc_type hypre -pc_hypre_type pilut -ksp_monitor<br>
> <br>
> If I use a single core, GMRES (+ pilut or euclid) converges. However, when using multiple cores the next error appears after some number of iterations:<br>
> <br>
> [0]PETSC ERROR: Scalar value must be same on all processes, argument # 3<br>
> <br>
> relative to the function VecMAXPY. I attached a screenshot with more detailed output. The same happens when using euclid. Can you please give me some insight on this?<br>
> <br>
> best regards<br>
> Alfredo<br>
> <Screenshot from 2020-08-24 17-57-52.png><br>
<br>
</blockquote></div><br clear="all"><div><br></div>-- <br><div dir="ltr"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>-- Norbert Wiener</div><div><br></div><div><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>
</div></blockquote></div><br></div></div></div></blockquote></div>