<html>
<head>
<meta http-equiv="Content-Type" content="text/html; charset=utf-8">
</head>
<body style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;" class="">
Well yes naturally for the residual but adding -ksp_true_residual just gives
<div class=""><span class=""><br class="">
<font size="1" face="Courier" class="">  0 KSP unpreconditioned resid norm 3.583290589961e+00 true resid norm 3.583290589961e+00 ||r(i)||/||b|| 1.000000000000e+00<br class="">
  1 KSP unpreconditioned resid norm 0.000000000000e+00 true resid norm 3.583290589961e+00 ||r(i)||/||b|| 1.000000000000e+00<br class="">
Linear solve converged due to CONVERGED_ATOL iterations 1</font></span></div>
<div class=""><span class=""><font face="Courier New" size="1" class=""><br class="">
</font></span></div>
<div class="">Mark - if that helps - a Poisson equation is used for the pressure so the Helmholtz is the same as for the velocity in the interior.</div>
<div class=""><br class="">
</div>
<div class="">Thibaut</div>
<div class=""><span class=""><font face="Courier New" size="1" class=""><br class="">
</font>
<blockquote type="cite" class="">Le 31 oct. 2018 à 21:05, Mark Adams <<a href="mailto:mfadams@lbl.gov" class="">mfadams@lbl.gov</a>> a écrit :<br class="">
<br class="">
These are indefinite (bad) Helmholtz problems. Right?<br class="">
<br class="">
On Wed, Oct 31, 2018 at 2:38 PM Matthew Knepley <<a href="mailto:knepley@gmail.com" class="">knepley@gmail.com</a>> wrote:<br class="">
On Wed, Oct 31, 2018 at 2:13 PM Thibaut Appel <t.appel17@imperial.ac.uk> wrote:<br class="">
Hi Mark, Matthew,<br class="">
<br class="">
Thanks for taking the time.<br class="">
<br class="">
1) You're not suggesting having -fieldsplit_X_ksp_type fgmres for each field, are you?<br class="">
<br class="">
2) No, the matrix has pressure in one of the fields. Here it's a 2D problem (but we're also doing 3D), the unknowns are (p,u,v) and those are my 3 fields. We are dealing with subsonic/transsonic flows so it is convection dominated indeed.<br class="">
<br class="">
3) We are in frequency domain with respect to time, i.e. \partial{phi}/\partial{t} = -i*omega*phi.<br class="">
<br class="">
4) Hypre is unfortunately not an option since we are in complex arithmetic.<br class="">
<br class="">
<br class="">
<br class="">
<blockquote type="cite" class="">I'm not sure about "-fieldsplit_pc_type gamg" GAMG should work on one block, and hence be a subpc. I'm not up on fieldsplit syntax.<br class="">
</blockquote>
According to the online manual page this syntax applies the suffix to all the defined fields?<br class="">
<br class="">
<br class="">
<br class="">
<blockquote type="cite" class="">Mark is correct. I wanted you to change the smoother. He shows how to change it to Richardson (make sure you add the self-scale option), which is probably the best choice.<br class="">
<br class="">
  Thanks,<br class="">
<br class="">
     Matt<br class="">
</blockquote>
<br class="">
You did tell me to set it to GMRES if I'm not mistaken, that's why I tried "-fieldsplit_mg_levels_ksp_type gmres" (mentioned in the email). Also, it wasn't clear whether these should be applied to each block or the whole system, as the online manual pages + .pdf
 manual barely mention smoothers and how to manipulate MG objects with KSP/PC, this especially with PCFIELDSPLIT where examples are scarce.<br class="">
<br class="">
From what I can gather from your suggestions I tried (lines with X are repeated for X={0,1,2}) <br class="">
<br class="">
This looks good. How can an identically zero vector produce a 0 residual? You should always monitor with<br class="">
<br class="">
  -ksp_monitor_true_residual.<br class="">
<br class="">
   Thanks,<br class="">
 <br class="">
    Matt <br class="">
