[petsc-users] Issue on Block Jacobi Preconditioner Reuse
Mark Adams
mfadams at lbl.gov
Fri Oct 15 12:26:01 CDT 2021
You seem to be trying to use the setup of one solver (lu) for a different
solver (hypre).
You can't do that in general and PCSetType(subpc, PCHYPRE); will delete
the old solvers data.
On Thu, Oct 14, 2021 at 11:48 PM 王一甲 via petsc-users <
petsc-users at mcs.anl.gov> wrote:
> Hi!Everyone:
>
> Glad to join the mailing list.Recently I 've been working on a program
> using block jacobi preconditioner for a sequence solve of two linear system
> A_1x_1=b_1,A_2x_2=b_2.
>
> Since A1 and A2 have same nonzero pattern and their element values
> are quite close, we hope to reuse the preconditoners constructed when
> solving A_1x_1=b_1 ,
>
> however after calling KSPSetReusePreconditioner, though the iteration
> number of second ksp solve is very small but the time used did not decrease
> much, these are my code:
>
> ierr = KSPSetUp(ksp);CHKERRQ(ierr);
> ierr = PCBJacobiGetSubKSP(pc, &num_local, &idx_first_local,
> &subksp);CHKERRQ(ierr);
>
>
> for (i=0; i<num_local; i++)
> {
> ierr = KSPGetPC(subksp[i], &subpc);CHKERRQ(ierr);
> if (i==0)
> {
> ierr =
> KSPSetType(subksp[i],"gmres");CHKERRQ(ierr);
> ierr = PCSetType(subpc,
> PCLU);CHKERRQ(ierr);
> ierr = PCFactorSetMatSolverType(subpc,
> "mkl_pardiso");CHKERRQ(ierr);
> ierr = KSPSetOptionsPrefix(subksp[i],
> "Blk0_");CHKERRQ(ierr);
> ierr =
> KSPSetReusePreconditioner(subksp[i],PETSC_TRUE);CHKERRQ(ierr);
> ierr =
> KSPSetFromOptions(subksp[i]);CHKERRQ(ierr);
> }
> if (i==1)
> {
> ierr =
> KSPSetType(subksp[i],"gmres");CHKERRQ(ierr);
> ierr = PCSetType(subpc,
> PCHYPRE);CHKERRQ(ierr);
> ierr = PCHYPRESetType(pc,
> "boomeramg");CHKERRQ(ierr);
> ierr =
> KSPSetReusePreconditioner(subksp[i],PETSC_TRUE);CHKERRQ(ierr);
> ierr = KSPSetOptionsPrefix(subksp[i],
> "Blk1_");CHKERRQ(ierr);
> }
> }//nlocal
> }//isbjacobi
>
>
> //Solve the Linear System
> ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
> //Solve the Linear System
> t0 = MPI_Wtime();
> ierr = KSPSolve(ksp, b1, x1);CHKERRQ(ierr);
> t0 = MPI_Wtime()-t0;
> ierr = PetscPrintf(PETSC_COMM_SELF,"First KSP Solve time: %g
> s\n",t0);CHKERRQ(ierr);
> //Preconditioner Reuse
> ierr = KSPSetReusePreconditioner(ksp, PETSC_TRUE);CHKERRQ(ierr);
> ierr = KSPSetOperators(ksp,A2,A1);CHKERRQ(ierr);
> t0 = MPI_Wtime();
> ierr = KSPSolve(ksp,b2,x2);CHKERRQ(ierr);
> t0 = MPI_Wtime()-t0;
> ierr = PetscPrintf(PETSC_COMM_SELF,"Second KSP Solve time: %g
> s\n",t0);CHKERRQ(ierr);
>
>
> The total block number is 2, and the first block is solved using direct
> method and the second block is solved using hypre's boomeramg, so I hope to
> reuse the factor of the first block and the set up phase of boomeramg as
> preconditioner of the next solve, is there anything wrong with the reuse
> code?
>
>
> Best Wishes
>
>
> WANG Yijia
>
>
> 2021/10/15
>
>
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20211015/9474e3fd/attachment.html>
More information about the petsc-users
mailing list