Hello<br><br>I&#39;ve almost got the petsc preconditioner free memory function to work. It now run the actual free-ing part, but segfaults on when it uses GMRES to solve the linear system.  The relevant section of code that setups the ksp/pc/localpc objects is below (in fortran) as well as the c-function for freeing the PC memory. <br>
<br>!***************** Fortran Code *********************<br>ksp_solver_type = &#39;asm&#39;<br>ksp_subspace = 50<br>asm_overlap =1<br>local_pc_ilu_level = 1<br>local_pc_ordering = &#39;rcm&#39;<br><br>!Setup global_ksp type and get the global PC<br>
  call KSPSetType(ksp,ksp_solver_type,ierr)<br>  call KSPGMRESSetRestart(ksp, ksp_subspace,ierr)<br>  call KSPSetPreconditionerSide(ksp,PC_RIGHT,ierr);<br>   call KSPGetPC(ksp,global_pc,ierr);             <br><br> ! Setup global_pc<br>
  call PCSetType(global_pc,&quot;asm&quot;,ierr)<br>  call PCASMSetOverlap(global_pc,asm_overlap,ierr)<br>  call PCSetup(global_pc,ierr)<br>  call PCASMGetSubKSP(global_pc, nlocal,  first, local_ksp, ierr )<br><br>! Setup local_pc<br>
  call KSPGetPC(local_ksp, local_pc, ierr )<br>  call PCSetType(local_pc, &#39;ilu&#39;, ierr)<br>  call PCFactorSetLevels(local_pc, local_pc_ilu_level, ierr)<br>  call PCFactorSetMatOrderingtype(local_pc, local_pc_ordering, ierr )<br>
  call KSPSetType(local_ksp, KSPPREONLY, ierr)<br><br>  call KSPSetup(global_ksp,ierr)<br>  call KSPSetUpOnBlocks(global_ksp,ierr)<br><br>  ! Free Memory<br>  call PCASMFreeSpace(global_pc)<br> <br>! ************* C Code **************************<br>
<br>void pcasmfreespace_(PC *inpc)<br>    {<br>      PC pc = *inpc;<br>      PC_ASM         *osm = (PC_ASM*)pc-&gt;data;<br>      PetscErrorCode ierr;<br>     <br>      if (osm-&gt;pmat) {<br>        if (osm-&gt;n_local_true &gt; 0) {<br>
          ierr = MatDestroyMatrices(osm-&gt;n_local_true,&amp;osm-&gt;pmat);<br>          CHKERRQ(ierr);<br>        }<br>      osm-&gt;pmat = 0;<br>      }<br>      <br>      if (pc-&gt;pmat) {<br>       ierr = MatDestroy(pc-&gt;pmat);<br>
       CHKERRQ(ierr);<br>       pc-&gt;pmat = 0;<br>      }<br>   }<br>! *********************************************8<br><br>When I run it, I get an error is PCApplyBAorAB<br>[0]PETSC ERROR: [0] MatMult line 1877 src/mat/interface/matrix.c<br>
[0]PETSC ERROR: [0] PCApplyBAorAB line 540 src/ksp/pc/interface/precon.c<br>[0]PETSC ERROR: [0] GMREScycle line 132 src/ksp/ksp/impls/gmres/gmres.c<br>[0]PETSC ERROR: [0] KSPSolve_GMRES line 227 src/ksp/ksp/impls/gmres/gmres.c<br>
[0]PETSC ERROR: [0] KSPSolve line 308 src/ksp/ksp/interface/itfunc.c<br><br>It appears then that you must keep the un-factored preconditioner matrix in memory. Is this actually the case? Or have I done something silly here?<br>
<br>Thanks,<br>Gaetan<br><br><br>
<br>
On Aug 25, 2011, at 1:29 PM, Gaetan Kenway wrote:<br>
<br>
&gt; Hello<br>
&gt;<br>
&gt; I&#39;ve managed to get the c-function for freeing preconditioner 
memory written. The contents of my new &#39;pcasmfreespace.c&#39; is below:<br>
&gt;<br>
&gt; #include &lt;private/pcimpl.h&gt;     /*I &quot;petscpc.h&quot; I*/<br>
&gt;<br>
&gt; typedef struct {<br>
&gt;  PetscInt   n, n_local, n_local_true;<br>
&gt;  PetscInt   overlap;             /* overlap requested by user */<br>
&gt;  KSP        *ksp;                /* linear solvers for each block */<br>
&gt;  VecScatter *restriction;        /* mapping from global to subregion */<br>
&gt;  VecScatter *localization;       /* mapping from overlapping to non-overlapping subregion */<br>
&gt;  VecScatter *prolongation;       /* mapping from subregion to global */<br>
&gt;  Vec        *x,*y,*y_local;      /* work vectors */<br>
&gt;  IS         *is;                 /* index set that defines each overlapping subdomain */<br>
&gt;  IS         *is_local;           /* index set that defines each non-overlapping subdomain, may be NULL */<br>
&gt;  Mat        *mat,*pmat;          /* mat is not currently used */<br>
&gt;  PCASMType  type;                /* use reduced interpolation, restriction or both */<br>
&gt;  PetscInt  type_set;            /* if user set this value (so won&#39;t change it for symmetric problems) */<br>
&gt;  PetscInt  same_local_solves;   /* flag indicating whether all local solvers are same */<br>
&gt;  PetscInt  sort_indices;        /* flag to sort subdomain indices */<br>
&gt; } PC_ASM;<br>
&gt;<br>
&gt; void pcasmfreespace_(PC pc)<br>
&gt;     {<br>
&gt;<br>
&gt;       PC_ASM         *osm = (PC_ASM*)pc-&gt;data;<br>
&gt;       PetscErrorCode ierr;<br>
&gt;<br>
&gt;       if (osm-&gt;pmat) {<br>
&gt;     if (osm-&gt;n_local_true &gt; 0) {<br>
&gt;       ierr = MatDestroyMatrices(osm-&gt;n_<div id=":7w">local_true,&amp;osm-&gt;pmat);CHKERRQ(ierr);<br>
&gt;     }<br>
&gt;     osm-&gt;pmat = 0;<br>
&gt;       }<br>
&gt;<br>
&gt;       if (pc-&gt;pmat) {ierr = MatDestroy(pc-&gt;pmat);CHKERRQ(ierr); pc-&gt;pmat = 0;}<br>
&gt;<br>
&gt;     }<br>
&gt;<br>
&gt; Note the underscore as I&#39;m trying to call it from Fortran. When I call it from fortran, I use:<br>
&gt;<br>
&gt; call pcasmfreespace(global_pc)<br>
&gt;<br>
&gt; This calls, the function, ok, but (according to Valgrind) I have an invalid read on the line containing:<br>
&gt;<br>
&gt; if (osm-&gt;pmat){<br>
&gt;<br>
&gt; I suspect this is something to do with passing the fortran 
&quot;pointer&quot; of the PC to c, or something along this lines. Is there 
anything else special you have to do to pass the &quot;fortran&quot; petsc objects
 to c?<br>
&gt;<br>
&gt; Thanks<br>
&gt;<br>
&gt; Gaetan<br>
&gt;<br>
&gt;<br>
&gt; Message: 2<br>
&gt; Date: Sun, 21 Aug 2011 17:22:28 -0500<br>
&gt; From: Barry Smith &lt;<a href="mailto:bsmith@mcs.anl.gov">bsmith@mcs.anl.gov</a>&gt;<br>
&gt; Subject: Re: [petsc-users] Freeing Preconditioner Matrix Space<br>
&gt; To: PETSc users list &lt;<a href="mailto:petsc-users@mcs.anl.gov">petsc-users@mcs.anl.gov</a>&gt;<br>
&gt; Message-ID: &lt;<a href="mailto:953121EF-B6AC-42EE-87BE-D4402C121652@mcs.anl.gov">953121EF-B6AC-42EE-87BE-D4402C121652@mcs.anl.gov</a>&gt;<br>
&gt; Content-Type: text/plain; charset=us-ascii<br>
&gt;<br>
&gt;<br>
&gt;  You don&#39;t need to put that in the PETSc source. Just built it in 
the same directory you build your application and link it in like any of
 your application code. You will need to stick<br>
&gt; #include &lt;private/pcimpl.h&gt;     /*I &quot;petscpc.h&quot; I*/<br>
&gt;<br>
&gt; typedef struct {<br>
&gt;  PetscInt   n, n_local, n_local_true;<br>
&gt;  PetscInt   overlap;             /* overlap requested by user */<br>
&gt;  KSP        *ksp;                /* linear solvers for each block */<br>
&gt;  VecScatter *restriction;        /* mapping from global to subregion */<br>
&gt;  VecScatter *localization;       /* mapping from overlapping to non-overlapping subregion */<br>
&gt;  VecScatter *prolongation;       /* mapping from subregion to global */<br>
&gt;  Vec        *x,*y,*y_local;      /* work vectors */<br>
&gt;  IS         *is;                 /* index set that defines each overlapping subdomain */<br>
&gt;  IS         *is_local;           /* index set that defines each non-overlapping subdomain, may be NULL */<br>
&gt;  Mat        *mat,*pmat;          /* mat is not currently used */<br>
&gt;  PCASMType  type;                /* use reduced interpolation, restriction or both */<br>
&gt;  PetscBool  type_set;            /* if user set this value (so won&#39;t change it for symmetric problems) */<br>
&gt;  PetscBool  same_local_solves;   /* flag indicating whether all local solvers are same */<br>
&gt;  PetscBool  sort_indices;        /* flag to sort subdomain indices */<br>
&gt; } PC_ASM;<br>
&gt;<br>
&gt; before the subroutine.<br>
&gt;<br>
&gt;  Barry<br>
&gt;<br>
&gt; On Aug 21, 2011, at 3:24 PM, Gaetan Kenway wrote:<br>
&gt;<br>
&gt; &gt; Hello<br>
&gt; &gt;<br>
&gt; &gt; I am attempting to implement a &quot;hack&quot; that was posted on the 
list a while back. I&#39;m working with the adjoint linear system solver for
 a CFD solver. I&#39;m using the ASM (or Block Jacobi) preconditioner with 
ILU(p) on each of the sub-domains. I use a different Preconditioner 
matrix (Pmat) than the actual jacobian. What I want to do is destroy the
 Pmat memory after the ILU factorization is performed. The hack that was
 posted is copied below:<br>
&gt; &gt;  PCASMFreeSpace(PC pc)<br>
&gt; &gt;    {<br>
&gt; &gt;   PC_ASM         *osm = (PC_ASM*)pc-&gt;data;<br>
&gt; &gt;   PetscErrorCode ierr;<br>
&gt; &gt;<br>
&gt; &gt;      if (osm-&gt;pmat) {<br>
&gt; &gt;     if (osm-&gt;n_local_true &gt; 0) {<br>
&gt; &gt;<br>
&gt; &gt;       ierr = MatDestroyMatrices(osm-&gt;n_<br>
&gt; local_true,&amp;osm-&gt;pmat);CHKERRQ(ierr);<br>
&gt; &gt;     }<br>
&gt; &gt;   osm-&gt;pmat = 0;<br>
&gt; &gt;   }<br>
&gt; &gt;   if (pc-&gt;pmat) {ierr = MatDestroy(pc-&gt;pmat);CHKERRQ(ierr); pc-&gt;pmat = 0;}<br>
&gt; &gt;<br>
&gt; &gt;   return 0;<br>
&gt; &gt;  }<br>
&gt; &gt;<br>
&gt; &gt; However, I&#39;ve had no luck actually getting the function 
compiled into petsc. There are no erorrs reported with i type &quot;make&quot; in 
the asm directory, but when I try to use the function in my application 
it can&#39;t find the symbol while linking. Where does it go in the asm.c 
file? Does it use &quot;static PetscErrorCode&quot; or &quot;PetscErrorCode 
PETSCKSP_DLLEXPORT&quot;? Does it have to be added to the .h include files? 
What has to be done for it work with Fortran?<br>
&gt; &gt;<br>
&gt; &gt; Any suggestions would be greatly appreciated.  This represents
 a significant chunk of my application&#39;s memory (10%-30%) and as such 
its too much to ignore.   Also is there any chance something like this 
would make it into an actual PETSc release?<br>
&gt; &gt;<br>
&gt; &gt; Gaetan<br>
&gt; &gt;<br>
&gt;<br>
&gt;<br>
</div><br>