Processing math: 100%

Nonlinear problems

We want to solve a nonlinear PDE.

A simple scalar PDE

We consider the simple PDE Δu+5u2=1 in Ω on the unit square Ω=(0,1)2. We note that this PDE can also be formulated as a nonlinear minimization problem (cf. unit 3.6 ).

In [1]:
import netgen.gui
%gui tk
import tkinter
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)=Ωuv+5u2v1v dx(=0  vH10)

In [2]:
V = H1(mesh, order=3, dirichlet=[1,2,3,4])
u = V.TrialFunction()
v = V.TestFunction()

a = BilinearForm(V)
a += SymbolicBFI( grad(u) * grad(v) + 5*u*u*v- 1 * v)

We use Newton's method and make the loop:

  • Given an initial guess u0
  • loop over i=0,.. until convergence:
    • Compute linearization: A(ui)+δA(ui)Δui=0:
      • fi=A(ui)
      • Bi=δA(ui)
      • Solve BiΔui=fi
    • Update ui+1=ui+Δui
    • Evaluate stopping criteria
  • Evaluate E(ui+1)

As a stopping criteria we take A(ui),Δui=A(ui),A(ui)(Bi)1<ε.

In [3]:
u = GridFunction(V)
Draw(u,mesh,"u")

res = u.vec.CreateVector()
du = u.vec.CreateVector()

for it in range(25):
    print ("Iteration",it)
    a.Apply(u.vec, res)
    a.AssembleLinearization(u.vec)
    #a.Assemble()
    
    du.data = a.mat.Inverse(V.FreeDofs()) * res
    u.vec.data -= du
    
    #stopping criteria
    stopcritval = sqrt(abs(InnerProduct(du,res)))
    print ("<A u",it,", A u",it,">_{-1}^0.5 = ", 
           stopcritval)
    if stopcritval < 1e-13:
        break
    
Iteration 0
<A u 0 , A u 0 >_{-1}^0.5 =  0.18745205801776055
Iteration 1
<A u 1 , A u 1 >_{-1}^0.5 =  0.0025712160474074893
Iteration 2
<A u 2 , A u 2 >_{-1}^0.5 =  5.324479829649632e-07
Iteration 3
<A u 3 , A u 3 >_{-1}^0.5 =  2.3040755730987054e-14