This page was generated from unit-3.7-nonlinear/nonlinear.ipynb.

3.7 Nonlinear problems

We want to solve a nonlinear PDE.

A simple scalar PDE

We consider the simple PDE

\[- \Delta u + 3 u^3 = 1 \text{ in } \Omega\]

on the unit square \(\Omega = (0,1)^2\).

We note that this PDE can also be formulated as a nonlinear minimization problem (cf. 3.8).

[1]:
from netgen import gui
from netgen.geom2d import unit_square
from ngsolve import *

mesh = Mesh (unit_square.GenerateMesh(maxh=0.3))

In NGSolve we can solve the PDE conveniently using the linearization feature of SymbolicBFI.

The BilinearForm (which is not bilinear!) needed in the weak formulation is

\[A(u,v) = \int_{\Omega} \nabla u \nabla v + 3 u^3 v - 1 v ~ dx \quad ( = 0 ~ \forall~v \in H^1_0)\]
[2]:
V = H1(mesh, order=3, dirichlet=[1,2,3,4])
u,v = V.TnT()
a = BilinearForm(V)
a += (grad(u) * grad(v) + 3*u**3*v- 1 * v)*dx

Newton’s method

We use Newton’s method and make the loop:

  • Given an initial guess \(u^0\)

  • loop over \(i=0,..\) until convergence:

  • Compute linearization: \(A u^i + \delta A(u^i) \Delta u^{i} = 0\):

    • \(f^i = A u^i\)

    • \(B^i = \delta A(u^i)\)

    • Solve \(B^i \Delta u^i = -f^i\)

  • Update \(u^{i+1} = u^i + \Delta u^{i}\)

  • Evaluate stopping criteria

As a stopping criteria we take \(\langle A u^i,\Delta u^i \rangle = \langle A u^i, A u^i \rangle_{(B^i)^{-1}}< \varepsilon\).

[3]:
def SimpleNewtonSolve(gfu,a,tol=1e-13,maxits=25):
    res = gfu.vec.CreateVector()
    du = gfu.vec.CreateVector()
    fes = gfu.space
    for it in range(maxits):
        print ("Iteration {:3}  ".format(it),end="")
        a.Apply(gfu.vec, res)
        a.AssembleLinearization(gfu.vec)
        du.data = a.mat.Inverse(fes.FreeDofs()) * res
        gfu.vec.data -= du

        #stopping criteria
        stopcritval = sqrt(abs(InnerProduct(du,res)))
        print ("<A u",it,", A u",it,">_{-1}^0.5 = ", stopcritval)
        if stopcritval < tol:
            break
[4]:
gfu = GridFunction(V)
Draw(gfu,mesh,"u")
SimpleNewtonSolve(gfu,a)
Iteration   0  <A u 0 , A u 0 >_{-1}^0.5 =  0.1874380387122074
Iteration   1  <A u 1 , A u 1 >_{-1}^0.5 =  9.417946030048514e-05
Iteration   2  <A u 2 , A u 2 >_{-1}^0.5 =  8.541642343381769e-11
Iteration   3  <A u 3 , A u 3 >_{-1}^0.5 =  3.6594922142680475e-17

There are also some solvers shipped with NGSolve now:

[5]:
from ngsolve.solvers import *
help(Newton)
Help on function Newton in module ngsolve.nonlinearsolvers:

Newton(a, u, freedofs=None, maxit=100, maxerr=1e-11, inverse='umfpack', dirichletvalues=None, dampfactor=1, printing=True, callback=None)
    Newton's method for solving non-linear problems of the form A(u)=0.

    Parameters
    ----------
    a : BilinearForm
      The BilinearForm of the non-linear variational problem. It does not have to be assembled.

    u : GridFunction
      The GridFunction where the solution is saved. The values are used as initial guess for Newton's method.

    freedofs : BitArray
      The FreeDofs on which the assembled matrix is inverted. If argument is 'None' then the FreeDofs of the underlying FESpace is used.

    maxit : int
      Number of maximal iteration for Newton. If the maximal number is reached before the maximal error Newton might no converge and a warning is displayed.

    maxerr : float
      The maximal error which Newton should reach before it stops. The error is computed by the square root of the inner product of the residuum and the correction.

    inverse : string
      A string of the sparse direct solver which should be solved for inverting the assembled Newton matrix.

    dampfactor : float
      Set the damping factor for Newton's method. If dampfactor is 1 then no damping is done. If value is < 1 then the damping is done by the formula 'min(1,dampfactor*numit)' for the correction, where 'numit' denotes the Newton iteration.

    printing : bool
      Set if Newton's method should print informations about the actual iteration like the error.

    Returns
    -------
    (int, int)
      List of two integers. The first one is 0 if Newton's method did converge, -1 otherwise. The second one gives the number of Newton iterations needed.

