[petsc-users] how to use schur fieldsplit preconditioner twice ~
Liang ZHENG
liangzh.cug at gmail.com
Wed Mar 9 15:23:02 CST 2016
Hi, all,
I wrote a code as following:
void petsc_solve_thm2d(type_matrix_petsc *matrix)
{
IS isg[4],iss,isf;
MSGS("[petsc] Set blocked preconditioner index");
PetscScalar np_Sv=matrix->np_col/3;
PetscScalar np_Sp=matrix->np_col/6;
PetscScalar np_Fv=matrix->np_col/3;
PetscScalar np_Fp=matrix->np_col/6;
ISCreateStride(PETSC_COMM_SELF,np_Sv,0,1,&isg[0]);
ISCreateStride(PETSC_COMM_SELF,np_Sp,np_Sv,1,&isg[1]);
ISCreateStride(PETSC_COMM_SELF,np_Fv,np_Sv+np_Sp,1,&isg[2]);
ISCreateStride(PETSC_COMM_SELF,np_Fp,np_Sv+np_Sp+np_Fv,1,&isg[3]);
IS is_solid[]={isg[0],isg[1]};
IS is_fluid[]={isg[2],isg[3]};
ierr = ISConcatenate(PETSC_COMM_SELF,2,is_solid,&iss);CHKERRV(ierr);
ierr = ISConcatenate(PETSC_COMM_SELF,2,is_fluid,&isf);CHKERRV(ierr);
MSGS("[petsc] Create solver");
ierr = KSPCreate(PETSC_COMM_SELF,&matrix->ksp);CHKERRV(ierr);
ierr = KSPSetOptionsPrefix(matrix->ksp,"thm_");CHKERRV(ierr);
ierr = KSPSetOperators(matrix->ksp,matrix->LHS,matrix->LHS);CHKERRV(ierr);
ierr = KSPSetFromOptions(matrix->ksp);CHKERRV(ierr);
ierr = KSPGetPC(matrix->ksp,&matrix->pc);CHKERRV(ierr);
// ierr = PCFieldSplitSetIS(matrix->pc,"Sv",isg[0]);CHKERRV(ierr);
// ierr = PCFieldSplitSetIS(matrix->pc,"Sp",isg[1]);CHKERRV(ierr);
// ierr = PCFieldSplitSetIS(matrix->pc,"Fv",isg[2]);CHKERRV(ierr);
// ierr = PCFieldSplitSetIS(matrix->pc,"Fp",isg[3]);CHKERRV(ierr);
ierr = PCFieldSplitSetIS(matrix->pc,"S",iss);CHKERRV(ierr);
ierr = PCFieldSplitSetIS(matrix->pc,"F",isf);CHKERRV(ierr);
//ierr = PCFieldSplitSetBlockSize(matrix->pc,4);CHKERRV(ierr);
MSGS("[petsc] Solving");
ierr = KSPSolve(matrix->ksp,matrix->RHS,matrix->SOL);CHKERRV(ierr);
MSGS("[petsc] Destroy solver");
ierr = KSPDestroy(&matrix->ksp);CHKERRV(ierr);
ierr = ISDestroy(&isg);CHKERRV(ierr);
}
then I can call petsc in the command line as:
./thm2d_nonuniform -thm_ksp_monitor -thm_ksp_view -thm_ksp_type fgmres
-thm_ksp_rtol 1e-10 -thm_pc_type fieldsplit -thm_pc_fieldsplit_type schur
-thm_pc_fieldsplit_S_ksp_type fgmres -thm_fieldsplit_S_pc_type lu
-thm_pc_fieldsplit_F_ksp_type fgmres -thm_fieldsplit_F_pc_type lu
-thm_fieldsplit_S_pc_factor_mat_solver_package umfpack
this command is working. but I also want to set -thm_fieldsplit_S_pc_type
using fieldsplit schur preconditioner. Is there a solution to do something
like:
-thm_fieldsplit_S_pc_type fieldsplit -thm_fieldsplit_S_pc_field_type schur
-thm_fieldsplit_S_fieldsplit_Sv_ksp_type xxx
-thm_fieldsplit_S_fieldsplit_Sp_ksp_type xxx
-thm_fieldsplit_S_fieldsplit_Sv_pc_type xxx -
thm_fieldsplit_S_fieldsplit_Sp_pc_type xxx
Thanks a lot ~
Cheers,
Larry
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20160309/82f86fee/attachment.html>
More information about the petsc-users
mailing list