scikit-fem is a lightweight Python 3.7+ library for performing finite element assembly.

Overview

scikit-fem is a lightweight Python 3.7+ library for performing finite element assembly. Its main purpose is the transformation of bilinear forms into sparse matrices and linear forms into vectors.

Features:

  • minimal dependencies, no compiled code
  • meshes: 1D, tri, quad, tet, hex
  • elements: 1D and quad (any order), tri (order < 5), tet and hex (order < 3)
  • special elements: H(div), H(curl), MINI, Crouzeix-Raviart, Argyris, Morley, ...
  • conforming adaptive mesh refinement

Installation

The most recent release can be installed simply by

pip install scikit-fem[all]

Specifying [all] includes meshio for loading and saving to external formats, and matplotlib for visualization. The minimal dependencies are numpy and scipy. You can also try the library in browser through Google Colab.

Examples

Solve the Poisson problem (see also ex01.py):

from skfem import *
from skfem.helpers import dot, grad

# create the mesh
m = MeshTri().refined(4)
# or, with your own points and elements:
# m = MeshTri(points, elements)

e = ElementTriP1()
basis = Basis(m, e)  # shorthand for CellBasis

@BilinearForm
def laplace(u, v, _):
    return dot(grad(u), grad(v))

@LinearForm
def rhs(v, _):
    return 1.0 * v

A = asm(laplace, basis)
b = asm(rhs, basis)
# or:
# A = laplace.assemble(basis)
# b = rhs.assemble(basis)

# Dirichlet boundary conditions
A, b = enforce(A, b, D=m.boundary_nodes())

# solve the linear system
x = solve(A, b)

# plot the solution
from skfem.visuals.matplotlib import plot, savefig
plot(m, x, shading='gouraud', colorbar=True)
savefig('solution.png')

Meshes can be initialized manually, loaded from external files using meshio, or created with the help of special constructors:

import numpy as np
from skfem import MeshLine, MeshTri, MeshTet

mesh = MeshLine(np.array([0., .5, 1.]))
mesh = MeshTri(
    np.array([[0., 0.],
              [1., 0.],
              [0., 1.]]).T,
    np.array([[0, 1, 2]]).T,
)
mesh = MeshTri.load("docs/examples/meshes/square.msh")  # requires meshio
mesh = MeshTet.init_tensor(*((np.linspace(0, 1, 60),) * 3))

We support many common finite elements. Below the stiffness matrix is assembled using second-order tetrahedra:

from skfem import Basis, ElementTetP2

basis = Basis(mesh, ElementTetP2())  # quadratic tetrahedron
A = laplace.assemble(basis)  # type: scipy.sparse.csr_matrix

More examples can be found in the gallery.

Benchmark

The following benchmark (docs/examples/performance.py) demonstrates the time spent on finite element assembly in comparison to the time spent on linear solve. The given numbers were calculated using a ThinkPad X1 Carbon laptop (7th gen). Note that the timings are only illustrative as they depend on, e.g., the type of element used, the number of quadrature points used, the type of linear solver, and the complexity of the forms. This benchmark solves the Laplace equation using linear tetrahedral elements and the default direct sparse solver of scipy.sparse.linalg.spsolve.

Degrees-of-freedom Assembly (s) Linear solve (s)
4096 0.04805 0.04241
8000 0.09804 0.16269
15625 0.20347 0.87741
32768 0.46399 5.98163
64000 1.00143 36.47855
125000 2.05274 nan
262144 4.48825 nan
512000 8.82814 nan
1030301 18.25461 nan

Documentation

The project is documented using Sphinx under docs/. Built version can be found from Read the Docs. Here are direct links to additional resources:

Getting help

If you encounter an issue you can use GitHub issue tracker. If you cannot find help from the documentation, you can use the GitHub Discussions to ask questions. Try to provide a snippet of code which fails and include also the version of the library you are using. The version can be found as follows:

import skfem; print(skfem.__version__)

Dependencies

The minimal dependencies for installing scikit-fem are numpy and scipy. In addition, many examples use matplotlib for visualization and meshio for loading/saving different mesh file formats. Some examples demonstrate the use of other external packages; see requirements.txt for a list of test dependencies.

Testing

The tests are run by GitHub Actions. The Makefile in the repository root has targets for running the testing container locally using docker. For example, make test_py38 runs the tests using py38 branch from kinnala/scikit-fem-docker-action. The releases are tested in kinnala/scikit-fem-release-tests.

Licensing

The contents of skfem/ and the PyPI package scikit-fem are licensed under the 3-clause BSD license. Some examples under docs/examples/ or snippets in the documentation may have a different license.

Acknowledgements

This project was started while working under a grant from the Finnish Cultural Foundation. Versions 2.0.0+ were prepared while working in a project funded by the Academy of Finland. The approach used in the finite element assembly has been inspired by the work of A. Hannukainen and M. Juntunen.

Contributing

We are happy to welcome any contributions to the library. Reasonable projects for first timers include:

  • Reporting a bug
  • Writing an example
  • Improving the tests
  • Finding typos in the documentation.

By contributing code to scikit-fem, you are agreeing to release it under BSD-3-Clause, see LICENSE.md.

Citing the library

We appreciate if you cite the following article in your scientific works:

@article{skfem2020,
  doi = {10.21105/joss.02369},
  year = {2020},
  volume = {5},
  number = {52},
  pages = {2369},
  author = {Tom Gustafsson and G. D. McBain},
  title = {scikit-fem: A {P}ython package for finite element assembly},
  journal = {Journal of Open Source Software}
}

