[petsc-users] Use of PCFieldSplitSetSchurPre

Mark Adams mfadams at lbl.gov
Sat Aug 16 13:50:33 CDT 2014


I am creating an Schur complement matrix explicitly (solver%schur) and
setting it with code like:

call SNESGetKSP(solver%snes,innerksp,ierr)
call KSPGetPC(innerksp,spc,ierr)
call
PCFieldSplitSetSchurPre(spc,PC_FIELDSPLIT_SCHUR_PRE_USER,solver%schur,ierr)

I would think that the "inner" solver would not be used but it is getting
used.  It looks to me like I am not getting this hooked up correctly and my
Schur complement is not being used.  I've attached an output with view data
and converged_reason data.  I use the "inner" parameters for my Schur
construction (code appended), so you see it at the beginning of this
output.  But it is also in each solve, which seems wrong to me.

Any thoughts?
Thanks,
Mark


        call MatGetSize(solver%Bmat,M,N,ierr)
        call VecDuplicate(Xsub(1),v1,ierr)
        call VecDuplicate(Xsub(1),v2,ierr)
        call KSPCreate(solver%comm,innerksp,ierr)
        call
KSPSetOptionsPrefix(innerksp,'fsa_fieldsplit_lambda_inner_',ierr)
        call KSPSetFromOptions(innerksp,ierr)
        call KSPSetOperators(innerksp,  solver%Amat,  solver%Amat, ierr )
        call KSPSetUp(innerksp,ierr)
        call MatDuplicate(solver%Dmat,MAT_COPY_VALUES,solver%schur,ierr)
        call MatGetOwnershipRange(solver%Cmat,low,high,ierr)
        do k=0,N-1
           call MatGetColumnVector(solver%Bmat,v1,k,ierr)
           call KSPSolve(innerksp,v1,v2,ierr)
           do j=0,N-1
              call VecZeroEntries(v1,ierr)
              if (j.ge.low .and.j.lt.high) then
                 call MatGetRow(solver%Cmat,j,ncols,cols,vals,ierr)
                 call VecSetValues(v1,ncols,cols,vals,INSERT_VALUES,ierr)
                 !write(*,*) (vals(i), i=1,ncols)
                 call MatRestoreRow(solver%Cmat,j,ncols,cols,vals,ierr)
              end if
              call VecAssemblyBegin(v1,ierr)
              call VecAssemblyEnd(v1,ierr)
              call VecTDot(v1,v2,dot,ierr)
              dot = -dot
              call
MatSetValues(solver%schur,ione,j,ione,k,dot,INSERT_VALUES,ierr) ! D - C
(A-1) B
           end do
        end do
        call MatAssemblyBegin(solver%schur,MAT_FINAL_ASSEMBLY,ierr)
        call MatAssemblyEnd(solver%schur,MAT_FINAL_ASSEMBLY,ierr)
        call VecDestroy(v1,ierr)
        call VecDestroy(v2,ierr)
        call KSPDestroy(innerksp,ierr)

        call SNESGetKSP(solver%snes,innerksp,ierr)
        call KSPGetPC(innerksp,spc,ierr)
        call
PCFieldSplitSetSchurPre(spc,PC_FIELDSPLIT_SCHUR_PRE_USER,solver%schur,ierr)
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20140816/b469b723/attachment-0001.html>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: out
Type: application/octet-stream
Size: 85715 bytes
Desc: not available
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20140816/b469b723/attachment-0001.obj>


More information about the petsc-users mailing list