Skip to content

[WIP][solver] Implement block preconditioner for arrowhead like matrices#1872

Open
maxfirmbach wants to merge 1 commit into4C-multiphysics:mainfrom
maxfirmbach:Block-preconditioner-for-augmented-systems
Open

[WIP][solver] Implement block preconditioner for arrowhead like matrices#1872
maxfirmbach wants to merge 1 commit into4C-multiphysics:mainfrom
maxfirmbach:Block-preconditioner-for-augmented-systems

Conversation

@maxfirmbach
Copy link
Copy Markdown
Contributor

Description and Context

This PR aims to implement a block preconditioner for "arrowhead" like block matrices involving the calculation of a proper Schur complement. This is still WIP as there are a few open questions on the implementation part. This PR presents a first implementation with a respective working test case.

@maxfirmbach maxfirmbach force-pushed the Block-preconditioner-for-augmented-systems branch from 7ee12dd to 482db04 Compare March 18, 2026 17:15
Comment on lines +97 to +104
/**
* Generic block preconditioner for arrowhead matrices of the form
* | K^S -M^T |
* | K^B D^T |
* | -M D |
* based on a Gauss-Seidel procedure using a Schur complement.
*/
class ArrowHeadPreconditionerFactory : public Teko::GaussSeidelPreconditionerFactory
Copy link
Copy Markdown
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@maxfirmbach Are you focusing on systems with zero block on the main diagonal? If not, the BGS(AMG) FSI preconditioner by Gee et al. (2011) also relies on such an arrowhead layout and, thus, could serve as another test case here.

Copy link
Copy Markdown
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@mayrmt Currently yes, but not necessarily in the near future.

Copy link
Copy Markdown
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Sure! FSI can also due Schur complements with the fluid block on the main diagonal. It is not part of Gee et al. (2011), but I have written it in my PhD thesis. There is also a paper by Ulrich Langer from 2016, that shows Schur complement preconditioners for FSI.

Copy link
Copy Markdown
Member

@mayrmt mayrmt left a comment

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@maxfirmbach I assume that such code would normally go into Teko directly, if we were able to update Trilinos, right?

Comment on lines +78 to +82
if (A->matrix(A->rows() - 1, A->cols() - 1).norm_inf() < 1e-12)
{
maps.push_back(
Core::Utils::shared_ptr_from_ref(A->matrix(A->rows() - 1, A->cols() - 1).domain_map()));
}
Copy link
Copy Markdown
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I'd hope that we can still extract the map from the operator domain map of the blocked operator, since the last blocked column of the blocked operator should be set properly ... at least I hope...

std::string precStr = "None";
if (pl.isParameter(inverse_type)) invStr = pl.get<std::string>(inverse_type);
if (pl.isParameter(preconditioner_type)) precStr = pl.get<std::string>(preconditioner_type);
solveType_ = Teko::GS_UseLowerTriangle;
Copy link
Copy Markdown
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

This hard-codes BGS. Is that intended?


// inserting inverse factory into vector
std::string invStr2 = pl.get<std::string>(fieldName);
if (position > (int)inverses.size())
Copy link
Copy Markdown
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Suggested change
if (position > (int)inverses.size())
if (position > static_cast<int>(inverses.size()))


// inserting preconditioner factory into vector
std::string precStr2 = pl.get<std::string>(fieldName);
if (position > (int)preconditioners.size())
Copy link
Copy Markdown
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Suggested change
if (position > (int)preconditioners.size())
if (position > static_cast<int>(preconditioners.size()))

@maxfirmbach
Copy link
Copy Markdown
Contributor Author

maxfirmbach commented Mar 20, 2026

@maxfirmbach I assume that such code would normally go into Teko directly, if we were able to update Trilinos, right?

@mayrmt Not necessarily, I think the current implementation is really trimmed to work for a specific block matrix type / layout.

Edit: The implementation would look a bit different though I guess.

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment

Projects

None yet

Development

Successfully merging this pull request may close these issues.

2 participants