[6]:
gfu.vec[:]=0
Newton(a,gfu,freedofs=gfu.space.FreeDofs(),maxit=100,maxerr=1e-11,inverse="umfpack",dampfactor=1,printing=True)
Newton iteration  0
err =  0.1874380387122074
Newton iteration  1
err =  9.417946030052061e-05
Newton iteration  2
err =  8.541642208644915e-11
Newton iteration  3
err =  4.2513011097456555e-17
[6]:
(0, 4)

A trivial problem:

\[5 u^2 = 1, \qquad u \in \mathbb{R}.\]
[7]:
V = NumberSpace(mesh)
u,v = V.TnT()
a = BilinearForm(V)
a += ( 5*u*u*v - 1 * v)*dx
gfu = GridFunction(V)
gfu.vec[:] = 1
SimpleNewtonSolve(gfu,a)

print("\nscalar solution", gfu.vec[0], "(exact: ", sqrt(0.2), ")")
Iteration   0  <A u 0 , A u 0 >_{-1}^0.5 =  1.2649110640673518
Iteration   1  <A u 1 , A u 1 >_{-1}^0.5 =  0.3265986323710904
Iteration   2  <A u 2 , A u 2 >_{-1}^0.5 =  0.04114755998989123
Iteration   3  <A u 3 , A u 3 >_{-1}^0.5 =  0.0008574269268691778
Iteration   4  <A u 4 , A u 4 >_{-1}^0.5 =  3.8832745226099997e-07
Iteration   5  <A u 5 , A u 5 >_{-1}^0.5 =  7.979879233426059e-14

scalar solution 0.4472135954999579 (exact:  0.4472135954999579 )

Another example: Stationary Navier-Stokes:

Find \(\mathbf{u} \in \mathbf{V}\), \(p \in Q\), \(\lambda \in \mathbb{R}\) so that

\begin{align} \int_{\Omega} \nu \nabla \mathbf{u} : \nabla \mathbf{v} + (\mathbf{u} \cdot \nabla) \mathbf{u} \cdot \mathbf{v}& - \int_{\Omega} \operatorname{div}(\mathbf{v}) p & &= \int \mathbf{f} \cdot \mathbf{v} && \forall \mathbf{v} \in \mathbf{V}, \\ - \int_{\Omega} \operatorname{div}(\mathbf{u}) q & & + \int_{\Omega} \lambda q &= 0 && \forall q \in Q, \\ & \int_{\Omega} \mu p & &= 0 && \forall \mu \in \mathbb{R}. \end{align}
[8]:
mesh = Mesh (unit_square.GenerateMesh(maxh=0.05)); nu = Parameter(1)
V = VectorH1(mesh,order=3,dirichlet="bottom|right|top|left")
Q = H1(mesh,order=2);
N = NumberSpace(mesh);
X = FESpace([V,Q,N])
(u,p,lam), (v,q,mu) = X.TnT()
a = BilinearForm(X)
a += (nu*InnerProduct(grad(u),grad(v))+InnerProduct(grad(u)*u,v)
      -div(u)*q-div(v)*p-lam*q-mu*p)*dx
[9]:
gfu = GridFunction(X)
gfu.components[0].Set(CoefficientFunction((4*x*(1-x),0)),
                      definedon=mesh.Boundaries("top"))