Changelog

The format is based on Keep a Changelog, and this project adheres to Semantic Versioning with respect to documented and/or tested features.

Unreleased

[5.0.0] - 2021-11-21

  • Changed: meshio is now an optional dependency
  • Changed: ElementComposite uses DiscreteField() to represent zero
  • Changed: Better string __repr__ for Basis and DofsView
  • Added: Support more argument types in Basis.get_dofs
  • Added: Version information in skfem.__version__
  • Added: Preserve Mesh.boundaries during uniform refinement of MeshLine1, MeshTri1 and MeshQuad1
  • Fixed: Refinement of quadratic meshes will now fall back to the refinement algorithm of first-order meshes instead of crashing
  • Fixed: Edge cases in the adaptive refine of MeshTet1 that failed to produce a valid mesh
  • Deprecated: Basis.find_dofs in favor of Basis.get_dofs
  • Deprecated: Merging DofsView objects via + and | in favor of using np.hstack

[4.0.1] - 2021-10-15

  • Fixed: MappingIsoparametric can now be pickled

[4.0.0] - 2021-09-27

  • Added: Mesh.save/Mesh.load now exports/imports Mesh.subdomains and Mesh.boundaries
  • Added: Mesh.load now optionally writes any mesh data to a list passed via the keyword argument out, e.g., out=data where data = ['point_data']
  • Added: Mesh.load (and skfem.io.meshio.from_file) now supports the additional keyword argument force_meshio_type for loading mesh files that have multiple element types written in the same file, one element type at a time
  • Added: asm will now accept a list of bases, assemble the same form using all of the bases and sum the result (useful for jump terms and mixed meshes, see Example 41)
  • Added: Mesh.with_boundaries now allows the definition of internal boundaries/interfaces via the flag boundaries_only=False
  • Added: MeshTri1DG, MeshQuad1DG, MeshHex1DG, MeshLine1DG; new mesh types for describing meshes with a discontinuous topology, e.g., periodic meshes (see Example 42)
  • Added: ElementHexDG for transforming hexahedral H1 elements to DG/L2 elements.
  • Added: ElementTriP1DG, ElementQuad1DG, ElementHex1DG, ElementLineP1DG; shorthands for ElementTriDG(ElementTriP1()) etc.
  • Added: ElementTriSkeletonP0 and ElementTriSkeletonP1 for defining Lagrange multipliers on the skeleton mesh (see Example 40)
  • Added: TrilinearForm for assembling a sparse 3-tensor, e.g., when dealing with unknown material data
  • Added: MeshTri.oriented for CCW oriented triangular meshes which can be useful for debugging or interfacing to external tools
  • Added: partial support for MeshWedge1 and ElementWedge1, the lowest order wedge mesh and element
  • Added: ElementTriP3, cubic triangular Lagrange element
  • Added: ElementTriP4, quartic triangular Lagrange element
  • Added: ElementTri15ParamPlate, 15-parameter nonconforming triangular element for plates
  • Added: ElementTriBDM1, the lowest order Brezzi-Douglas-Marini element
  • Added: Mesh.draw().show() will now visualize any mesh interactively (requires vedo)
  • Added: Adaptive refinement for MeshTet1
  • Fixed: MappingIsoparametric is now about 2x faster for large meshes thanks to additional caching
  • Fixed: MeshHex2.save did not work properly
  • Fixed: Mesh.load ignores unparseable cell_sets inserted by meshio in MSH 4.1
  • Changed: Mesh string representation is now more informative
  • Changed: Form.assemble no more allows keyword arguments with list or dict type: from now on only DiscreteField or 1d/2d ndarray objects are allowed and 1d ndarray is passed automatically to Basis.interpolate for convenience
  • Changed: MeshLine is now a function which initializes MeshLine1 and not an alias to MeshLine1
  • Changed: FacetBasis is now a shorthand for BoundaryFacetBasis and no longer initializes InteriorFacetBasis or MortarFacetBasis if the keyword argument side is passed to the constructor
  • Removed: the deprecated Mesh.define_boundary method
  • Removed: the unused Mesh.validate attribute

[3.2.0] - 2021-08-02

  • Added: ElementTriCCR and ElementTetCCR, conforming Crouzeix-Raviart finite elements
  • Fixed: Mesh.mirrored returned a wrong mesh when a point other than the origin was used
  • Fixed: MeshLine constructor accepted only NumPy arrays and not plain Python lists
  • Fixed: Mesh.element_finder (and CellBasis.probes, CellBasis.interpolator) was not working properly for a small number of elements (<5) or a large number of input points (>1000)
  • Fixed: MeshTet and MeshTri.element_finder are now more robust against degenerate elements
  • Fixed: Mesh.element_finder (and CellBasis.probes, CellBasis.interpolator) raises exception if the query point is outside of the domain

[3.1.0] - 2021-06-18

  • Added: Basis, a shorthand for CellBasis
  • Added: CellBasis, a new preferred name for InteriorBasis
  • Added: BoundaryFacetBasis, a new preferred name for ExteriorFacetBasis
  • Added: utils.penalize, an alternative to condense and enforce for essential boundary conditions
  • Added: InteriorBasis.point_source, with ex38
  • Added: ElementTetDG, similar to ElementTriDG for tetrahedral meshes
  • Fixed: MeshLine1.element_finder

