In NGSolve, CoefficientFunctions
are representations of functions defined on the computational domain $\Omega$. Examples are expressions of coordinate variables $x, y, z$ and functions that are constant on subdomains. Much of the magic behind the seamless integration of NGSolve with python lies in CoefficientFunctions
. This tutorial introduces you to them.
import netgen.gui
%gui tk
from netgen.geom2d import unit_square
from ngsolve import *
mesh = Mesh (unit_square.GenerateMesh(maxh=0.2))
myfunc = x*(1-x)
myfunc # You have just created a CoefficientFunction
x # This is a built-in CoefficientFunction
Use the mesh
to visualize the function.
Draw(myfunc, mesh, "firstfun")
mip = mesh(0.2, 0.2)
myfunc(mip)
Note that myfunc(0.2,0.3)
does not work: You need to give points in the form of MappedIntegrationPoint
s like mip
above. The mesh
knows how to produce them.
pts = [(0.1*i, 0.2) for i in range(10)]
vals = [myfunc(mesh(*p)) for p in pts]
for p,v in zip(pts, vals):
print("point=(%3.2f,%3.2f), value=%6.5f"
%(p[0], p[1], v))
We may plot the restriction of the CoefficientFunction
on a line using matplotlib.
import matplotlib.pyplot as plt
px = [0.01*i for i in range(100)]
vals = [myfunc(mesh(p,0.5)) for p in px]
plt.plot(px,vals)
plt.xlabel('x')
plt.show()
CoefficientFunction
¶We may Set
a GridFunction
using a CoefficientFunction
:
fes = H1(mesh, order=1)
u = GridFunction(fes)
u.Set(myfunc)
Draw(u) # Cf.: Draw(myfunc, mesh, "firstfun")
The Set
method interpolates myfunc
to obtain the grid function u
.
Set
does an Oswald-type interpolation as follows:
myfunc
in $L^2$ on each mesh element;CoefficientFunction
¶We can numerically integrate the function using the mesh:
Integrate(myfunc, mesh)
There is no facility to directly differentiate a CoefficientFunction
. But you can interpolate it into a GridFunction
and then differentiate the GridFunction
.
u.Set(myfunc)
gradu = grad(u)
gradu[0]
Draw(gradu[0], mesh, 'dx_firstfun')
Obviously the accuracy of this process can be improved for smooth functions by using higher order finite element spaces.
CoefficientFunctions
¶Above, gradu
provided an example of a vector-valued coefficient function. To visualize it, click on Visual
menu in GUI and check Draw Surface Vectors
.
Draw(gradu, mesh, "grad_firstfun")
You can also define vector coefficient expressions directly:
vecfun = CoefficientFunction((-y, sin(x)))
Draw(vecfun, mesh, "vecfun")
Internally, coefficient functions are implemented as an expression tree made from building blocks like x
, y
, sin
, etc., and arithmetic operations.
E.g., the expression tree for myfunc = x*(1-x)
looks like this:
print(myfunc)