NGSolve / ngsolve

Netgen/NGSolve is a high performance multiphysics finite element software. It is widely used to analyze models from solid mechanics, fluid dynamics and electromagnetics. Due to its flexible Python interface new physical equations and solution algorithms can be implemented easily.
https://ngsolve.org/
GNU Lesser General Public License v2.1
421 stars 77 forks source link

CF(0) * TestFunction does not have TestFunction #73

Open CAarset opened 1 month ago

CAarset commented 1 month ago

Ran into an unexpected issue applying a PDE solve to a coefficient-function-based source: The code

v = fes.TestFunction() LinearForm(f v dx)

works for any CoefficientFunction f except f = CoefficientFunction(0), in which case it fails with error "Linearform must have TestFunction". The reason is apparent from

print(CoefficientFunction(0) * VV)

returning ZeroCoefficientFunction, while

print(CoefficientFunction(1) * VV)

returns

coef binary operation '*', real coef unary operation ' ', real coef 1, real coef test-function diffop = Id, real

While this makes some sense from a compression perspective (0 times anything is 0), it makes applying PDE solves to a CF prone to unexpected failure, and notionally it also removes the obvious way to define a 0-LinearForm. Note that interpolating CoefficientFunction(0) to a GridFunction on the fes does work and createas a 0-LinearForm as expected, but is a step I'd rather bypass.

Minimal non-working example attached. CF0.txt

Horep commented 3 weeks ago

As a workaround, you can wrap your float with Parameter. This prevents the error from occurring, i.e. write f0 = CoefficientFunction(Parameter(0)). The other way to produce a 0-LinearForm is to write LinearForm(fes) which will assemble to a vector of zeroes.