[petsc-users] GAMG speed

Matthew Knepley knepley at gmail.com
Tue Aug 13 20:34:55 CDT 2013


On Tue, Aug 13, 2013 at 8:25 PM, Michele Rosso <mrosso at uci.edu> wrote:

>  Matt,
>
> thank you! I will try to reduce the number of levels and see how it goes.
> I asked about the speed since CG + Block Jacobi with ICC in each block
> runs faster then CG + MG, so I thought I was missing something.
>

What is the operator? If its the mass matrix,  we expect this.


> Could you please tell me how to get rid of Chebichev?
>

I believe its

  -mg_levels_ksp_type richardson

but you can check.

  matt


> Best,
> Michele
>
>  On 08/13/2013 05:51 M, Matthew Knepley wrote:
>
> On Tue, Aug 13, 2013 at 7:05 PM, Michele Rosso <mrosso at uci.edu> wrote:
>
>>  Hi Matt,
>>
>> I attached the output of the commands you suggested.
>> The options I used are:
>>
>> -log_summary -ksp_monitor -ksp_view -ksp_converged_reason -pc_type mg
>> -pc_mg_galerkin -pc_mg_levels 5 -options_left
>>
>
>  The convergence is great. I notice that your coarse solve takes no time.
> You could probably use fewer levels for
> this problem. For this problem there is no easy things left I think. We
> are currently debating how you can squeeze
> something extra out of the smoother. Here you could probably get rid of
> Chebychev and use only SOR.
>
>      Matt
>
>
>>  and here are the lines of codes where I setup the solution process:
>>
>>     call DMDACreate3d( PETSC_COMM_WORLD
>> ,                                  &
>>                     & DMDA_BOUNDARY_PERIODIC ,
>> DMDA_BOUNDARY_PERIODIC,     &
>>                     & DMDA_BOUNDARY_PERIODIC ,
>> DMDA_STENCIL_STAR,          &
>>                     & N_Z , N_Y , N_X , N_B3 , N_B2 , 1_ip,  1_ip , 1_ip
>> , &
>>                     & NNZ ,NNY , NNX, da , ierr)
>>
>>
>>     ! Create Global Vectors
>>     call DMCreateGlobalVector(da,b,ierr)
>>     call VecDuplicate(b,x,ierr)
>>
>>     ! Set initial guess for first use of the module to 0
>>     call VecSet(x,0.0_rp,ierr)
>>
>>     ! Create matrix
>>     call DMCreateMatrix(da,MATAIJ,A,ierr)
>>
>>     ! Create solver
>>     call KSPCreate(PETSC_COMM_WORLD,ksp,ierr)
>>     call KSPSetDM(ksp,da,ierr)
>>     call KSPSetDMActive(ksp,PETSC_FALSE,ierr)
>>     call KSPSetOperators(ksp,A,A,SAME_NONZERO_PATTERN,ierr)
>>     call KSPSetType(ksp,KSPCG,ierr)
>>     call KSPSetNormType(ksp,KSP_NORM_UNPRECONDITIONED,ierr)
>>     call KSPSetInitialGuessNonzero(ksp,PETSC_TRUE,ierr)
>>     call KSPSetTolerances(ksp, tol ,PETSC_DEFAULT_DOUBLE_PRECISION,&
>>         & PETSC_DEFAULT_DOUBLE_PRECISION,PETSC_DEFAULT_INTEGER,ierr)
>>
>>     ! Nullspace removal
>>     call MatNullSpaceCreate(
>> PETSC_COMM_WORLD,PETSC_TRUE,PETSC_NULL_INTEGER,&
>>                            & PETSC_NULL_INTEGER,nullspace,ierr)
>>     call KSPSetNullspace(ksp,nullspace,ierr)
>>     call MatNullSpaceDestroy(nullspace,ierr)
>>
>>     ! To allow using option from command line
>>     call KSPSetFromOptions(ksp,ierr)
>>
>>
>> Hope I did not omit anything useful.
>> Thank you for your time.
>>
>> Best,
>> Michele
>>
>>
>>
>>
>>  On 08/13/2013 04:26 PM, Matthew Knepley wrote:
>>
>> On Tue, Aug 13, 2013 at 6:09 PM, Michele Rosso <mrosso at uci.edu> wrote:
>>
>>>  Hi Karli,
>>>
>>> thank you for your hint: now it works.
>>> Now I would like to speed up the solution: I was counting on increasing
>>> the number of levels/the number of processors used, but now I see I cannot
>>> do that.
>>> Do you have any hint to achieve better speed?
>>> Thanks!
>>>
>>
>>  "Better speed" is not very helpful for us, and thus we cannot offer
>> much help. You could
>>
>>   1) Send the output of -log_summary -ksp_monitor -ksp_view
>>
>>   2) Describe the operator succintly
>>
>>      Matt
>>
>>
>>>  Best,
>>> Michele
>>>
>>>
>>>
>>
>
>
>  --
> What most experimenters take for granted before they begin their
> experiments is infinitely more interesting than any results to which their
> experiments lead.
> -- Norbert Wiener
>
>
>


-- 
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20130813/bde3655f/attachment-0001.html>


More information about the petsc-users mailing list