Since 2.9 release

General changes

Updating to newer version of UFL

Currently we are compatible with the 2022 version of UFL but planning to move to the 2024 or newer versions soon. This requires some minor changes to the usage of UFL shown in the tutorial:

  • Replace `ufl.SpatialCoordinate(ufl.triangle)` with `ufl.SpatialCoordinate(dune.ufl.domain(2))`

  • Replace `ufl.atan_2` with `ufl.atan2`

There is still an issue with dune.ufl.Constant which needs to be resolved before we can update the package requirements of dune-fem.

New clone method on spaces

This new method allows to build a new space from a given space with some parameters changed, e.g., using a different dimRange, order or storage.

Dirichlet BCs for operators with different range/domain space

It is now possible to add Dirichlet boundary conditions to operators with different range and domain spaces. Boundary conditions have to be compatible with the range space and rows for dofs on the boundary are set to zero in the assembled matrix. The boundary values provided in the dune.ufl.DirichletBC structure are ignored. Applying the operator $L$ to a $u$ in the domain space simply results in a function $v$ in the range space with $v_i=0$ for all degrees of freedom on the boundary.

Changes to solver parameter keys

The parameter keys used to customize the (non)-linear solvers during the construction of a scheme have been changed. Keys with prefix newton.linear have been replaced with linear while parameters with only a newton prefix are now prefixed with nonlinear. So

{"newton.tolerance":1e-8,
 "newton.verbose":True,
 "newton.linear.preconditioning.method":"jacobi",
 "newton.linear.tolerance":1e-10}

should be replace with

{"nonlinear.tolerance":1e-8,
 "nonlinear.verbose":True,
 "linear.preconditioning.method": "jacobi",
 "linear.tolerance":1e-10}

In addition the parameter newton.linear.tolerance.strategy has been replaced with nonlinear.forcing.

Warning

while all the changes are indendet to be backwards compatible with suitable deprecation warnings printed if old keys are still used, it is possible that some corner cases were missed.

Added a solve method to the linearization of schemes

The object returned by the linear method on a scheme now has a solve method.

Updated dune.fem.scheme.linearized

This function takes a scheme and returns a scheme implementing the first two terms of its Taylor expansion. Like the object returned by scheme.linear this has a solve method.

Add operator __len__ to spaces.

This allows to use the typical Python built-in function len(space).

Add DirichletIndices method to operators and schemes

This method returns the indices of dofs on the Dirichlet boundary. By default all indices are returned but a boundary id can be provided to only return a subset of indices.

Extension of the jacobian method on schemes and operators

In addition to the existing jacobian method taking a linear operator object a new overload is available taking a further rhs argument which returns the functional part of the operator as though it was on the right hand side. So if the ufl form is a(u,v)=l(v) then this method assembles both the system matrix A and the right hand side functional (rhs needs to be a discrete function over the range space of the operator).

Set non default quadratures in quadrature rules

At the moment only dx(metadata={"quadrature_rule":"lumped"}) is available that sets a quadrature that can be used for lower order Lagrange elements to use lumped mass matrices. More general cases will be added soon.

Assembly function

The function dune.fem.assemble takes a ufl form and returns their assembled versions.

  • If the first argument has arity zero (i.e. has no trial or test function) then the intral will be computed, e.g., dune.fem.assemble(dot.ufl(x,x)*ufl.dx, gridView=gv) where x is a ufl.SpatialCoordinate. Note that the grid view needs to be provided since it can not be determined from the ufl form. It can be omitted if the form contains for example a grid function, e.g., dune.fem.assemble(inner(grad(u_h),grad(u_h))*dx). Note that forms are always scalar valued. For vector valued integration use dune.fem.integrate.

  • If the first function has arity one (i.e. contains a test function) then the functional we be assembled returning a discrete function: dune.fem.assemble(f*v*dx). Again if possible the grid view will be extracted from the form.

  • If the first argument has arity two (i.e. contains both test and trial function) then the function assembles the system matrix and (if the form contains a arity one part) the right hand side functional: dune.fem.assemble(u*v*dx - f*v*dx) returns the linear operator for u*v*dx and the discrete function for f*v*dx. If the form does not contain an arity one part, e.g., dune.fem.assemble(u*v*dx) returns only the linear operator.

Improved Pickling support

Discrete functions can now be directly pickled and will be reconstructed with the grid on load. The required JIT modules will also be generated on load if they do not exist on the system. For this to work one needs to use the dune.common.pickle module as discussed in the Section on pickling.

Paraview Reader

A paraview reader is added to the dune.fem package which can directly work with the pickled files. UFL transformations can be applied after loading the discrete functions in paraview. For paraview to find the reader a environment variable needs to be set. The correct command can be found running python -m dune.fem reader or the correct env variable can be set by running

