[petsc-users] How to use multigrid?

Jed Brown jedbrown at mcs.anl.gov
Sat Nov 3 10:31:52 CDT 2012


1. *Always* send -log_summary when asking about performance.
2. AMG setup costs more, the solve should be faster, especially for large
problems.
3. 30k degrees of freedom is not large.


On Sat, Nov 3, 2012 at 10:27 AM, w_ang_temp <w_ang_temp at 163.com> wrote:

> Hello,
>     I have tried AMG, but there are some problems. I use the command:
>     mpiexec -n 4 ./ex4f -ksp_type gmres -pc_type gamg
> -pc_gamg_agg_nsmooths 1 -ksp_gmres_restart 170 -ksp_rtol 1.0e-15
> -ksp_converged_reason.
>     The matrix has a size of 30000. However, compared with -pc_type asm,
> the amg need more time:asm needs 4.9s, amg needs 13.7s. I did several tests
> and got the same conclusion. When it begins, the screen shows the
> information:
> [0]PCSetData_AGG bs=1 MM=7601. I do not know the meaning. And if there is
> some
> parameters that affect the performance of AMG?
>     Besides, I want to confirm a conception. In my view, AMG itself can be
> a solver
> like gmres. It can also be used as a preconditioner like jacobi and is
> used by combining
> with other solver. Is it right? If it is right, how use AMG solver?
>     My codes are attached.
>     Thanks.
>                                              Jim
>
> -----------codes------------------------------------------
>       call MatCreate(PETSC_COMM_WORLD,A,ierr)
>       call MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n,ierr)
>       call MatSetType(A, MATMPIAIJ,ierr)
>       call MatSetFromOptions(A,ierr)
> !premalloc
>       !find the max non-zero numbers of all rows
>        maxnonzero=0
>        do 19,II=1,m
>         !no-zero numbers of this row
>         maxnonzeroII=NROWIN(II+1)-NROWIN(II)
>         if (maxnonzeroII>maxnonzero) then
>         maxnonzero=maxnonzeroII
>         endif
>   19  continue
>      &nbs p; call
> MatMPIAIJSetPreallocation(A,maxnonzero,PETSC_NULL_INTEGER,
>      &  maxnonzero,PETSC_NULL_INTEGER,ierr)
>       call MatGetOwnershipRange(A,Istart,Iend,ierr)
> !set values per row
>        do 10,II=Istart+1,Iend
>         !no-zero numbers of this row
>         rowNum=NROWIN(II+1)-NROWIN(II)
>
>         allocate(nColPerRow(rowNum))
>         allocate(valuePerRow(rowNum))
>
>         kValStart=NROWIN(II)+1-1
>         kValEnd=NROWIN(II)+rowNum-1
>
>         !column index
>         nColPerRow=NNZIJ(kValStart:kValEnd)-1
>         valuePerRow=VALUE(kValStart:kValEnd)
> &nb sp;       nRow=II-1
>         call MatSetValues(A,ione,nRow,rowNum,nColPerRow,valuePerRow,
>      &       INSERT_VALUES,ierr)
>        deallocate(nColPerRow)
>        deallocate(valuePerRow)
>   10  continue
>       call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
>       call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)
>       call VecCreateMPI(PETSC_COMM_WORLD,PETSC_DECIDE,m,b,ierr)
>       call VecSetFromOptions(b,ierr)
>       call VecDuplicate(b,u,ierr)
>       call VecDuplicate(b,x,ierr)
>       call KSPCreate(PETSC_COMM_WORLD,ksp,ierr)
>       call KSPSetOperators(ksp,A,A,DIF FERENT_NONZERO_PATTERN,ierr)
>       call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-my_ksp_monitor',  &
>      &                    flg,ierr)
>       if (flg) then
>         call KSPMonitorSet(ksp,MyKSPMonitor,PETSC_NULL_OBJECT,          &
>      &                     PETSC_NULL_FUNCTION,ierr)
>       endif
>       call KSPSetFromOptions(ksp,ierr)
>       call PetscOptionsHasName(PETSC_NULL_CHARACTER,                    &
>      &     '-my_ksp_convergence',flg,ierr)
>       if (flg) then
>         call KSPSetConvergenceTest(ksp,MyKSPConverged,                  &
>      &          PETSC_NULL_OBJECT,PETSC_NULL_FUNCTION,ierr)
>       endif
>       !Assing values to 'b'
>       bTemp=Iend-Istart
>       ioneb=1
>
>       do 12,II=Istart,Iend-1
>             voneb=F(II+1)
>          call VecSetValues(b,ioneb,II,voneb,INSERT_VALUES,ierr)
>   12    continue
>       call VecAssemblyBegin(b,ierr)
>       call VecAssemblyEnd(b,ierr)
>       call KSPSolve(ksp,b,x,ierr)
> ------codes-------------------------------------------------
>
>
>
>
>
>
>
>
>
> >At 2012-11-01 22:00:28,"Jed Brown" <jedbrown at mcs.anl.gov> wrote:
>
> >Yes, it's faster to understand this error message than to have
> "mysteriously slow performance".
>
> >** Preallocation routines now automatically set
> MAT_NEW_NONZERO_ALLOCATION_ERR, if you intentionally preallocate less than
> necessary then *>*use
> MatSetOption(mat,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE) to disable the
> error generation.
> *
> >http://www.mcs.anl.gov/petsc/documentation/changes/33.html
>
> >On Thu, Nov 1, 2012 at 8:57 AM, w_ang_temp <w_ang_temp at 163.com> wrote:
>
>> >Do you mean that the two versions have a different in this point? If I
>> use the new version, I have to
>> >make some modifications on my codes?
>>
>
>
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20121103/c4ad41d3/attachment.html>


More information about the petsc-users mailing list