[petsc-users] KSPMonitorSingularValue

Barry Smith bsmith at mcs.anl.gov
Fri Oct 21 08:19:00 CDT 2011


On Oct 21, 2011, at 2:48 AM, Klaij, Christiaan wrote:

>> If you just want to see the monitor, why not use the command line or
>> PetscOptionsSetValue()?
> 
> I have three ksp's: one for the momentum eqs (GMRES) and one for
> the pressure eq (CG) inside the SIMPLE preconditioner and one for
> the matrix-free coupled mass-momentum system (FGMRES).  The
> command line shows me results for all but right now I'm only
> interested in monitoring FGMRES on the coupled system.

  http://www.mcs.anl.gov/petsc/petsc-2/snapshots/petsc-current/docs/manualpages/KSP/KSPSetOptionsPrefix.html

   Give each ksp a unique prefix: for example -momentum   -pressure and -mass-momentum  Then use -mass-momentum_ksp_monitor_singular_value at the command line.

    Barry

> 
> 
> 
> 
>> This is the code for setting up the monitor, you can call the part inside
>> the if statement yourself if you like.
>> 
>>    ierr = PetscOptionsString("-ksp_monitor_singular_value","Monitor
>> singular
>> values","KSPMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
>>    if (flg) {
>>      ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE);CHKERRQ(ierr);
>>      ierr =
>> PetscViewerASCIIOpen(((PetscObject)ksp)->comm,monfilename,&monviewer);CHKERRQ(ierr);
>>      ierr =
>> KSPMonitorSet(ksp,KSPMonitorSingularValue,monviewer,(PetscErrorCode
>> (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
>>    }
> 
> I'm still confused whether it is supposed to work with FGMRES,
> the manual states only CG and GMRES. (Besides, I'm using fortran)
> 
> 
> 
>> How are you applying the action of the linear operator? If you use finite
>> differencing, it could be inaccurate. Is this incompressible or a low-Mach
>> compressible formulation? Try -ksp_monitor_true_residual, if the true
>> residual drifts from the unpreconditioned residual computed by FGMRES, the
>> Krylov space could be losing orthogonality. You can try
>> -ksp_gmres_modifiedgramschmidt. Are you losing a lot of progress in
>> restarts?
> 
> It's incompressible Navier-Stokes. No finite differencing, the
> action is computed directly without approximations. It's right
> preconditioning, so preconditioned and true residual should be
> the same. I don't get any progress, the residual is stagnating
> from the very first iteration way before any restart.
> 
> Regarding modified Gram Schmidt, I tried to set it as follows:
> 
> call KSPGMRESSetOrthogonalization(ksp,KSPGMRESModifiedGramSchmidtOrthogonalization,ierr)
> 
> But my compiler tells me:
> 
> This name does not have a type, and must have an explicit type.   [KSPGMRESMODIFIEDGRAMSCHMIDTORTHOGONALIZATIO]
> 
> (petsc-3.1-p7, fortran with "use petscksp" and #include "finclude/petsckspdef.h")
> 
> 
> dr. ir. Christiaan Klaij
> CFD Researcher
> Research & Development
> E mailto:C.Klaij at marin.nl
> T +31 317 49 33 44
> 
> MARIN
> 2, Haagsteeg, P.O. Box 28, 6700 AA Wageningen, The Netherlands
> T +31 317 49 39 11, F +31 317 49 32 45, I www.marin.nl
> 



More information about the petsc-users mailing list