[petsc-dev] sor smoothers

Jed Brown jedbrown at mcs.anl.gov
Mon Aug 12 22:23:27 CDT 2013


Barry Smith <bsmith at mcs.anl.gov> writes:

>   Jed,
>
>    I don't understand your computations. I run
>
> ./ex29 -ksp_type chebyshev -ksp_max_it 2 -ksp_monitor -pc_type eisenstat -log_summary -ksp_view -ksp_chebyshev_eigenvalues 1.0,2.0 -ksp_norm_type none
>
>
> 85	  ierr =   MatSOR(eis->A,eis->b[pc->presolvedone-1],eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_FORWARD_SWEEP),0.0,1,1,b);CHKERRQ(ierr);
> Breakpoint 2, MatSOR (mat=0x7f91ac1a7e70, b=0x7f91ac1cf470, omega=1, flag=SOR_EISENSTAT, shift=0, its=1, lits=1, x=0x7f91ac1db070) at matrix.c:3681
> Breakpoint 2, MatSOR (mat=0x7f91ac1a7e70, b=0x7f91ac1d5270, omega=1, flag=SOR_EISENSTAT, shift=0, its=1, lits=1, 
> 102	  ierr = MatSOR(eis->A,eis->b[pc->presolvedone],eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_BACKWARD_SWEEP),0.0,1,1,x);CHKERRQ(ierr);
>
> so get .5 + 1 + 1 + .5 = 3 work units (where a work unit is a matmult).  

I miscounted the middle two MatSOR.  A normal SSOR with nonzero initial
guess is two work units, but Eisenstat is indeed only one work unit.

> Now I run with 
>
> ./ex29 -ksp_type chebyshev -ksp_max_it 2 -ksp_monitor -pc_type sor  -log_summary -ksp_view -ksp_chebyshev_eigenvalues 1.0,2.0 -ksp_norm_type none 
>
> #1  0x000000010232ca57 in PCApply_SOR (pc=0x7ffb4697a470, x=0x7ffb469dc470, y=0x7ffb469d0870) at sor.c:35
> 35	  ierr = MatSOR(pc->pmat,x,jac->omega,(MatSORType)flag,jac->fshift,jac->its,jac->lits,y);CHKERRQ(ierr);
> (gdb) p flag
> $1 = 28
>
> Breakpoint 2, MatMult (mat=0x7ffb469a9270, x=0x7ffb469d0870, y=0x7ffb469dc470) at matrix.c:2212
> Breakpoint 1, MatSOR (mat=0x7ffb469a9270, b=0x7ffb469dc470, omega=1, flag=28, shift=0, its=1, lits=1, x=0x7ffb469d6670) at 
> Breakpoint 2, MatMult (mat=0x7ffb469a9270, x=0x7ffb469d6670, y=0x7ffb469dc470) at matrix.c:2212
> Breakpoint 1, MatSOR (mat=0x7ffb469a9270, b=0x7ffb469dc470, omega=1, flag=28, shift=0, its=1, lits=1, x=0x7ffb469a1670) at 
>
> so get (.5 + .5) + 1 + (.5 + .5) + 1 + (.5 + .5) = 5 work units
>
> Note that for the first iteration with a zero initial guess the down
> AND UP sweep together can be done in one work unit; the values in the
> downsweep are saved and used to save work in the up sweep.

Yup, my mistake.  I was counting for all those SSORs having a nonzero
initial guess.

>
> MatMult                2 1.0 1.1801e-02 1.0 1.16e+07 1.0 0.0e+00 0.0e+00 0.0e+00  0 10  0  0  0   0 10  0  0  0   981
> MatSOR                 3 1.0 4.6818e-02 1.0 1.78e+07 1.0 0.0e+00 0.0e+00 0.0e+00  0 16  0  0  0   0 16  0  0  0   380
>
> Thus we see that we save all of the MatMult time which is 2 units of the 5 units needed with SOR in terms of flops computed so 40% of the work but only 20% of the time.
>
> On the post-smooth of the multigrid there is a nonzero initial guess eisenstat does 
>
>   if (nonzero) {
>     ierr = VecCopy(x,eis->b[pc->presolvedone-1]);CHKERRQ(ierr);
>     ierr = MatSOR(eis->A,eis->b[pc->presolvedone-1],eis->omega,SOR_APPLY_UPPER,0.0,1,1,x);CHKERRQ(ierr);
>
> so an extra .5 work unit 
>
> while Chebychev does the matrix vector product to get the initial residual so 
>
> Eisenstat is 3 units + .5 unit + 1 unit = 4.5 units
> SOR           5  units               + 1 unit = 6 units 
>
> so for combined pre and post smooth Eisenstat/SOR = 7.5/11 work units

I think that is right, and indeed, that looks like enough benefit to
justify converting the matrix format.
-------------- next part --------------
A non-text attachment was scrubbed...
Name: not available
Type: application/pgp-signature
Size: 835 bytes
Desc: not available
URL: <http://lists.mcs.anl.gov/pipermail/petsc-dev/attachments/20130812/e223c20c/attachment.sig>


More information about the petsc-dev mailing list