If it is sequential, which it probably should be, then you can you MatLUFactorSymbolic(), MatLUFactorNumeric() and MatMatSolve() where you put a bunch of your right hand side vectors into a dense array; not all million of them but maybe 10 to 100 at a time.
Barry > On Aug 10, 2016, at 10:18 PM, Harshad Ranadive <harshadranad...@gmail.com> > wrote: > > Hi Barry > > The matrix A is mostly tridiagonal > > 1 α 0 ......... 0 > > α 1 α 0 .......0 > > > 0 α 1 α 0 ....0 > > > .................... > 0..............α 1 > > In some cases (periodic boundaries) there would be an 'α' in right-top-corner > and left-bottom corner. > > I am not using multigrid approach. I just implemented an implicit filtering > approach (instead of an explicit existing one) which requires the solution of > the above system. > > Thanks > Harshad > > On Thu, Aug 11, 2016 at 1:07 PM, Barry Smith <bsm...@mcs.anl.gov> wrote: > > Effectively utilizing multiple right hand sides with the same system can > result in roughly 2 or at absolute most 3 times improvement in solve time. > A great improvement but when you have a million right hand sides not a giant > improvement. > > The first step is to get the best (most efficient) preconditioner for you > problem. Since you have many right hand sides it obviously pays to spend more > time building the preconditioner so that each solve is faster. If you provide > more information on your linear system we might have suggestions. CFD so is > your linear system a Poisson problem? Are you using geometric or algebraic > multigrid with PETSc? It not a Poisson problem how can you describe the > linear system? > > Barry > > > > > On Aug 10, 2016, at 9:54 PM, Harshad Ranadive <harshadranad...@gmail.com> > > wrote: > > > > Hi All, > > > > I have currently added the PETSc library with our CFD solver. > > > > In this I need to use KSPSolve(...) multiple time for the same matrix A. I > > have read that PETSc does not support passing multiple RHS vectors in the > > form of a matrix and the only solution to this is calling KSPSolve multiple > > times as in example given here: > > http://www.mcs.anl.gov/petsc/petsc-current/src/ksp/ksp/examples/tutorials/ex16.c.html > > > > I have followed this technique, but I find that the performance of the code > > is very slow now. I basically have a mesh size of 8-10 Million and I need > > to solve the matrix A very large number of times. I have checked that the > > statement KSPSolve(..) is taking close to 90% of my computation time. > > > > I am setting up the matrix A, KSPCreate, KSPSetup etc just once at the > > start. Only the following statements are executed in a repeated loop > > > > Loop begin: (say million times !!) > > > > loop over vector length > > VecSetValues( ....) > > end > > > > VecAssemblyBegin( ... ) > > VecAssemblyEnd (...) > > > > KSPSolve (...) > > > > VecGetValues > > > > Loop end. > > > > Is there an efficient way of doing this rather than using KSPSolve multiple > > times? > > > > Please note my matrix A never changes during the time steps or across the > > mesh ... So essentially if I can get the inverse once would it be good > > enough? It has been recommended in the FAQ that matrix inverse should be > > avoided but would it be okay to use in my case? > > > > Also could someone please provide an example of how to use MatLUFactor and > > MatCholeskyFactor() to find the matrix inverse... the arguments below were > > not clear to me. > > IS row > > IS col > > const MatFactorInfo *info > > > > Apologies for a long email and thanks to anyone for help. > > > > Regards > > Harshad > > > > > > > > > > > >