-ksp_view_pre -ksp_monitor -ksp_converged_reason \<br class="">
-ksp_type fgmres -ksp_rtol 1.0e-8 \<br class="">
-pc_type fieldsplit \<br class="">
-pc_fieldsplit_type multiplicative \<br class="">
-pc_fieldsplit_block_size 3 \<br class="">
-pc_fieldsplit_0_fields 0 \<br class="">
-pc_fieldsplit_1_fields 1 \<br class="">
-pc_fieldsplit_2_fields 2 \<br class="">
-fieldsplit_X_pc_type gamg \<br class="">
-fieldsplit_X_ksp_type gmres \<br class="">
-fieldsplit_X_ksp_rtol 1e-10 \<br class="">
-fieldsplit_X_mg_levels_ksp_type richardson \<br class="">
-fieldsplit_X_mg_levels_pc_type sor \<br class="">
-fieldsplit_X_pc_gamg_agg_nsmooths 0 \<br class="">
-fieldsplit_X_mg_levels_ksp_richardson_self_scale \<br class="">
-log_view<br class="">
<br class="">
which yields <br class="">
<br class="">
KSP Object: 1 MPI processes<br class="">
  type: fgmres<br class="">
    restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement<br class="">
    happy breakdown tolerance 1e-30<br class="">
  maximum iterations=10000, initial guess is zero<br class="">
  tolerances:  relative=1e-08, absolute=1e-50, divergence=10000.<br class="">
  left preconditioning<br class="">
  using DEFAULT norm type for convergence test<br class="">
PC Object: 1 MPI processes<br class="">
  type: fieldsplit<br class="">
  PC has not been set up so information may be incomplete<br class="">
    FieldSplit with MULTIPLICATIVE composition: total splits = 3, blocksize = 3<br class="">
    Solver info for each split is in the following KSP objects:<br class="">
  Split number 0 Fields  0<br class="">
  KSP Object: (fieldsplit_0_) 1 MPI processes<br class="">
    type: preonly<br class="">
    maximum iterations=10000, initial guess is zero<br class="">
    tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.<br class="">
    left preconditioning<br class="">
    using DEFAULT norm type for convergence test<br class="">
  PC Object: (fieldsplit_0_) 1 MPI processes<br class="">
    type not yet set<br class="">
    PC has not been set up so information may be incomplete<br class="">
  Split number 1 Fields  1<br class="">
  KSP Object: (fieldsplit_1_) 1 MPI processes<br class="">
    type: preonly<br class="">
    maximum iterations=10000, initial guess is zero<br class="">
    tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.<br class="">
    left preconditioning<br class="">
    using DEFAULT norm type for convergence test<br class="">
  PC Object: (fieldsplit_1_) 1 MPI processes<br class="">
    type not yet set<br class="">
    PC has not been set up so information may be incomplete<br class="">
  Split number 2 Fields  2<br class="">
  KSP Object: (fieldsplit_2_) 1 MPI processes<br class="">
    type: preonly<br class="">
    maximum iterations=10000, initial guess is zero<br class="">
    tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.<br class="">
    left preconditioning<br class="">
    using DEFAULT norm type for convergence test<br class="">
  PC Object: (fieldsplit_2_) 1 MPI processes<br class="">
    type not yet set<br class="">
    PC has not been set up so information may be incomplete<br class="">
  linear system matrix = precond matrix:<br class="">
  Mat Object: 1 MPI processes<br class="">
    type: seqaij<br class="">
    rows=52500, cols=52500<br class="">
    total: nonzeros=1127079, allocated nonzeros=1128624<br class="">
    total number of mallocs used during MatSetValues calls =0<br class="">
      not using I-node routines<br class="">
  0 KSP Residual norm 3.583290589961e+00 <br class="">
  1 KSP Residual norm 0.000000000000e+00 <br class="">
