Symbolic Integrators

Symbolic integrators can be used to define arbitrary (multi-)physical problems. The variational formulation of the (non-)linear problem can be implemented in a very natural way. Examples of the usage of these integrators can i.e. be found in Navier Stokes Equation, Magnetic Field or Nonlinear elasticity.

The finite element space provides placeholder coefficient functions by the methods TestFunction and TrialFunction. They have implemented canonical derivatives and traces of the finite element space and insert the basis functions of the FESpace during assembling of the system matrices.

Linear Problems

For linear problems we use the function Assemble of the BilinearForm to assemble the matrix and vector. For example for the \(H_\text{curl}\) linear problem

\[\begin{split}\int_\Omega \mu^{-1} \nabla \times u \cdot \nabla \times v + 10^{-6} u \cdot v \, dx = \int_C \begin{pmatrix} y \\ -x \\ 0 \end{pmatrix} \cdot v \, dx\end{split}\]

from example Magnetic Fields we have to define the space


fes = HCurl(mesh, order=4, dirichlet="outer", nograds = True)

and the BilinearForm


a = BilinearForm(fes, symmetric=True)
a += nu*curl(u)*curl(v)*dx + 1e-6*nu*u*v*dx

a.Assemble()

as well as the LinearForm

The argument of the symbolic integrator must be a coefficient function depending linearily on the test and trial function.

BilinearForm.Assemble(self: ngsolve.comp.BilinearForm, reallocate: bool = False)ngsolve.comp.BilinearForm

Assemble the bilinear form.

Parameters:

reallocatebool

input reallocate

Nonlinear Problems

If your left hand side of the variational formulation is nonlinear there are multiple ways to get a discretisation, depending on what you want.

Apply

The function Apply applies the formulation to the given BaseVector. You can get a BaseVector form your GridFunction with GridFunction.vec. The output vector can be created with BaseVector.CreateVector.

BilinearForm.Apply(self: ngsolve.comp.BilinearForm, x: ngsolve.la.BaseVector, y: ngsolve.la.BaseVector)None

Applies a (non-)linear variational formulation to x and stores the result in y.

Parameters:

xngsolve.BaseVector

input vector

yngsolve.BaseVector

output vector

AssembleLinearization

For a variational formulation

\[\int_\Omega f(u) v \, dx\]

the method AssembleLinearization computes

\[\int_\Omega f'(u_\text{lin}) u v \, dx\]

with automatic differentiation of \(f(u)\) and an input BaseVector \(u_\text{in}\).

BilinearForm.AssembleLinearization(self: ngsolve.comp.BilinearForm, ulin: ngsolve.la.BaseVector, reallocate: bool = False)None

Computes linearization of the bilinear form at given vecor.

Parameters:

ulinngsolve.la.BaseVector

input vector

Assemble

You can do your own linearization as well using Assemble and a GridFunction as a CoefficientFunction in your integrator. Let gfu_old be this gridfunction then

a = BilinearForm(fes)
a += SymbolicBFI(gfu_old * u * v)

will be a linearization for

\[\int_\Omega u^2 v \, dx\]

Every time you call Assemble the bilinearform is updated with the new values of the GridFunction.

Symbolic Energy

SymbolicEnergy can be used to solve a minimization problem. In this tutorial we show how to solve the nonlinear problem

\[\min_{u \in V} 0.05 \nabla u + u^4 - 100u\]

For this we use SymbolicEnergy:


a = BilinearForm (V, symmetric=False)
a += Variation( (0.05*grad(u)*grad(u) + u*u*u*u - 100*u)*dx )

from the GridFunction we can create new BaseVector:

With this we can use AssembleLinearization to do a Newton iteration to solve the problem:


for it in range(20):
    print ("Newton iteration", it)
    print ("energy = ", a.Energy(u.vec))

    a.Apply (u.vec, res)
    a.AssembleLinearization (u.vec)
    inv = a.mat.Inverse(V.FreeDofs())
    w.data = inv * res
    print ("w*r =", InnerProduct(w,res))

    u.vec.data -= w
    Redraw()