[3.0.0] - 2021-04-19

  • Added: Completely rewritten Mesh base class which is "immutable" and uses Element classes to define the ordering of nodes; better support for high-order and other more general mesh types in the future
  • Added: New quadratic mesh types: MeshTri2, MeshQuad2, MeshTet2 and MeshHex2
  • Added: InteriorBasis.probes; like InteriorBasis.interpolator but returns a matrix that operates on solution vectors to interpolate them at the given points
  • Added: More overloads for DiscreteField, e.g., multiplication, summation and subtraction are now explicitly supported inside the form definitions
  • Added: MeshHex.to_meshtet for splitting hexahedra into tetrahedra
  • Added: MeshHex.element_finder for interpolating finite element solutions on hexahedral meshes via InteriorBasis.interpolator
  • Added: Mesh.with_boundaries, a functional replacement to Mesh.define_boundary, i.e. defining boundaries via Boolean lambda function
  • Added: Mesh.with_subdomains for defining subdomains via Boolean lambda function
  • Added: skfem.utils.projection, a replacement of skfem.utils.project with a different, more intuitive order of arguments
  • Added: skfem.utils.enforce for setting essential boundary conditions by changing matrix rows to zero and diagonals to one.
  • Deprecated: skfem.utils.project in favor of skfem.utils.projection
  • Deprecated: Mesh.define_boundary in favor of Mesh.with_boundaries
  • Removed: Mesh.{refine,scale,translate}; the replacements are Mesh.{refined,scaled,translated}
  • Removed: skfem.models.helpers; available as skfem.helpers
  • Removed: DiscreteField.{f,df,ddf,hod}; available as DiscreteField.{value,grad,hess,grad3,...}
  • Removed: Python 3.6 support
  • Changed: Mesh.refined no more attempts to fix the indexing of Mesh.boundaries after refine
  • Changed: skfem.utils.solve now uses scipy.sparse.eigs instead of scipy.sparse.eigsh by default; the old behavior can be retained by explicitly passing solver=solver_scipy_eigs_sym()
  • Fixed: High memory usage and other small fixes in skfem.visuals.matplotlib related to 1D plotting

[2.5.0] - 2021-02-13

  • Deprecated: side keyword argument to FacetBasis in favor of the more explicit InteriorFacetBasis and MortarFacetBasis.
  • Added: InteriorFacetBasis for integrating over the interior facets, e.g., evaluating error estimators with jumps and implementing DG methods.
  • Added: MortarFacetBasis for integrating over the mortar mesh.
  • Added: InteriorBasis.with_element for reinitializing an equivalent basis that uses a different element.
  • Added: Form.partial for applying functools.partial to the form function wrapped by Form.
  • Fixed: Include explicit Python 3.9 support.

[2.4.0] - 2021-01-20

  • Deprecated: List and tuple keyword argument types to asm.
  • Deprecated: Mesh2D.mirror in favor of the more general Mesh.mirrored.
  • Deprecated: Mesh.refine, Mesh.scale and Mesh.translate in favor of Mesh.refined, Mesh.scaled and Mesh.translated.
  • Added: Mesh.refined, Mesh.scaled, and Mesh.translated. The new methods return a copy instead of modifying self.
  • Added: Mesh.mirrored for mirroring a mesh using a normal and a point.
  • Added: Functional now supports forms that evaluate to vectors or other tensors.
  • Added: ElementHex0, piecewise constant element for hexahedral meshes.
  • Added: FacetBasis.trace for restricting existing solutions to lower dimensional meshes on boundaries or interfaces.
  • Fixed: MeshLine.refined now correctly performs adaptive refinement of one-dimensional meshes.

[2.3.0] - 2020-11-24

  • Added: ElementLineP0, one-dimensional piecewise constant element.
  • Added: skfem.helpers.curl now calculates the rotated gradient for two-dimensional elements.
  • Added: MeshTet.init_ball for meshing a ball.
  • Fixed: ElementQuad0 was not compatible with FacetBasis.

[2.2.3] - 2020-10-16

  • Fixed: Remove an unnecessary dependency.

[2.2.2] - 2020-10-15

  • Fixed: Make the preconditioner in TestEx32 more robust.

[2.2.1] - 2020-10-15

  • Fixed: Remove tests from the PyPI distribution.

[2.2.0] - 2020-10-14

  • Deprecated: L2_projection will be replaced by project.
  • Deprecated: derivative will be replaced by project.
  • Added: MeshTet.element_finder and MeshLine.element_finder for using InteriorBasis.interpolator.
  • Added: ElementTriCR, the nonconforming Crouzeix-Raviart element for Stokes flow.
  • Added: ElementTetCR, tetrahedral nonconforming Crouzeix-Raviart element.
  • Added: ElementTriHermite, an extension of ElementLineHermite to triangular meshes.
  • Fixed: Fix Mesh.validate for unsigned Mesh.t.

[2.1.1] - 2020-10-01

  • Fixed: Further optimizations to Mesh3D.boundary_edges: tested to run on a laptop with over 10 million elements.

[2.1.0] - 2020-09-30

  • Added: ElementHex2, a triquadratic hexahedral element.
  • Added: MeshTri.init_circle, constructor for a circle mesh.
  • Fixed: Mesh3D.boundary_edges (and, consequently, Basis.find_dofs) was slow and used lots of memory due to an exhaustive search of all edges.

[2.0.0] - 2020-08-21

  • Deprecated: project will only support functions like lambda x: x[0] instead of lambda x, y, z: x in the future.
  • Added: Support for complex-valued forms: BilinearForm and LinearForm now take an optional argument dtype which defaults to np.float64 but can be also np.complex64.
  • Added: Dofs.__or__ and Dofs.__add__, for merging degree-of-freedom sets (i.e. Dofs objects) using | and + operators.
  • Added: Dofs.drop and Dofs.keep, for further filtering the degree-of-freedom sets
  • Removed: Support for old-style decorators bilinear_form, linear_form, and functional (deprecated since 1.0.0).
  • Fixed: FacetBasis did not initialize with ElementQuadP.