export PV_PLUGIN_PATH=`python -m dune.fem readerpath`

See merge request 581.

Sampler along Boundary

In addition to a lineSample and pointSample the dune.fem.utility module now also provides a boundarySample function

boundarySample(gridFunction, boundaryId=None, boundaryDomain=None, order=-1)

See the help text of dune.fem.utility.boundarySample for details.

Point sampler returning (entity,local coordinate) for given global coordinate

Improved Dirichlet Boundary Indexing in UFL Forms

To restrict a domain integral to a part of the boundary given by a boundary id one can now use ufl.ds(id) or multiply the integral with dune.ufl.BoundaryId(space) (see the section on Mixing different types of boundary conditions.

Note: we have also unified the default boundary ids when using a dune.grid.cartesianDomain to always be the same as when using dune.grid.structuredGrid. This was not the case with some grid implementations.

Eisenstatwalker in Newton method

Set parameter nonlinear.forcing to eisenstatwalker during scheme construction.

Storage setup in scheme

When creating a scheme one can now select the storage for the underlying linear solver backend. This is demonstrated in the Monolithic solver (Lagrange and DG formulation) where a scheme with storage istl is created while the discrete spaces have storage numpy.

Methods on space to return local mapper and for evaluating basis functions

New methods have been added to the space to return the basis functions evaluated at some given point. Also a new method has been added to return the global indices of the dofs attached to a given element. See section on discrete spaces.

Logging of parameter values

See section on internal solvers.

Deprecation of the rhs argument for scheme.solve

The rhs argument has been replaced with rightHandSide and has to be provided as a named argument.

Warning

schemes with Dirichlet constraints did not handle a given rhs argument correctly. The returned target was always equal to the Dirichlet conditions on the boundary independent of the values in the provided rhs. This has now been corrected so that target=rhs+dirichlet at the boundary which is consistent with the evaluation of the operator (operator(u,w) results in w=u-dirichlet on the boundary). To obtain the old behavior add scheme.setConstraints(0,rhs) before calling scheme.solve(target=target,rightHandSide=rhs). No changes are needed for schemes without Dirichlet boundary constraints. See merge request 661.

Deprecation of dune.fem.function.integrate

This function is deprecated. Use dune.fem.integrate instead. The ufl expression to integrate is now the first argument. If possible the grid view and order will be determined from the ufl expression if possible but can be passed in as extra arguments similar to the old function: dune.fem.function.integrate(grid,expression,order=None) with dune.fem.integrate(expression,gridView=None,order=None).

Deprecation of dune.fem.function.uflFunction and dune.fem.function.cppFunction

To construct different versions of grid function use dune.fem.function.gridFunction. The signature of the new function is def gridFunction(expr=None,gridView=None,*,name=None,order=None, fctName=None, view=None, **kwargs).

  • If expr is a gridView or is None (and the gridView is not None) then the function acts like the decorator from 2.9 version. Note that the old view parameter is also deprecated and replaced by gridView.

    Warning

    @gridFunction(gv,"name",5) will not work anymore since name and order need to be named arguments. So this will lead to an error instead of a deprecation warning and needs to be replaced by @gridFunction(gv,name="name",order=5).

  • If expr is a ufl expression then this functions acts like the uflFunction from 2.9 with the old ufl argument being replaced by expr. In addition gridView and order can be omitted in case where the grid view can be determined from the ufl expression (e.g. the expression contains a grid function). The order will be determined heuristically from the expression if omitted.

  • If expr is a string then the function acts as like the old cppFunction. If fctName is not provided the name argument is assumed to also be the name of the cpp function to export.

Deprecation of dune.ufl.expression2GF

Use the extended version of dune.fem.function.gridFunction instead.

Deprecation of dune.fem.operator.linear

Use the new method linear on the schemes or operators. Note that the new method has no ubar argument (use the jacobian method to assemble around a non zero value). The method has a new bool argument assemble which is True by default. Passing in assemble=False returns an empty linear operator which can be passed into jacobian to be filled.

Deprecation of dune.fem.space.combinedSpace

This has been renamed to dune.fem.space.compositeSpace and will be deprecated soon.

DUNE-VEM

Dirichlet constraints for conforming $C^1$ space

The conforming $C^1$ VEM space can now be used with all types of Dirichlet boundary conditions, i.e., only setting the values, only setting the normal derivative, or setting both.

Warning

setting only value or derivative dofs only works on axis aligned boundaries.

Intersection integrals now use the edge projections

Integrals over the intersection (ds or dS in UFL) now correctly use the edge projections instead of the element projections, so an integral of the jump of a VEM function in a continuous VEM space will now be zero as expected.

https://zenodo.org/badge/DOI/10.5281/zenodo.3706994.svg