[Trilinos-Users] MueLu setting block size warning (Tpetra+Belos+MueLu)

Brian Granzow brian.granzow at gmail.com
Tue Apr 25 13:49:02 EDT 2017


Hello!

I have an application that uses Tpetra + MueLu + Belos to solve linear
systems obtained via finite element assembly. To create the MueLu
preconditioner, I am using the function:
MueLu::CreateTpetraPreconditioner:

https://github.com/trilinos/Trilinos/blob/master/packages/muelu/adapters/tpetra/MueLu_CreateTpetraPreconditioner.hpp#L49

I am passing this function the required Tpetra::Operator (inA) and
input Teuchos::ParameterList (inParamList), as well as the optional
coordinate Tpetra::MultiVector (inCoords). For a problem with 3
degrees of freedom per node (linear elasticity), I specify the MueLu
input parameter "number of equations" to be equal to 3, and I see the
warning below:

```
******* WARNING *******
Setting matrix block size to 3 (value of the parameter in the list)
instead of 1 (provided matrix).
You may want to check "number of equations" (or "PDE equations" for
factory style list) parameter.
```

I am curious if anyone has experienced this before, and if the fix is
immediately obvious.

---

For more detail, the Tpetra::Operator (inA) that is passed to
MueLu::CreateTpetraPreconditioner by dynamically casting an existing
Tpetra::CrsMatrix to a Tpetra::Operator. The coordinate
Tpetra::MultiVector is constructed from a Tpetra::Map that describes
the parallel local to global distribution of nodes (in this case mesh
vertices). The distinction here being that the number of rows of the
coordinate MultiVector is smaller than the number of rows of my
Tpetra::CrsMatrix by a factor of 3 (since there are 3 degrees of
freedom per node). In fact, I've used this exact same approach in a
small application in Albany where I did not see the above warning:

https://github.com/gahansen/Albany/blob/master/src/CTM/CTM_LinearSolver.cpp#L50

My current application is nearly identical to this Albany application,
so I'm unsure where I am going wrong.  I am pretty confident that the
Tpetra::Operator is formed correctly as I am able to solve my FEM
problem with correct results when I simply omit passing the coordinate
vector to MueLu::CreateTpetraPreconditioner and leave the "number of
equations" parameter as the default value.

Thanks for the consideration,
Brian Granzow


More information about the Trilinos-Users mailing list