<html>
<head>
<meta http-equiv="Content-Type" content="text/html; charset=utf-8">
</head>
<body style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;" class="">
Roman,<br class="">
<br class="">
When I tried what you suggested, I get the following error:<br class="">
<br class="">
 <font face="Courier" class="">      call PCSetOperators(pc,PETSC_NULL_MAT,my_matmul_k,ierr);CHKERRA(ierr) <br class="">
                                            1<br class="">
Error: Invalid procedure argument at (1)<br class="">
make: *** [/Users/hall/git/Mustang/Mustang/build_dir/slepc_eigen2.o] Error 1<br class="">
<br class="">
</font><br class="">
I assume this is because PCSetOperators is expecting a MAT type rather than a subroutine.  I modified your code snippet using a shell matrix, with the result:<br class="">
<br class="">
call EPSSetWhichEigenpairs(eps,EPS_TARGET_REAL,ierr);CHKERRA(ierr)<br class="">
call EPSSetTarget(eps,sigma,ierr);CHKERRA(ierr)<br class="">
call EPSGetST(eps,st,ierr);CHKERRA(ierr)<br class="">
call STSetType(st,STSINVERT,ierr);CHKERRA(ierr)  <br class="">
call STGetKSP(st,ksp,ierr);CHKERRA(ierr)      <br class="">
call KSPGetPC(ksp,pc,ierr);CHKERRA(ierr)      <br class="">
call KSPSetType(ksp,KSPBCGS,ierr);CHKERRA(ierr)    <br class="">
call PCSetType(pc,PCNONE,ierr);CHKERRA(ierr)     <br class="">
call EPSSetUp(eps,ierr);CHKERRA(ierr)      <br class="">
<br class="">
<div class="">! begin modified code</div>
<div class=""><br class="">
call MatCreateShell(PETSC_COMM_WORLD,n,n,n,n,PETSC_NULL_INTEGER,P,ierr) ; CHKERRA(ierr)<br class="">
call MatShellSetOperation(P,MATOP_MULT,my_matmul_k,ierr); CHKERRA(ierr)<br class="">
call PCSetOperators(pc,PETSC_NULL_MAT,P,ierr);CHKERRA(ierr)</div>
<div class=""><br class="">
! end modified code<br class="">
<br class="">
</div>
<div class="">call PCSetType(pc,PCMAT,ierr);CHKERRA(ierr)<br class="">
<br class="">
<br class="">
<br class="">
</div>
<div class="">Is this the correct thing to do here?</div>
<div class="">Again, thank you for your help.<br class="">
<div class="">Kenneth<br class="">
<br class="">
<br class="">
<blockquote type="cite" class="">On Aug 7, 2018, at 11:13 AM, Jose E. Roman <<a href="mailto:jroman@dsic.upv.es" class="">jroman@dsic.upv.es</a>> wrote:<br class="">
<br class="">
Currently there is no straightforward way to do this in SLEPc. For a preconditioned eigensolver such as JD you could use this:<br class="">
<a href="https://urldefense.proofpoint.com/v2/url?u=http-" class="">https://urldefense.proofpoint.com/v2/url?u=http-</a>3A__slepc.upv.es_documentation_current_docs_manualpages_ST_STPrecondSetMatForPC.html&d=DwIFaQ&c=imBPVzF25OnBgGmVOlcsiEgHoG1i6YHLR0Sj_gZ4adc&r=HM_IVChvU7XXVyLZB9xFCbw5_4LZNuzGEpVw5dUCyyo&m=L-kD1fTipDH_5BEC-T3FTfZKFQBUO6_ymvkpAiFfe38&s=kZTO6rmsUygc5qVZm5IkYMfjEMPVXloxoE0xi16ccUQ&e=<br class="">
in combination with pc=PCMAT.<br class="">
<br class="">
For a Krylov solver you will have to do something like this:<br class="">
<br class="">
     call EPSSetWhichEigenpairs(eps,EPS_TARGET_REAL,ierr);CHKERRA(ierr)<br class="">
     call EPSSetTarget(eps,sigma,ierr);CHKERRA(ierr)<br class="">
     call EPSGetST(eps,st,ierr);CHKERRA(ierr)<br class="">
     call STSetType(st,STSINVERT,ierr);CHKERRA(ierr)<br class="">
     call STGetKSP(st,ksp,ierr);CHKERRA(ierr)<br class="">
     call KSPGetPC(ksp,pc,ierr);CHKERRA(ierr)<br class="">
     call KSPSetType(ksp,KSPBCGS,ierr);CHKERRA(ierr)<br class="">
     call PCSetType(pc,PCNONE,ierr);CHKERRA(ierr)  ! dummy preconditioner<br class="">
