Trilinos block matrix/vector + gmres

26 views
Skip to first unread message

Praveen C

unread,
Feb 28, 2026, 8:29:06 AM (3 days ago) Feb 28
to Deal.II Googlegroup
Dear all

We are exploring different linear solvers for a mixed H(div) method for poisson equation.

We want to try both block solvers based on Schur form and also gmres, direct, etc.

So we declare matrices and vectors of block form.

But then we are not able to use TrilinosWrappers::SolverGMRES with these blocked variables.

Is this a fundamental limitation or something else ?

I have attached a MWE which shows compiler error when I call GMRES.

thanks
praveen
main.cc

Wolfgang Bangerth

unread,
Feb 28, 2026, 11:52:29 AM (3 days ago) Feb 28
to dea...@googlegroups.com
Hi Praveen,
yes, when you use TrilinosWrappers::SolverGMRES, then it wants a single
Trilinos matrix which a block matrix is not. But you can use
dealii::SolverGMRES with any kind of linear operator, which a block matrix is.
So just use the latter.

Best
W.

Praveen C

unread,
Feb 28, 2026, 10:51:58 PM (2 days ago) Feb 28
to dea...@googlegroups.com
Thanks, I did not know we could use Trilinos objects with deal.II solver.

I wanted to try Trilinos because the deal.II solver did not have AMG preconditioner. Also I could not use ILU preconditioner with dealii blockmatrix/vector and dealii::solverGMRES.

petsc:gmres seems to support more general objects, so I will try that one. Petsc also has AMG.

thanks
praveen

--
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see https://groups.google.com/d/forum/dealii?hl=en
---
You received this message because you are subscribed to the Google Groups "deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email to dealii+un...@googlegroups.com.
To view this discussion visit https://groups.google.com/d/msgid/dealii/e753b8d7-df9b-4cd2-96b2-66c3651be488%40colostate.edu.

Praveen C

unread,
Mar 2, 2026, 3:00:25 AM (yesterday) Mar 2
to Praveen C, Deal. II Googlegroup
Hello all

Here I am trying PetscWrappers::BlockSparseMatrix/Vector and PETScWrappers::SolverGMRES and this wont work, because gmres solve function expects MatrixBase type but I have BlockMatrixBase type.

GMRES should still work on block matrices since it only needs vmult, so maybe this is just a matter of implementing a wrapper for block matrix case.

If this function takes BlockMatrixBase


I am hoping gmres should still work.

Could you somebody with more experience comment on this, before I do some experiments ?

thanks
praveen

James Shi

unread,
Mar 2, 2026, 3:07:32 AM (yesterday) Mar 2
to dea...@googlegroups.com
Hi Praveen,

I've been using PETSc's BlockVector and Trilinos' BlockVector with dealii::SolverGMRES<almost_any_vector_type> for some time, even the unstable Trilinos' Tpetra vectors works just fine.

You can try using dealii::SolverGMRES<> as Prof. Wolfgang Bangerth suggested.

Best,
James Shi

Praveen C <cpra...@gmail.com> 于2026年3月2日周一 16:00写道:

Praveen C

unread,
Mar 2, 2026, 3:26:28 AM (yesterday) Mar 2
to Deal. II Googlegroup



I've been using PETSc's BlockVector and Trilinos' BlockVector with dealii::SolverGMRES<almost_any_vector_type> for some time, even the unstable Trilinos' Tpetra vectors works just fine.

You can try using dealii::SolverGMRES<> as Prof. Wolfgang Bangerth suggested.


Thank you for sharing your experience.

Is it possible to use AMG preconditioner from Petsc or Trilinos with dealii::SolverGMRES ?

thanks
praveen

James Shi

unread,
Mar 2, 2026, 3:48:14 AM (24 hours ago) Mar 2
to dea...@googlegroups.com
Hi Praveen,

I haven’t tried plugging PETSc/Trilinos AMG “directly” into deal.II’s SolverGMRES.

For the non-block PETSc/Trilinos vector types, I prefer to use PETSc’s KSP or Trilinos’ Belos solvers directly, with their AMG preconditioners. For the block PETSc/Trilinos vector types, I implement a block preconditioner for dealii::SolverGMRES; inside that preconditioner, I still call KSP or Belos directly.

Regarding AMG: I don’t know whether BoomerAMG (via PETSc/HYPRE) works well with block matrices—I’ve only used it with regular matrices. Trilinos’ MueLu is also complicated; based on recent reports it seems to support block matrices, but it requires striding settings that are not very well documented. I’m also concerned that Trilinos’ definition of a “block matrix” may not match deal.II’s  BlockSparseMatrix.

Best,
James



Praveen C <pra...@tifrbng.res.in> 于2026年3月2日周一 16:26写道:
--
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see https://groups.google.com/d/forum/dealii?hl=en
---
You received this message because you are subscribed to the Google Groups "deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email to dealii+un...@googlegroups.com.

Wolfgang Bangerth

unread,
Mar 2, 2026, 4:47:42 AM (23 hours ago) Mar 2
to dea...@googlegroups.com
On 3/2/26 00:26, Praveen C wrote:
>
> Is it possible to use AMG preconditioner from Petsc or Trilinos with
> dealii::SolverGMRES ?

The dealii::Solver* classes all just require the matrix and the preconditioner
to be linear operators, i.e., classes with a vmult() function. In your case,
these may be block operators that internally call AMG functions from PETSc or
Trilinos. You may be interested in taking a look here, for example:

https://dealii.org/developer/doxygen/deal.II/step_22.html#step_22-BlockSchurcomplementpreconditioner

Best
W.

Praveen C

unread,
Mar 2, 2026, 10:22:51 PM (5 hours ago) Mar 2
to dea...@googlegroups.com
Thank you, we are also implementing block solvers/preconditioner as in these examples you suggested.

We also want to try GMRES on the entire block matrix with say ILU preconditioner, but does not seem to be possible, since all available ILU do  not take block matrices.

If we assemble into PETSc::BlockSparseMatrix, there is also no way to copy that into PETSc::SparseMatrix. So it looks like we have to switch to SparseMatrix during assembly itself for this purpose.

I also asked on Petsc list,


and they may make it possible to use block matrices with ILU, but it does not work at present.

best
praveen

--
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see https://groups.google.com/d/forum/dealii?hl=en
---
You received this message because you are subscribed to the Google Groups "deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email to dealii+un...@googlegroups.com.

Wolfgang Bangerth

unread,
12:26 AM (3 hours ago) 12:26 AM
to dea...@googlegroups.com
On 3/2/26 19:22, Praveen C wrote:
>
> We also want to try GMRES on the entire block matrix with say ILU
> preconditioner, but does not seem to be possible, since all available ILU do
> not take block matrices.
>
> If we assemble into PETSc::BlockSparseMatrix, there is also no way to copy
> that into PETSc::SparseMatrix. So it looks like we have to switch to
> SparseMatrix during assembly itself for this purpose.

I think you misunderstand what the purpose of block matrices is, then. It is
not *necessary* to use block matrices just because you have a system of
equations. If you want to use solvers and preconditioners that work on the
*entire* matrix, then by all means assemble everything into a single matrix.
The point of block matrices and vectors is so that you *can* decompose the
matrix into blocks because you want to use a block preconditioner. But if
that's not what you want to do, then just don't use block objects.

Best
W.
Reply all
Reply to author
Forward
0 new messages