[1.2.0] - 2020-07-07

  • Added: MeshQuad._splitquads aliased as MeshQuad.to_meshtri.
  • Added: Mesh.__add__, for merging meshes using + operator: duplicated nodes are joined.
  • Added: ElementHexS2, a 20-node quadratic hexahedral serendipity element.
  • Added: ElementLineMini, MINI-element for one-dimensional mesh.
  • Fixed: Mesh3D.boundary_edges was broken in case of hexahedral meshes.
  • Fixed: skfem.utils.project did not work for ElementGlobal.

[1.1.0] - 2020-05-18

  • Added: ElementTetMini, MINI-element for tetrahedral mesh.
  • Fixed: Mesh3D.boundary_edges incorrectly returned all edges where both nodes are on the boundary.

[1.0.0] - 2020-04-22

  • Deprecated: Old-style form constructors bilinear_form, linear_form, and functional.
  • Changed: Basis.interpolate returns DiscreteField objects instead of ndarray tuples.
  • Changed: Basis.interpolate works now properly for vectorial and high-order elements by interpolating all components and higher order derivatives.
  • Changed: Form.assemble accepts now any keyword arguments (with type DiscreteField) that are passed over to the forms.
  • Changed: Renamed skfem.importers to skfem.io.
  • Changed: Renamed skfem.models.helpers to skfem.helpers.
  • Changed: skfem.utils.solve will now expand also the solutions of eigenvalue problems.
  • Added: New-style form constructors BilinearForm, LinearForm, and Functional.
  • Added: skfem.io.json for serialization of meshes to/from json-files.
  • Added: ElementLinePp, p-th order one-dimensional elements.
  • Added: ElementQuadP, p-th order quadrilateral elements.
  • Added: ElementQuadDG for transforming quadrilateral H1 elements to DG elements.
  • Added: ElementQuadBFS, Bogner-Fox-Schmit element for biharmonic problems.
  • Added: ElementTriMini, MINI-element for Stokes problems.
  • Added: ElementComposite for using multiple elements in one bilinear form.
  • Added: ElementQuadS2, quadratic Serendipity element.
  • Added: ElementLineHermite, cubic Hermite element for Euler-Bernoulli beams.
  • Added: Mesh.define_boundary for defining named boundaries.
  • Added: Basis.find_dofs for finding degree-of-freedom indices.
  • Added: Mesh.from_basis for defining high-order meshes.
  • Added: Basis.split for splitting multicomponent solutions.
  • Added: MortarMapping with basic support for mortar methods in 2D.
  • Added: Basis constructors now accept quadrature keyword argument for specifying a custom quadrature rule.