Linear solve converged due to CONVERGED_ATOL iterations 1<br class="">
<br class="">
so something must not be set correctly. The solution is identically zero everywhere.<br class="">
<br class="">
Is that option list what you meant? If you could let me know what should be corrected.<br class="">
<br class="">
<br class="">
<br class="">
Thanks for your support,<br class="">
<br class="">
<br class="">
<br class="">
Thibaut<br class="">
<br class="">
<br class="">
<br class="">
On 31/10/2018 16:43, Mark Adams wrote:<br class="">
<blockquote type="cite" class=""><br class="">
<br class="">
On Tue, Oct 30, 2018 at 5:23 PM Appel, Thibaut via petsc-users <petsc-users@mcs.anl.gov> wrote:<br class="">
Dear users,<br class="">
<br class="">
Following a suggestion from Matthew Knepley I’ve been trying to apply fieldsplit/gamg for my set of PDEs but I’m still encountering issues despite various tests. pc_gamg simply won’t start.<br class="">
Note that direct solvers always yield the correct, physical result.<br class="">
Removing the fieldsplit to focus on the gamg bit and trying to solve the linear system on a modest size problem still gives, with<br class="">
<br class="">
'-ksp_monitor -ksp_rtol 1.0e-10 -ksp_gmres_restart 300 -ksp_type gmres -pc_type gamg'<br class="">
<br class="">
[3]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------<br class="">
[3]PETSC ERROR: Petsc has generated inconsistent data<br class="">
[3]PETSC ERROR: Have un-symmetric graph (apparently). Use '-(null)pc_gamg_sym_graph true' to symetrize the graph or '-(null)pc_gamg_threshold -1' if the matrix is structurally symmetric.<br class="">
<br class="">
And since then, after adding '-pc_gamg_sym_graph true' I have been getting<br class="">
[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------<br class="">
[0]PETSC ERROR: Petsc has generated inconsistent data<br class="">
[0]PETSC ERROR: Eigen estimator failed: DIVERGED_NANORINF at iteration<br class="">
<br class="">
-ksp_chebyshev_esteig_noisy 0/1 does not change anything<br class="">
<br class="">
Knowing that Chebyshev eigen estimator needs a positive spectrum I tried ‘-mg_levels_ksp_type gmres’ but iterations would just go on endlessly.<br class="">
<br class="">
This is OK, but you need to use '-ksp_type fgmres' (this could be why it is failing ...). <br class="">
<br class="">
It looks like your matrix is 1) just the velocity field and 2) very unsymmetric (eg, convection dominated). I would start with ‘-mg_levels_ksp_type richardson -mg_levels_pc_type sor’.<br class="">
<br class="">
I would also start with unsmoothed aggregation: '-pc_gamg_nsmooths 0' <br class="">
 <br class="">
<br class="">
It seems that I have indeed eigenvalues of rather high magnitude in the spectrum of my operator without being able to determine the reason.<br class="">
The eigenvectors look like small artifacts at the wall-inflow or wall-outflow corners with zero anywhere else but I do not know how to interpret this.<br class="">
Equations are time-harmonic linearized Navier-Stokes to which a forcing is applied, there’s no time-marching.<br class="">
<br class="">
You mean you are in frequency domain?<br class="">
 <br class="">
<br class="">
Matrix is formed with a MPIAIJ type. The formulation is incompressible, in complex arithmetic and the 2D physical domain is mapped to a logically rectangular,<br class="">
<br class="">
This kind of messes up the null space that AMG depends on but AMG theory is gone for NS anyway.<br class="">
 <br class="">
regular collocated grid with a high-order finite difference method.<br class="">
I determine the ownership of the rows/degrees of freedom of the matrix with PetscSplitOwnership and I’m not using DMDA.<br class="">
<br class="">
Our iterative solvers are probably not going to work well on this but you should test hypre also (-pc_type hypre -pc_hypre_type boomeramg). You need to configure PETSc to download hypre.<br class="">
<br class="">
Mark<br class="">
 <br class="">
<br class="">
The Fortran application code is memory-leak free and has undergone a strict verification/validation procedure for different variations of the PDEs.<br class="">
<br class="">
If there’s any problem with the matrix what could help for the diagnostic? At this point I’m running out of ideas so I would really appreciate additional suggestions and discussions.<br class="">
<br class="">
Thanks for your continued support,<br class="">
<br class="">
<br class="">
Thibaut<br class="">
</blockquote>
<br class="">
<br class="">
-- <br 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<br class="">
<br class="">
https://www.cse.buffalo.edu/~knepley/<br class="">
</blockquote>
<br class="">
</span></div>
</body>
</html>