Working with CoefficientFunctions¶
A CoefficientFunction is a function which can be evaluated on a mesh, and may be used to provide the coefficient or a right-hand-side to the variational formulation. Because typical finite element procedures iterate over elements, and map integration points from a reference element to a physical element, the evaluation of a CoefficientFunction requires a mapped integration point. The mapped integration point contains the coordinate, as well as the Jacobian, but also element number and region index. This allows the efficient evaluation of a wide class of functions needed for finite element computations.
- class ngsolve.CoefficientFunction¶
A CoefficientFunction (CF) is some function defined on a mesh. Examples are coordinates x, y, z, domain-wise constants, solution-fields, ... CFs can be combined by mathematical operations (+,-,sin(), ...) to form new CFs Parameters:
val : can be one of the following:
- scalar (float or complex):
Creates a constant CoefficientFunction with value val
- tuple of scalars or CoefficientFunctions:
Creates a vector or matrix valued CoefficientFunction, use dims=(h,w) for matrix valued CF
- list of scalars or CoefficientFunctions:
Creates a domain-wise CF, use with generator expressions and mesh.GetMaterials() and mesh.GetBoundaries()
Basic Coefficient Functions¶
Python objects implicitly convertible to float or Complex are implicitly converted to CoefficientFunction if needed. The Cartesian coordinates x, y, and z are pre-defined coordinate coefficient functions:
f = LinearForm(V)
f += SymbolicLFI( x * V.TestFunction())
A CoefficientFunction initialized with a list stores one value per region, depending where the CoefficientFunction is used this is a domain or boundary region. Here one usually uses generator expressions in combination with mesh.GetMaterials()
and mesh.GetBoundaries()
:
alpha_d = {"air" : 1, "box" : 100}
alpha = CoefficientFunction([alpha_d[mat] for mat in mesh.GetMaterials()])
A CoefficientFunction initialized with a tuple gives a vector coefficient function. Components of the coefficient function can be accessed by the bracket operator:
vec_cf = CoefficientFunction( (fx,fy) ) # Creates a vector CF
fx_again = vec_cf[0] # Creates a scalar CF from the first component
A matrix valued CF can be created with the additional dims argument:
mat_cf = CoefficientFunction((f11,f12,f21,f22),dims=(2,2))
ProxyFunctions¶
The finite element spaces provide proxy placeholder functions for Symbolic Integrators when assembling the system matrices the proxy functions are replaced with the FEM basis functions.
- FESpace.TrialFunction(self: ngsolve.comp.FESpace) → object¶
Return a proxy to be used as a trialfunction in Symbolic Integrators
- FESpace.TestFunction(self: ngsolve.comp.FESpace) → object¶
Return a proxy to be used as a testfunction for Symbolic Integrators
Operations on CoefficientFunctions¶
CoefficientFunctions can be combined with algebraic operations (+,-,*,/,**). Math functions are available as well for CoefficientFunctions: sin, cos, tan, exp, log, atan, sqrt, Conj.
Special CoefficientFunctions¶
NGSolve provides special coefficient functions needed for finite element computations. Special coefficient functions are mesh dependent.
You can get the normal vector on an interface with
- specialcf.normal = <bound method PyCapsule.normal of <ngsolve.fem.SpecialCFCreator object>>¶
this can be used in discontinuous Galerkin methods.
The tangential vector can be obtained in the same way:
- specialcf.tangential = <bound method PyCapsule.tangential of <ngsolve.fem.SpecialCFCreator object>>¶
The local mesh size can be obtained by a coefficient function as well:
- specialcf.mesh_size = <ngsolve.fem.CoefficientFunction object>¶
This has application i.e. in hybrid DG methods
.
Additional CoefficientFunctions¶
IfPos¶
The function IfPos provides CoefficientFunctions depending on some condition:
- ngsolve.IfPos(c1: ngfem::CoefficientFunction, then_obj: object, else_obj: object) → ngfem::CoefficientFunction¶
Returns new CoefficientFunction with values then_obj if c1 is positive and else_obj else.
Parameters:
- c1ngsolve.CoefficientFunction
Indicator function
- then_objobject
Values of new CF if c1 is positive, object must be implicitly convertible to ngsolve.CoefficientFunction. See help(
CoefficientFunction
) for information.- else_objobject
Values of new CF if c1 is not positive, object must be implicitly convertible to ngsolve.CoefficientFunction. See help(
CoefficientFunction
) for information.
Parameter CoefficientFunction¶
If you want a CoefficientFunction with a variable parameter, instead of using a new BilinearForm every time you can use a Parameter CF. The parameter can be modified with the Set method and when assembling the BilinearForm again, the updated parameter is used.
- class ngsolve.Parameter¶
CoefficientFunction with a modifiable value
Parameters:
- valuefloat
Parameter value
- Get(self: ngsolve.fem.Parameter) → float¶
return parameter value
- Set(self: ngsolve.fem.Parameter, value: float) → None¶
Modify parameter value.
Parameters:
- valuedouble
input scalar
BSpline CoefficientFunction¶
You can create a BSpline as a CoefficientFunction as well. BSplines are differentiable and integrable:
- class ngsolve.BSpline¶
BSpline of arbitrary order
Parameters:
- orderint
order of the BSpline
- knotslist
list of float
- valslist
list of float
- Differentiate(self: ngsolve.fem.BSpline) → ngsolve.fem.BSpline¶
Differentiate the BSpline
- Integrate(self: ngsolve.fem.BSpline) → ngsolve.fem.BSpline¶
Integrate the BSpline
Compiling a CoefficientFunctions¶
- CoefficientFunction.Compile(self: ngsolve.fem.CoefficientFunction, realcompile: bool = False, maxderiv: int = 2, wait: bool = False) → ngsolve.fem.CoefficientFunction¶
Compile list of individual steps, experimental improvement for deep trees
Parameters:
- realcompilebool
True -> Compile to C++ code
- maxderivint
input maximal derivative
- waitbool
True -> Waits until the previous Compile call is finished before start compiling
Evaluating CoefficientFunctions¶
Sometimes it can be useful to evaluate a CoefficientFunction. Since some CoefficientFunctions are only defined on the mesh (like a GridFunction) this can only be done with information from the mesh. For this we request a mapped integration point from the mesh and plug it into the CoefficientFunction:
>>> mip = mesh(0.2,0.4)
>>> cf = x*x*y
>>> cf(mip)
0.016000000000000004