Comments
  • Add support for facet orientation in FacetBasis

    Add support for facet orientation in FacetBasis

    This PR now acknowledges the existence of oriented boundaries/interfaces (i.e. sets of facets).

    We can load Gmsh file with an oriented interface in the middle:

    In [1]: from skfem import *                                                     
    In [2]: m = MeshTri.load('docs/examples/meshes/interface.msh')                  
    In [3]: m.boundaries                                                            
    Out[3]: {'interfacee': OrientedBoundary([ 7, 11, 50, 55, 60, 65, 70, 75])}
    

    Here OrientedBoundary is a subclass of ndarray with the additional attribute ori which is either 0 or 1 for each facet. We can visualize the orientation:

    In [4]: m.draw(boundaries=True).show() 
    

    Figure_1

    The orientation is taken into account in FacetBasis, e.g.,

    In [5]: fb = FacetBasis(m, ElementTriP1(), facets='interfacee')  
    In [7]: from skfem.helpers import dot                                           
    In [8]: Functional(lambda w: dot(w.y.grad, w.n)).assemble(fb, y=m.p[1])         
    Out[8]: 1.0
    

    Obviously there can be multiple oriented facet sets:

    In [13]: m = MeshTri.load('docs/examples/meshes/internal.msh')                  
    In [14]: m.draw(boundaries=True).show()
    

    Figure_2

    There is a facility for creating oriented facet sets around subdomains in Mesh.facets_around:

    In [2]: from skfem import *                                                     
    In [3]: m = MeshTri.load('docs/examples/meshes/internal.msh')                   
    In [4]: M = m.with_boundaries({'left': m.facets_around(lambda x: x[0] < 0)})    
    In [5]: M.draw(boundaries=True).show()
    

    Figure_3

    We can also now orient facet sets constructed with Mesh.facets_satisfying:

    In [1]: from skfem import *
    In [2]: m = MeshTri().refined(4)
    In [3]: M = m.with_boundaries({'mid': m.facets_satisfying(lambda x: x[0] == 0.5,
       ...:  normal=[1, 0])})
    In [4]: M.draw(boundaries=True).show()
    

    Figure_1

    Fixes #821 and fixes #870.

    Missing work: boundary orientations are invalidated by save/load cycle and by refine.

    opened by kinnala 107
  • weak forms for hyperelasticity

    weak forms for hyperelasticity

    @bhaveshshrimali (from #438):

    Thanks! I was starting to take a look at implementing the hyperelasticity demo. The most natural examples to follow seem to be linear elasticity, nonlinear poisson and laplace with inhomogeneous bcs.

    What would be the best place to look for functions like log, inv, det (like in UFL) when writing the weak forms? For instance the stress would have an expression that looks like the following in UFL:

    def firstPKStress(u):
        F = Identity(len(u)) + grad(u) 
        J = det(F)
        return mu * F - mu * inv(F).T + J * (J-1) * inv(F).T
    

    I can see the helper functions eye and grad, which should help in defining F as eye(1, u.shape[0]) + grad(u), but how about the determinant (det) and inverse (inv) ?

    opened by gdmcbain 71
  • Examples providing your own mesh and retriving the basis functions?

    Examples providing your own mesh and retriving the basis functions?

    Are there any examples for supplying your own mesh (for example from scipy Delaunay) and extracting basis functions of some degree?

    For example you might evaluate the basis functions at a set of points for fitting against data as a kind of smoothing spline.

    question 
    opened by cottrell 58
  • skfem.ivp

    skfem.ivp

    From #529 at 2020-12-24:

    Note that we don't have any helper functions for time integration yet. That would be a welcome improvement though, e.g., skfem.ivp or something.

    discussion 
    opened by gdmcbain 54
  • Mesh.save subdomains and boundaries

    Mesh.save subdomains and boundaries

    Should Mesh.save save the .subdomains and .boundaries attributes? In a way that they might be reread by Mesh.load.

    If I/O is to continue to rely on meshio.write and meshio.read (as it should), I think that saving boundaries requires saving cells of type bnd_type (e.g. 'line' for MeshQuad); currently skfem.Mesh.save only saves cells of type Mesh.meshio_type.

    (I'm currently parsing a polyMesh from OpenFOAM and constructing an skfem.Mesh. My parser is very slow so I want to be able to save the mesh while experimenting with the solver.)

    Perhaps instead for the moment I should just construct a meshio.Mesh and meshio.write that, leaving scikit-fem out of it, but this might be something to think about for other cases, e.g. when the mesh has been refined or adapted inside scikit-fem.

    feature request 
    opened by gdmcbain 44
  • Singular matrix when mesh contains points not used by cells.

    Singular matrix when mesh contains points not used by cells.

    Expected behavior is that these points would be ignored (as is the case in other packages, e.g. SfePy, dolfinx). Observed behavior is assembly of a singular matrix and failure to solve.

    Minimum code to reproduce is below, with one line to add/not add the offending point. You can inspect mesh.points and mesh.cells_dict['lines'] and mesh.cells_dict['triangles'] to see this point is in the points array but not used by any line or trinagle.

    import gmsh
    import meshio
    import skfem
    from skfem.helpers import dot, grad
    from skfem.io.meshio import from_meshio
    
    lc0 = .25
    p_center = [.1234, .5678]
    gmsh.initialize()
    try:
        gmsh.model.add('mesh')
        # Define unit square area
        # Corners
        bottom_left = gmsh.model.geo.addPoint(0, 0, 0, lc0)
        bottom_right = gmsh.model.geo.addPoint(1, 0, 0, lc0)
        top_right = gmsh.model.geo.addPoint(1, 1, 0, lc0)
        top_left = gmsh.model.geo.addPoint(0, 1, 0, lc0)
        # Edges
        bottom = gmsh.model.geo.addLine(bottom_left, bottom_right)
        right = gmsh.model.geo.addLine(bottom_right, top_right)
        top = gmsh.model.geo.addLine(top_right, top_left)
        left = gmsh.model.geo.addLine(top_left, bottom_left)
        # Surface
        perimeter = gmsh.model.geo.addCurveLoop([bottom, right, top, left])
        surface = gmsh.model.geo.addPlaneSurface([perimeter])
    #
    #
    # comment out the following line to remove the offending point
    # and clear the warnings/failure 
    #
        gmsh.model.geo.add_point(p_center[0], p_center[1], 0, lc0)
    #
    #
    #
        gmsh.model.geo.synchronize()
        gmsh.model.mesh.generate(dim=2)
        gmsh.write('demo.msh')
    finally:
        gmsh.finalize()
    
    mesh = meshio.read('demo.msh')
    fem_mesh = from_meshio(mesh)
    basis = skfem.Basis(fem_mesh, skfem.ElementTriP1())
    
    @skfem.BilinearForm
    def laplace(u, v, _):
        return dot(grad(u), grad(v))
    
    @skfem.LinearForm
    def rhs(v, _):
        return 1.0 * v
    
    A = laplace.assemble(basis)
    b = rhs.assemble(basis)
    
    D = fem_mesh.boundary_nodes()
    A,b = skfem.enforce(A, b, D=D)  # warning generated here:
    # site-packages\scipy\sparse\_index.py:125: SparseEfficiencyWarning: Changing the
    # sparsity structure of a csr_matrix is expensive. lil_matrix is more efficient.
    #    self._set_arrayXarray(i, j, x)
    #
    
    x = skfem.solve(A, b)  # another warning generated here:
    # site-packages\scipy\sparse\linalg\dsolve\linsolve.py:206: MatrixRankWarning:
    # Matrix is exactly singular
    #    warn("Matrix is exactly singular", MatrixRankWarning)
    
    opened by gatling-nrl 39
  • Experimental autodiff support

    Experimental autodiff support

    Obstacle problem with penalty?

    from skfem import *
    from skfem.experimental.autodiff import NonlinearForm
    from skfem.experimental.autodiff.helpers import *
    import autograd.numpy as np
    
    
    m = MeshTri().refined(5)
    basis = Basis(m, ElementTriP1())
    
    eps = 1e-3
    f0 = -1.
    
    def g(x):
        return np.sin(np.pi * x[0]) - 1.05
    
    
    @NonlinearForm
    def F(u, v, w):
        return (dot(grad(u), grad(v))
                + 1. / w.h ** 2 * (u[0] - g(w.x)) * (u[0] - g(w.x) <= 0) * v[0]
                - f0 * v[0])
    
    
    x = basis.zeros()
    D = basis.get_dofs()
    
    for itr in range(50):
    
        J, f = F.assemble(x, basis)
    
        xprev = x.copy()
        x += 0.99 * solve(*condense(J, -f, D=D))
    
        res = np.linalg.norm(x - xprev)
        print(res)
        if res < 1e-8:
            break
    
    
    basis.plot(x, shading='gouraud', colorbar=True).show()
    

    Minimal surface problem?

    from skfem import *
    from skfem.experimental.autodiff import NonlinearForm
    from skfem.experimental.autodiff.helpers import grad, dot
    import numpy as np
    import autograd.numpy as anp
    
    m = MeshTri().refined(5)
    
    
    @NonlinearForm(hessian=True)
    def energy(u, _):
        return anp.sqrt(1. + dot(grad(u), grad(u)))
    
    
    basis = Basis(m, ElementTriP1())
    
    D = m.boundary_nodes()
    x = np.sin(np.pi * m.p[0])
    
    for itr in range(100):
        J, F = energy.assemble(x, basis)
        x_prev = x.copy()
        x += 0.7 * solve(*condense(J, -F, D=D))
        res = np.linalg.norm(x - x_prev)
        if res < 1e-8:
            break
        if __name__ == "__main__":
            print(res)
    
    if __name__ == "__main__":
        from skfem.visuals.matplotlib import plot3, show
        plot3(m, x)
        show()
    

    Lid driven cavity with Navier-Stokes?

    from skfem import *
    from skfem.experimental.autodiff import NonlinearForm
    from skfem.experimental.autodiff.helpers import *
    import numpy as np
    
    
    m = MeshTri().refined(3)
    basis = Basis(m, ElementVector(ElementTriP2()) * ElementTriP1())
    
    
    @NonlinearForm
    def F(u, p, v, q, w):
        return (ddot(grad(u), grad(v))
                - dot(mul(grad(u), u[0]), v[0])
                - 1e-2 * p[0] * q[0]
                - div(u) * q[0]
                - div(v) * p[0])
    
    
    x = basis.zeros()
    D = basis.get_dofs().all(['u^1^1', 'u^2^1'])
    x[basis.get_dofs('top').all(['u^1^1'])] = 100.
    
    for itr in range(50):
    
        J, f = F.assemble(x, basis)
    
        xprev = x.copy()
        x += solve(*condense(J, -f, D=D))
    
        res = np.linalg.norm(x - xprev)
        print(res)
        if res < 1e-8:
            break
    
    
    (u, ubasis), (p, pbasis) = basis.split(x)
    
    pbasis.plot(p).show()
    ubasis.plot(u).show()
    
    opened by kinnala 37
  • 'DiscreteField' object has no attribute 'reshape'

    'DiscreteField' object has no attribute 'reshape'

    This leads to confusing behavior in forms.

    def some_functional():
        def _form(w):
            # arrays are shaped to group by facet, but we don't care for that here so reshape(2,-1)
            ii = w['ind_p0'].reshape(2,-1).T
            gg = grad(w['ind_p1']).reshape(2,-1).T
            pp = w.x.reshape(2,-1).T  
            nn = w.n.reshape(2,-1).T
            ...
            return w.x[0] # must return *something* of the right shape
        return skfem.Functional(_form)
    
    skfem.asm(
        some_functional(),
        skfem.InteriorFacetBasis(mesh, skfem.ElementTriP1(), side=0),
        ind_p1=ind_p1,
        ind_p0=basis_p1.project(basis_p0.interpolate(ind_p0)),
    )
    

    gg, pp, nn all work, but ii returns 'DiscreteField' object has no attribute 'reshape'.

    The "work around" or maybe the intended usage appears to be

    ii = w['ind_p0'].value.reshape(2,-1).T
    

    and this might be okay, except that

    gg = grad(w['ind_p1']).reshape(2,-1).T
    

    strongly implies (imo) that w['ind_p1'] is the "value". Moreover w[ind_p0]**2 works, and returns (w[ind_p0].value)**2 making me think maybe w[ind_p0] supports any ufunc.

    I think reshape should work.

    opened by gatling-nrl 37
  • Matrix Transforms, Displacements, and Resultant Forces

    Matrix Transforms, Displacements, and Resultant Forces

    Hello All,

    I don't have an issue with scikit-fem so much as I have a general question of use.

    I'd like to take a beam, constrain one end then displace the other end, but not just in translation-- I'd like to perform a 4x4 matrix transform on that end consisting of both rotation and translation then determine the resultant force.

    Every package I've investigated so far only implements the translation without the rotation: Solidworks, Fusion 360, CalculiX, etc.

    Is this something that scikit-fem can do? Would you be so kind as to show me how?

    Thanks!

    question 
    opened by mbatesole 35
  • Using `scikit-fem` on pip: solver remark

    Using `scikit-fem` on pip: solver remark

    Although this is tangentially related to scikit-fem since the primary objective of scikit-fem is assembly and not solving the resulting systems, it maybe important to note that when working with a default pip installation on Windows, there is a significant performance degradation in solving, since the default sparse solver is SuperLU which is super slow. I think the conda installation may not suffer from this since the default solver is UMFPACK (?) which can be faster. There maybe further improvements possible when numpy is compiled against an optimized BLAS but that is something most users wouldn't care for.

    I managed to get Pardiso working on Windows and with pypardiso things are significantly faster. See the timings on example 36 with nthreads=8 (and Mesh.refined(3)) and a slight addition at https://github.com/kinnala/scikit-fem/blob/c07ff9f6d6a51e62eac2777351b20b4549027ada/skfem/utils.py#L102

    namely simply calling pypardiso.spsolve

    def solver_direct_pypardiso(**kwargs) -> LinearSolver:
        """Linear solver provided by Pardiso."""
    
        def solver(A, b):
            return pypardiso.spsolve(A, b)
    
        return solver
    

    Timings

    Pardiso

    Solve time: 0.3992440700531006
    1, norm_du: 45.481410024640475, norm_dp: 14.934126744945504
    Solve time: 0.31940722465515137
    2, norm_du: 16.15926762709032, norm_dp: 20.93859584370273
    Solve time: 0.3229987621307373
    3, norm_du: 3.2471798205287667, norm_dp: 15.053533831714226
    Solve time: 0.35599684715270996
    4, norm_du: 0.4538801110428914, norm_dp: 2.6708267761548887
    Solve time: 0.3192098140716553
    5, norm_du: 0.027865264212750027, norm_dp: 0.1630714233226417
    Solve time: 0.328704833984375
    6, norm_du: 0.002041324778791658, norm_dp: 0.016283691403898768
    Solve time: 0.3122262954711914
    7, norm_du: 1.6033937733037654e-05, norm_dp: 0.0003318753774227752
    Solve time: 0.33296942710876465
    8, norm_du: 1.6887613490335165e-09, norm_dp: 4.315368227582554e-08
    Solve time: 0.35301685333251953
    9, norm_du: 4.4086266413504256e-15, norm_dp: 1.7099154835473912e-14
    Total time: 98.29319977760315
    
    

    SuperLU

    Solve time: 19.864436626434326
    1, norm_du: 45.48141002464045, norm_dp: 14.934126744945653
    Solve time: 17.160115242004395
    2, norm_du: 16.159267627090234, norm_dp: 20.93859584370295
    Solve time: 18.916175842285156
    3, norm_du: 3.247179820528735, norm_dp: 15.053533831714496
    Solve time: 17.143075704574585
    4, norm_du: 0.45388011104290815, norm_dp: 2.670826776154936
    Solve time: 17.078906536102295
    5, norm_du: 0.027865264212752802, norm_dp: 0.16307142332264873
    Solve time: 17.101752281188965
    6, norm_du: 0.0020413247787918446, norm_dp: 0.016283691403900905
    Solve time: 17.19295907020569
    7, norm_du: 1.6033937733110564e-05, norm_dp: 0.00033187537742298715
    Solve time: 17.220698833465576
    8, norm_du: 1.688761318114545e-09, norm_dp: 4.3153681931898745e-08
    Solve time: 16.51562809944153
    9, norm_du: 4.3255823949057516e-15, norm_dp: 1.5848762808515266e-14
    Total time: 251.57028770446777
    
    
    opened by bhaveshshrimali 34
  • periodic boundary conditions

    periodic boundary conditions

    Any ideas on periodic boundary conditions? There's a brief mention in the JOSS paper.

    I'm getting increasingly interested in problems of propagation, for which these often appear in textbook examples, so I'll soon look into throwing something together.

    feature request 
    opened by gdmcbain 34
  • Interpolator with ElementVector

    Interpolator with ElementVector

    import numpy as np
    from skfem import *
    
    mesh = MeshTri()
    basis = Basis(mesh, ElementTriP1())
    basis.interpolator(basis.zeros())(np.array([[[0, 0]]]).swapaxes(0, -1))
    

    works, but

    import numpy as np
    from skfem import *
    
    mesh = MeshTri()
    basis = Basis(mesh, ElementVector(ElementTriP1()))
    basis.interpolator(basis.zeros())(np.array([[[0, 0]]]).swapaxes(0, -1))
    

    doesn't.

    ValueError: row, column, and data array must all be the same length
    
    opened by HelgeGehring 3
  • facet_basis in high order meshes

    facet_basis in high order meshes

    Discussed in https://github.com/kinnala/scikit-fem/discussions/964

    Originally posted by Abelsylowlie November 11, 2022 Hi, all. I try to solve a problem with a boundary term on high order meshes and desire to know if facet_basis is compatible with high order meshes.

    >>> from skfem import MeshTri2, ElementTriP2, FacetBasis
    >>> FacetBasis(MeshTri2.init_circle(), ElementTriP2())
    

    raises

      File in <module>
        FacetBasis(MeshTri2.init_circle(), ElementTriP2())
    
      File ~\anaconda3\lib\site-packages\skfem\assembly\basis\facet_basis.py:91 in __init__
        Y = self.mapping.invF(x, tind=self.tind)
    
      File ~\anaconda3\lib\site-packages\skfem\mapping\mapping_isoparametric.py:153 in invF
        raise Exception(("Newton iteration didn't converge "
    
    Exception: Newton iteration didn't converge up to TOL=1e-12
    ```</div>
    opened by kinnala 0
  • MeshDG not compatible with ElementGlobal

    MeshDG not compatible with ElementGlobal

    Here is minimal example:

    m = MeshLine(np.linspace(0, 1, 10))
    mp = MeshLine1DG.periodic(m, [9], [0])
    basis = Basis(mp, ElementLineHermite())
    

    leads to

    LinAlgError: Singular matrix
    
    opened by kinnala 0
  • Release tests on Windows

    Release tests on Windows

    The release tests are breaking on Windows because of temporary files. Could be that delete=False must be added such as here: https://github.com/appveyor/ci/issues/2547

    opened by kinnala 0
  • Use outward normals on boundary for ElementGlobal DOFs

    Use outward normals on boundary for ElementGlobal DOFs

    I think it would be good if we use outward vectors to initialize the DOFs. I tried to study the codes in element_global.py but I'm not sure how to fix it.

    Let's open another issue for that. It's not completely obvious how to do this because ElementGlobal uses its own ordering independent of Mapping. The reason why it's done like this currently is to make sure that normals for the internal DOFs between two elements point to same directions. That's required for conformity, but boundary normals are not taken into account at all.

    Originally posted by @kinnala in https://github.com/kinnala/scikit-fem/issues/566#issuecomment-782877591

    bug 
    opened by kinnala 0
Releases(8.0.0)
Owner
Tom Gustafsson
A researcher working on computational methods.
Tom Gustafsson
Coursera Machine Learning - Python code

Coursera Machine Learning This repository contains python implementations of certain exercises from the course by Andrew Ng. For a number of assignmen

Jordi Warmenhoven 859 Dec 10, 2022
Predict the demand for electricity (R) - FRENCH

06.demand-electricity Predict the demand for electricity (R) - FRENCH Prédisez la demande en électricité Prérequis Pour effectuer ce projet, vous devr

1 Feb 13, 2022
A concept I came up which ditches the idea of "layers" in a neural network.

Dynet A concept I came up which ditches the idea of "layers" in a neural network. Install Copy Dynet.py to your project. Run the example Install matpl

Anik Patel 4 Dec 05, 2021
Laporan Proyek Machine Learning - Azhar Rizki Zulma

Laporan Proyek Machine Learning - Azhar Rizki Zulma Project Overview Domain proyek yang dipilih dalam proyek machine learning ini adalah mengenai hibu

Azhar Rizki Zulma 6 Mar 12, 2022
A simple application that calculates the probability distribution of a normal distribution

probability-density-function General info An application that calculates the probability density and cumulative distribution of a normal distribution

1 Oct 25, 2022
This is the code repository for Interpretable Machine Learning with Python, published by Packt.

Interpretable Machine Learning with Python, published by Packt

Packt 299 Jan 02, 2023
SageMaker Python SDK is an open source library for training and deploying machine learning models on Amazon SageMaker.

SageMaker Python SDK SageMaker Python SDK is an open source library for training and deploying machine learning models on Amazon SageMaker. With the S

Amazon Web Services 1.8k Jan 01, 2023
MICOM is a Python package for metabolic modeling of microbial communities

Welcome MICOM is a Python package for metabolic modeling of microbial communities currently developed in the Gibbons Lab at the Institute for Systems

57 Dec 21, 2022
ETNA is an easy-to-use time series forecasting framework.

ETNA is an easy-to-use time series forecasting framework. It includes built in toolkits for time series preprocessing, feature generation, a variety of predictive models with unified interface - from

Tinkoff.AI 674 Jan 07, 2023
A GitHub action that suggests type annotations for Python using machine learning.

Typilus: Suggest Python Type Annotations A GitHub action that suggests type annotations for Python using machine learning. This action makes suggestio

40 Sep 18, 2022
Sequence learning toolkit for Python

seqlearn seqlearn is a sequence classification toolkit for Python. It is designed to extend scikit-learn and offer as similar as possible an API. Comp

Lars 653 Dec 27, 2022
fastFM: A Library for Factorization Machines

Citing fastFM The library fastFM is an academic project. The time and resources spent developing fastFM are therefore justified by the number of citat

1k Dec 24, 2022
pandas, scikit-learn, xgboost and seaborn integration

pandas, scikit-learn and xgboost integration.

299 Dec 30, 2022
This project impelemented for midterm of the Machine Learning #Zoomcamp #Alexey Grigorev

MLProject_01 This project impelemented for midterm of the Machine Learning #Zoomcamp #Alexey Grigorev Context Dataset English question data set file F

Hadi Nakhi 1 Dec 18, 2021
Python implementation of Weng-Lin Bayesian ranking, a better, license-free alternative to TrueSkill

Python implementation of Weng-Lin Bayesian ranking, a better, license-free alternative to TrueSkill This is a port of the amazing openskill.js package

Open Debates Project 156 Dec 14, 2022
Machine Learning Algorithms

Machine-Learning-Algorithms In this project, the dataset was created through a survey opened on Google forms. The purpose of the form is to find the p

Göktuğ Ayar 3 Aug 10, 2022
PennyLane is a cross-platform Python library for differentiable programming of quantum computers

PennyLane is a cross-platform Python library for differentiable programming of quantum computers. Train a quantum computer the same way as a neural ne

PennyLaneAI 1.6k Jan 01, 2023
A statistical library designed to fill the void in Python's time series analysis capabilities, including the equivalent of R's auto.arima function.

pmdarima Pmdarima (originally pyramid-arima, for the anagram of 'py' + 'arima') is a statistical library designed to fill the void in Python's time se

alkaline-ml 1.3k Dec 22, 2022
Transform ML models into a native code with zero dependencies

m2cgen (Model 2 Code Generator) - is a lightweight library which provides an easy way to transpile trained statistical models into a native code

Bayes' Witnesses 2.3k Jan 03, 2023
Adaptive: parallel active learning of mathematical functions

adaptive Adaptive: parallel active learning of mathematical functions. adaptive is an open-source Python library designed to make adaptive parallel fu

741 Dec 27, 2022