[10]:
SimpleNewtonSolve(gfu,a)
Draw(gfu.components[1],mesh,"p")
Draw(gfu.components[0],mesh,"u")
Iteration   0  <A u 0 , A u 0 >_{-1}^0.5 =  2.783769725648835
Iteration   1  <A u 1 , A u 1 >_{-1}^0.5 =  0.007914999989203112
Iteration   2  <A u 2 , A u 2 >_{-1}^0.5 =  6.483100496209905e-08
Iteration   3  <A u 3 , A u 3 >_{-1}^0.5 =  8.902307290345756e-16
[11]:
nu.Set(0.01)
SimpleNewtonSolve(gfu,a)
Redraw()
Iteration   0  <A u 0 , A u 0 >_{-1}^0.5 =  0.08811930239376543
Iteration   1  <A u 1 , A u 1 >_{-1}^0.5 =  0.008271034539448304
Iteration   2  <A u 2 , A u 2 >_{-1}^0.5 =  0.00010638964361914529
Iteration   3  <A u 3 , A u 3 >_{-1}^0.5 =  1.9401697072005596e-08
Iteration   4  <A u 4 , A u 4 >_{-1}^0.5 =  5.631845363639007e-16
[12]:
nu.Set(0.001)
SimpleNewtonSolve(gfu,a)
Redraw()
Iteration   0  <A u 0 , A u 0 >_{-1}^0.5 =  0.07337554371730436
Iteration   1  <A u 1 , A u 1 >_{-1}^0.5 =  0.05991827549295081
Iteration   2  <A u 2 , A u 2 >_{-1}^0.5 =  0.03627869069808476
Iteration   3  <A u 3 , A u 3 >_{-1}^0.5 =  0.03185175712415947
Iteration   4  <A u 4 , A u 4 >_{-1}^0.5 =  0.030219933881877784
Iteration   5  <A u 5 , A u 5 >_{-1}^0.5 =  0.024153091283687263
Iteration   6  <A u 6 , A u 6 >_{-1}^0.5 =  0.05128559692462751
Iteration   7  <A u 7 , A u 7 >_{-1}^0.5 =  0.09875512946340462
Iteration   8  <A u 8 , A u 8 >_{-1}^0.5 =  0.10389113309636115
Iteration   9  <A u 9 , A u 9 >_{-1}^0.5 =  0.23380427560538825
Iteration  10  <A u 10 , A u 10 >_{-1}^0.5 =  0.28366440056541214
Iteration  11  <A u 11 , A u 11 >_{-1}^0.5 =  0.5566537810275759
Iteration  12  <A u 12 , A u 12 >_{-1}^0.5 =  0.7554637043469036
Iteration  13  <A u 13 , A u 13 >_{-1}^0.5 =  1.82911175770442
Iteration  14  <A u 14 , A u 14 >_{-1}^0.5 =  5.793437732769543
Iteration  15  <A u 15 , A u 15 >_{-1}^0.5 =  8.114588244578368
Iteration  16  <A u 16 , A u 16 >_{-1}^0.5 =  150.36658957337215
Iteration  17  <A u 17 , A u 17 >_{-1}^0.5 =  179.1165288844912
Iteration  18  <A u 18 , A u 18 >_{-1}^0.5 =  333.1018042447549
Iteration  19  <A u 19 , A u 19 >_{-1}^0.5 =  694.8536663937053
Iteration  20  <A u 20 , A u 20 >_{-1}^0.5 =  300.5623981210263
Iteration  21  <A u 21 , A u 21 >_{-1}^0.5 =  490.8427454975774
Iteration  22  <A u 22 , A u 22 >_{-1}^0.5 =  6055.570753342729
Iteration  23  <A u 23 , A u 23 >_{-1}^0.5 =  10720.213827432566
Iteration  24  <A u 24 , A u 24 >_{-1}^0.5 =  2688.688431138383
[13]:
nu.Set(0.001)
gfu.components[0].Set(CoefficientFunction((4*x*(1-x),0)),definedon=mesh.Boundaries("top"))
Newton(a,gfu,maxit=20,dampfactor=0.1)
Redraw()
Newton iteration  0
err =  10.30258116032596
Newton iteration  1
err =  9.272321489775477
Newton iteration  2
err =  7.417880759733821
Newton iteration  3
err =  5.19258707937858
Newton iteration  4
err =  3.115713918052335
Newton iteration  5
err =  1.5582105779835163
Newton iteration  6
err =  0.6247189584105143
Newton iteration  7
err =  0.18752548084660411
Newton iteration  8
err =  0.037813002580542925
Newton iteration  9
err =  0.0038348847243832996
Newton iteration  10
err =  1.8494024509496642e-05
Newton iteration  11
err =  1.2678518698189276e-08
Newton iteration  12
err =  4.922514640262317e-15