<br class="">
     ! provide the preconditioner matrix after setup<br class="">
     call EPSSetUp(eps,ierr);CHKERRA(ierr)<br class="">
     call PCSetOperators(pc,PETSC_NULL_MAT,my_matmul_K,ierr);CHKERRA(ierr)<br class="">
     call PCSetType(pc,PCMAT,ierr);CHKERRA(ierr)<br class="">
<br class="">
We may change the interface in future versions so that this can be done in a more natural way.<br class="">
<br class="">
Jose<br class="">
<br class="">
<br class="">
<blockquote type="cite" class="">El 7 ago 2018, a las 16:42, Kenneth Hall, SC.D. <kenneth.c.hall@duke.edu> escribió:<br class="">
<br class="">
Hello,<br class="">
<br class="">
I am new to SLEPc and Petsc, so I would appreciate any help I can get with my question. I have a large eigenvalue problem I am solving, based on the linearization of a nonlinear code.  The eigenvalue problem takes the form<br class="">
<br class="">
A.x = lambda.B.x<br class="">
<br class="">
I can generate the matrix B explicitly, but I can only form products of A.x, and not A itself, so I use a shell matrix.   I also have an iterative solver I could potentially use as a preconditioner. That is, I have an operator K^-1 where<br class="">
<br class="">
K^-1 \approx  (A - sigma.B)^-1<br class="">
<br class="">
But again, I don’t have K^-1 explicitly, rather I can form products of K^-1 with a vector, so I would need to use a shell matrix.<br class="">
<br class="">
I have had success with the un-preconditioned case, but can’t figure out how to do the preconditioned case.  Any suggestions would be greatly appreciated.  I work with Fortran.  A code snippet with the un-preconditioned set-up follows. <br class="">
<br class="">
Thanks in advance.<br class="">
Kenneth Hall<br class="">
<br class="">
<br class="">
===============================================================<br class="">
<br class="">
!<br class="">
!.... initialize slepc<br class="">
     call SlepcInitialize(PETSC_NULL_CHARACTER,ierr) ; if (ierr .ne. 0) then ; stop 'SlepcInitialize failed' ; end if<br class="">
!<br class="">
!.... get size and rank of MPI process<br class="">
     call MPI_Comm_size(PETSC_COMM_WORLD,sz  ,ierr)                                               ; CHKERRA(ierr)<br class="">
     call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)                                               ; CHKERRA(ierr)<br class="">
!<br class="">
     if (sz .ne. 1) then ; SETERRA(PETSC_COMM_SELF,1,'This is a uniprocessor example only!'); endif<br class="">
