[PETSC #15886] Re: matrix-free additive schwarz?

Barry Smith petsc-maint at mcs.anl.gov
Sun Mar 4 21:52:27 CST 2007


  Aaron,

   You can use PCASMSetLocalSubdomains(pc,n,is[]) where is[] is an
array of one or more IS index sets that define the subdomains (with
overlap) for that process. Usually n is 1, which means one subdomain
per processor. The routine PCASMCreateSubdomains2D() may be helpful
for you in figuring out how to compute the is. You must also call
PCASMSetOverlap(pc,0); so that PETSc will not try to increase the
overlap from what you already provided. From this information the ASM
PC figures out the size of the local vectors and the VecScatters
that are needed to map to and from the local vectors. It will also
create the appropriate local KSP solvers for the subdomains.


  BUT, BUT, BUT what operator (Mat) are you going to use on the 
overlapped region and if you are going to use matrix-free on the 
subdomain from what will you construct a preconditioner?

  The way PCASM currently works is it calls MatGetSubMatrices() on the
original global matrix to get the Mat that is used as the operator 
on the subdomain. So what you need to do is for your matrix-free matrix
(which I assume you create with MatCreateShell()) do
MatShellSetOperation(mat,MAT_GET_SUBMATRICES,(void(*)(void))YourSubMat);
Where YourSubMat() is a function that returns your subdomain matrix;
which, of course, could be matrix-free or you could assemble it
for the subdomain. Of course, if you use matrix-free here then what
will you use for a preconditioner for the subdomain?? You could
always build some other matrix to use to construct the preconditioner.

  Let us know if this helps,

   Barry

There is on caveat, the is's you pass in get sorted by the ASM PC
so when you provide your local matrix it must be provided in the sorted
order. We could possibly fix this limitation if you are working fully
with matrix-free.




On Sun, 4 Mar 2007, P. Aaron Lott wrote:

> 
> Hi Matt,
> 
> I hope I can explain a bit better what it is I'm looking for, and also try to
> understand what I would need to provide to petsc in order to use its ASM
> preconditioner. It seems that in essence, for a matrix-free ASM
> implementation, I would need to provide either the restriction operators for
> each subdomain, or the assembled system matrix so that petsc's routines could
> determine the restriction operators.
> 
> Currently, I have a  Newton-Krylov spectral element code for solving the 2D
> incompressible Navier-Stokes equations. I perform all my matrix-vector
> products in a matrix-free framework. I now need a preconditioner for the
> subsidiary Advection-Diffusion equation solve, and I'm wanting to use ASM. I'm
> hoping to apply this preconditioner in a matrix free setting as well. My code
> applies matrix-vector products element-wise, and then performs an assembly
> (gather-scatter) procedure to get the global solution.
> 
> If all the subdomains for ASM are extensions of the elements from the spectral
> element discretization, then challenge is to construct the the restriction
> operators for each element to produce the ASM subdomains. From your message,
> it seems like petsc isn't able to construct these restriction operators unless
> I provide an assembled global system matrix.  If this is what petsc needs, I
> could compute the assembled system matrix (although storing this requires a
> lot of memory), and pass it to petsc. I thought though, that it may be
> possible to compute the restriction operators without accessing the assembled
> system, but rather through local-global maps providing indirect addressing. If
> someone has done this before, in a completely matrix-free setting, the same
> code could be used, or modified to obtain the restriction operators for each
> subdomain. I would be particularly interested in the second option if it is
> available either in petsc, or some package that links with petsc. I suppose
> the third option is, as you mentioned, providing my own preconditioner with
> PCSHELL, but my original reason for wanting to use petsc, was for it to
> provide me with the preconditioner...
> 
> Thanks,
> 
> -Aaron
> 
> 
> On Mar 4, 2007, at 3:41 PM, Matthew Knepley wrote:
> 
> > On 3/4/07, P. Aaron Lott <palott at ipst.umd.edu> wrote:
> > > Hi,
> > > 
> > > I'm new to petsc, and I've been reading through the petsc manual
> > > today trying to figure out if petsc provides a matrix-free
> > > implementation of a Newton-Krylov-Schwarz solver. It looks like in
> > > terms of the Newton-Krylov part, everything can be done in a matrix-
> > > free setting, but I can't tell if there are routines available to
> > > allow for a matrix-free Additive Schwarz preconditioner, and if so,
> > > what needs to be provided by the user to implement it. I'm very keen
> > > on using petsc to provide an additive schwarz preconditioner for my
> > > spectral element CFD code, but I need a matrix-free implementation
> > > because I don't have access to the global system matrix. Does anyone
> > > have experience with applying matrix-free preconditioners in petsc?
> > 
> > Its not quite clear what you want. You can, of course, apply any
> > preconditioner
> > you want with PCSHELL. We do have an ASM preconditioner, however the extra
> > support this provides is to
> > 
> > 1) automatically figure out the overlapping domains using the matrix
> >     nonzero structure
> > 
> > 2) solves the equation restricted to these overlapping partitions
> > 
> > The first job cannot be done without a structure of some sort, which we do
> > not
> > have in the matrix-free case, and the second depends on the first.
> > 
> >  Matt
> > 
> > > Thanks,
> > > 
> > > -Aaron
> > > 
> > > 
> > > 
> > > P. Aaron Lott
> > > Ph.D. Candidate
> > > 4239 Computer and Space Sciences Building
> > > University of Maryland
> > > College Park, MD 20742-4015
> > > 
> > > palott at ipst.umd.edu
> > > http://www.lcv.umd.edu/~palott
> > > Office: 301.405.4894
> > > Fax:      301.314.0827
> > > 
> > > 
> > > 
> > 
> > 
> > -- 
> > One trouble is that despite this system, anyone who reads journals widely
> > and critically is forced to realize that there are scarcely any bars to
> > eventual
> > publication. There seems to be no study too fragmented, no hypothesis too
> > trivial, no literature citation too biased or too egotistical, no design too
> > warped, no methodology too bungled, no presentation of results too
> > inaccurate, too obscure, and too contradictory, no analysis too
> > self-serving,
> > no argument too circular, no conclusions too trifling or too unjustified,
> > and
> > no grammar and syntax too offensive for a paper to end up in print. --
> > Drummond Rennie
> 
> P. Aaron Lott
> Ph.D. Candidate
> 4239 Computer and Space Sciences Building
> University of Maryland
> College Park, MD 20742-4015
> 
> palott at ipst.umd.edu
> http://www.lcv.umd.edu/~palott
> Office: 301.405.4894
> Fax:      301.314.0827
> 
> 




More information about the petsc-users mailing list