[petsc-users] the specific preconditioner

Barry Smith bsmith at mcs.anl.gov
Fri Nov 21 15:38:38 CST 2014


> On Nov 21, 2014, at 2:53 PM, paul zhang <paulhuaizhang at gmail.com> wrote:
> 
> Hi Barry, 
> 
> Thanks for your checking. I know it for the first problem. 
> For the second, I am not quite sure how to use those flags though. I have compile all of my source code to a binary, say CFDsolver.  Then as I use it, I do 
> 
> $ mpirun -np 16 ../CFDsolver  case_file

  mpirun -np 16 ../CFDsolver case_file -ksp_type bcgs -pc_type jacobi -ksp_view 

for example

> 
> What am I supposed to do? 
> 
> Thanks,
> Paul
> 
> 
> 
> 
> 
> 
> 
> 
> 
> 
> On Fri, Nov 21, 2014 at 12:39 PM, Barry Smith <bsmith at mcs.anl.gov> wrote:
> 
>   You need to call KSPGetPC(ksp,&pc); before setting the pc
> 
>   Note that since you call KSPSetFromOptions(ksp); you can control the KSP's and PC's from the command line and don't need to edit the code and recompile it for each change to the KSP or PC.  For example
> 
>    -ksp_type bcgs -pc_type jacobi
> 
>   Barry
> 
> > On Nov 21, 2014, at 8:18 AM, huaibao zhang <paulhuaizhang at gmail.com> wrote:
> >
> > Hi Jed,
> >
> > I have a specific question regarding the attached my initialization function.  Previously I used  KSPSetType(ksp,KSPFGMRES) solver, and I would like to see if there is some improvement using the other choices like KSPBCGS. It seems hoverer, the PC         //PCSetType(pc,PCBJACOBI)  does  not work with me. And I am not quite sure why.
> >
> > By the way, can you please give some suggestions on my code? I am not quite sure if I am doing the things rightly.
> >
> >
> > Many thanks,
> > Paul
> >
> >
> >
> >
> > void petsc_init(void) {
> >
> >         vector<Cell>::iterator cit;
> >         vector<int>::iterator it;
> >
> >         //Create nonlinear solver context
> >         KSPCreate(PETSC_COMM_WORLD,&ksp);
> >
> >         VecCreateMPI(PETSC_COMM_WORLD,grid[gid].cellCount*nVars,grid[gid].globalCellCount*nVars,&RHS);
> >         VecSetFromOptions(RHS);
> >         VecDuplicate(RHS,&deltaX);
> >
> >
> >         VecSet(RHS,0.);
> >         VecSet(deltaX,0.);
> >
> >         vector<int> diagonal_nonzeros, off_diagonal_nonzeros;
> >         int nextCellCount;
> >
> >         // Calculate space necessary for matrix memory allocation
> >         for (cit=grid[gid].cell.begin();cit!=grid[gid].cell.end();cit++) {
> >                 nextCellCount=0;
> >                 for (it=(*cit).faces.begin();it!=(*cit).faces.end();it++) {
> >                         if (grid[gid].face[*it].bc==INTERNAL_FACE) {
> >                                 nextCellCount++;
> >                         }
> >                 }
> >                 for (int i=0;i<nVars;++i) {
> >                         diagonal_nonzeros.push_back( (nextCellCount+1)*nVars);
> >                         off_diagonal_nonzeros.push_back( ((*cit).ghosts.size())*nVars);
> >                 }
> >         }
> >
> >         MatCreateMPIAIJ(
> >                         PETSC_COMM_WORLD,
> >                         grid[gid].cellCount*nVars,
> >                         grid[gid].cellCount*nVars,
> >                         grid[gid].globalCellCount*nVars,
> >                         grid[gid].globalCellCount*nVars,
> >                         0,&diagonal_nonzeros[0],
> >                         0,&off_diagonal_nonzeros[0],
> >                         &impOP);
> >
> >         KSPSetOperators(ksp,impOP,impOP,SAME_NONZERO_PATTERN);
> >
> >         KSPSetTolerances(ksp,rtol,abstol,1.e15,maxits);
> >         KSPSetType(ksp, KSPBCGS);
> >         //PCSetType(pc,PCBJACOBI);
> >         KSPGMRESSetRestart(ksp,30);
> >         KSPSetFromOptions(ksp);
> >
> >         return;
> > }
> 
> 
> 
> 
> -- 
> Huaibao (Paul) Zhang
> Gas Surface Interactions Lab
> Department of Mechanical Engineering
> University of Kentucky,
> Lexington,
> KY, 40506-0503
> Office: 216 Ralph G. Anderson Building
> Web:gsil.engineering.uky.edu



More information about the petsc-users mailing list