!<br class="">
!.... get user options<br class="">
     call PetscOptionsGetInt (PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-n'   ,n ,flg,ierr)        ; CHKERRA(ierr)<br class="">
     call PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-terse',terse,ierr)        ; CHKERRA(ierr)<br class="">
!<br class="">
!.... write out option list<br class="">
     if (rank .eq. 0 ) then<br class="">
        write(*,*) 'n        = ',n<br class="">
     end if<br class="">
!<br class="">
!.... Register the matrix-vector subroutine for the operator that defines the eigensystem, Ax=kBx<br class="">
     call MatCreateShell(PETSC_COMM_WORLD,n,n,n,n,PETSC_NULL_INTEGER,A,ierr)                      ; CHKERRA(ierr)<br class="">
     call MatShellSetOperation(A,MATOP_MULT,my_matmul_A,ierr)                                     ; CHKERRA(ierr)<br class="">
!<br class="">
!.... set up B matrix<br class="">
     call MatCreate(PETSC_COMM_WORLD,B,ierr)                                                      ; CHKERRA(ierr)<br class="">
     call MatSetSizes(B,n,n,n,n,ierr)                                                             ; CHKERRA(ierr)<br class="">
     call MatSetFromOptions(B,ierr)                                                               ; CHKERRA(ierr)<br class="">
     call MatSetUp(B,ierr)                                                                        ; CHKERRA(ierr)<br class="">
     call My_Matrix_B(B,n)<br class="">
!<br class="">
!.... set up eigenvectors<br class="">
     call VecCreate(PETSC_COMM_WORLD,xr,ierr)                                                      ; CHKERRA(ierr)<br class="">
     call VecSetSizes(xr,n,n,ierr)                                                                 ; CHKERRA(ierr)<br class="">
     call VecSetFromOptions(xr,ierr)                                                               ; CHKERRA(ierr)<br class="">
     call VecSetUp(xr,ierr)                                                                        ; CHKERRA(ierr)<br class="">
!<br class="">
     call VecCreate(PETSC_COMM_WORLD,xi,ierr)                                                      ; CHKERRA(ierr)<br class="">
     call VecSetSizes(xi,n,n,ierr)                                                                 ; CHKERRA(ierr)<br class="">
     call VecSetFromOptions(xi,ierr)                                                               ; CHKERRA(ierr)<br class="">
     call VecSetUp(xi,ierr)                                                                        ; CHKERRA(ierr)<br class="">
!<br class="">
!.... Create eigensolver context<br class="">
     call EPSCreate(PETSC_COMM_WORLD,eps,ierr)                                                    ; CHKERRA(ierr)<br class="">
!<br class="">
!.... Set operators. In this case, it is a generalized eigenvalue problem<br class="">
     call EPSSetOperators(eps,A,B,ierr)                                                           ; CHKERRA(ierr)<br class="">
     call EPSSetProblemType(eps,EPS_GNHEP,ierr)                                                   ; CHKERRA(ierr)<br class="">
!<br class="">
!.... set default options<br class="">
     call EPSSetType(eps,EPSKRYLOVSCHUR,ierr)                                                     ; CHKERRA(ierr)<br class="">
!<br class="">
!.... set default eigenvalues to find<br class="">
     call EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL,ierr)                                        ; CHKERRA(ierr)<br class="">
!<br class="">
!.... number of requested eigenvalues<br class="">
     call EPSSetDimensions(eps,nev,PETSC_DEFAULT_INTEGER,PETSC_DEFAULT_INTEGER,ierr)              ; CHKERRA(ierr)<br class="">
!<br class="">
!.... set stopping criteria<br class="">
     call EPSSetTolerances(eps,tol,maxit,ierr)                                                    ; CHKERRA(ierr)<br class="">
!<br class="">
!.... set up monitors<br class="">
     call PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,vf,ierr)      ; CHKERRA(ierr)<br class="">
     call EPSMonitorSet(eps,EPSMONITORFIRST,vf,PetscViewerAndFormatDestroy,ierr)                  ; CHKERRA(ierr)<br class="">
     call SlepcConvMonitorCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,ctx,ierr)         ; CHKERRA(ierr)<br class="">
     call EPSMonitorSet(eps,EPSMONITORCONVERGED,ctx,SlepcConvMonitorDestroy,ierr)                 ; CHKERRA(ierr)<br class="">
!<br class="">
!.... Set solver parameters at runtime<br class="">
     call EPSSetFromOptions(eps,ierr)                                                             ; CHKERRA(ierr)<br class="">
!<br class="">
!.... Solve the eigensystem<br class="">
     call EPSSolve(eps,ierr)                                                                      ; CHKERRA(ierr)<br class="">
<br class="">
</blockquote>
<br class="">
</blockquote>
<br class="">
</div>
</div>
</body>
</html>