Skip to content

Commit

Permalink
Fixing section titles
Browse files Browse the repository at this point in the history
  • Loading branch information
Umberto Zerbinati committed Jun 4, 2024
1 parent 73538cd commit eb7506b
Show file tree
Hide file tree
Showing 2 changed files with 8 additions and 9 deletions.
4 changes: 2 additions & 2 deletions docs/src/PETScPC/poisson.py.rst
Original file line number Diff line number Diff line change
@@ -1,5 +1,5 @@
Preconditioning the Poisson problem
=====================================
Vertex Patch smoothing for p-Multigrid preconditioners
========================================================

In this tutorial, we explore using `PETSc PC` as a preconditioner inside NGSolve preconditioning infrastructure.
We will focus our attention on the Poisson problem, and we will consider different preconditioning strategies.
Expand Down
13 changes: 6 additions & 7 deletions docs/src/PETScPC/stokes.py.rst
Original file line number Diff line number Diff line change
@@ -1,5 +1,5 @@
Preconditioning the Stokes problem
===================================
Boffi-Lovadina Augmentation for Bernardi-Raugel discretization of the Stokes problem
======================================================================================

In this tutorial, we explore constructing preconditioners for saddle point problems using `PETSc PC`.
In particular, we will consider a Bernardi-Raugel inf-sup stable discretization of the Stokes problem, i.e.
Expand Down Expand Up @@ -33,8 +33,8 @@ Such a discretization can easily be constructed using NGSolve as follows: ::
else:
mesh = Mesh(ngm.Mesh.Receive(COMM_WORLD))
nu = Parameter(1.0)
V = VectorH1(mesh, order=4, dirichlet="wall|inlet|cyl", autoupdate=True)
Q = L2(mesh, order=2, autoupdate=True)
V = VectorH1(mesh, order=2, dirichlet="wall|inlet|cyl", autoupdate=True)
Q = L2(mesh, order=0, autoupdate=True)
u,v = V.TnT(); p,q = Q.TnT()
a = BilinearForm(nu*InnerProduct(Grad(u),Grad(v))*dx)
a.Assemble()
Expand Down Expand Up @@ -184,7 +184,7 @@ This formulation can easily be constructed by adding a new velocity block in the
gfu.Set(uin, definedon=mesh.Boundaries("inlet"))
sol = BlockVector( [gfu.vec, gfp.vec] )

print("-----------|Augmented LU|-----------")
print("-----------|Boffi--Lovadina Augmentation LU|-----------")
solvers.MinRes (mat=K, pre=C, rhs=rhs, sol=sol, tol=1e-10,
printrates=True, initialize=False)
Draw(gfu)
Expand Down Expand Up @@ -216,7 +216,7 @@ This is not ideal for large problems, and we can use a `Hypre` preconditioner fo
C = BlockMatrix( [ [twolvpre, None], [None, mGpre] ] )
gfu.vec.data[:] = 0; gfp.vec.data[:] = 0;
gfu.Set(uin, definedon=mesh.Boundaries("inlet"))
print("-----------|Augmented Two Level Additive Schwarz|-----------")
print("-----------|Boffi--Lovadina Augmentation Two Level Additive Schwarz|-----------")
solvers.MinRes (mat=K, pre=C, rhs=rhs, sol=sol, tol=1e-10,
printrates=True, initialize=False)
Draw(gfu)
Expand Down Expand Up @@ -345,4 +345,3 @@ We see that a purely algebraic approach based on the Sherman-Morrisson-Woodbory
- 100 (1.06e-03)
* - Augmentation Schermon-Morrisson-Woodbory
- 84 (1.16e-07)

0 comments on commit eb7506b

Please sign in to comment.