1.1 First NGSolve example#
Let us solve the Poisson problem of finding \(u\) satisfying
Quick steps to solution:#
1. Import NGSolve and Netgen Python modules:#
[ ]:
from ngsolve import *
from ngsolve.webgui import Draw
2. Generate an unstructured mesh#
[ ]:
mesh = Mesh(unit_square.GenerateMesh(maxh=0.2))
mesh.nv, mesh.ne # number of vertices & elements
Here we prescribed a maximal mesh-size of 0.2 using the maxh flag.
[ ]:
Draw(mesh);
3. Declare a finite element space:#
[ ]:
fes = H1(mesh, order=2, dirichlet="bottom|right")
fes.ndof # number of unknowns in this space
Python’s help system displays further documentation.
[ ]:
help(fes)
4. Declare test function, trial function, and grid function#
Test and trial function are symbolic objects - called
ProxyFunctions- that help you construct bilinear forms (and have no space to hold solutions).GridFunctions, on the other hand, represent functions in the finite element space and contains memory to hold coefficient vectors.
[ ]:
u = fes.TrialFunction() # symbolic object
v = fes.TestFunction() # symbolic object
gfu = GridFunction(fes) # solution
Alternately, you can get both the trial and test variables at once:
[ ]:
u, v = fes.TnT()
5. Define and assemble linear and bilinear forms:#
[ ]:
a = BilinearForm(fes)
a += grad(u)*grad(v)*dx
a.Assemble()
f = LinearForm(fes)
f += x*v*dx
f.Assemble();
Alternately, we can do one-liners:
[ ]:
a = BilinearForm(grad(u)*grad(v)*dx).Assemble()
f = LinearForm(x*v*dx).Assemble()
You can examine the linear system in more detail:
[ ]:
print(f.vec)
[ ]:
print(a.mat)
6. Solve the system:#
[ ]:
gfu.vec.data = \
a.mat.Inverse(freedofs=fes.FreeDofs()) * f.vec
Draw(gfu);
The Dirichlet boundary condition constrains some degrees of freedom. The argument fes.FreeDofs() indicates that only the remaining “free” degrees of freedom should participate in the linear solve.
You can examine the coefficient vector of solution if needed:
[ ]:
print(gfu.vec)
You can see the zeros coming from the zero boundary conditions.
Ways to interact with NGSolve#
A jupyter notebook (like this one) gives you one way to interact with NGSolve. When you have a complex sequence of tasks to perform, the notebook may not be adequate.
You can write an entire python module in a text editor and call python on the command line. (A script of the above is provided in
poisson.py.)python3 poisson.pyIf you want the Netgen GUI, then use
netgenon the command line:netgen poisson.pyYou can then ask for a python shell from the GUI’s menu options (Solve -> Python shell).