C++ C Makefile Fortran Shell CMake Other
Switch branches/tags
Permalink
Failed to load latest commit information.
.buildkite Update pipeline.yml May 12, 2017
.github github template Jan 28, 2017
.templates various changes: Mar 16, 2015
applications fixed script install #896 May 9, 2017
benchmarks add bench_perf_io_vector Feb 3, 2017
cmake move back to 0.104.0 [ci skip] May 21, 2017
contrib Merge pull request #897 from feelpp/feature/nlopt May 9, 2017
data Merge branch 'develop' of https://github.com/feelpp/feelpp into develop Feb 18, 2017
databases/materials dynamic viscosity is named mu Apr 6, 2015
doc use nlopt algo LD_MMA May 9, 2017
feel fix inconsistent override compilation warning May 17, 2017
ports Minor changes to cesga port in order to compile PETSc with MKL May 3, 2017
projects rm submodules from projects Jan 15, 2017
quickstart Merge remote-tracking branch 'origin/develop' into feature/cmake-depe… May 5, 2017
research rm cerema Jan 15, 2017
testsuite Merge remote-tracking branch 'origin/develop' into feature/nlopt May 9, 2017
toolboxes install headers of toolboxes/feel/modelvf/ May 16, 2017
tools Merge branch 'develop' May 21, 2017
.clang-format Feature/visibility (#783) Dec 28, 2016
.dir-locals.el various changes: Mar 16, 2015
.emacs-dirvars various changes: Mar 16, 2015
.gitattributes up git attributes + install example [ci skip] Jan 24, 2017
.gitignore Add Cotire Pch + cling interpreter script Mar 6, 2017
.gitmodules rm submodules from projects Jan 15, 2017
.travis.yml add ubuntu 17.04 to travis May 6, 2017
CMakeLists.txt fix feel++ installation with make install May 18, 2017
CMakeLists.txt.user rm Feb 6, 2014
CONTRIBUTING.adoc minor update [ci skip] Apr 17, 2017
COPYING.adoc update COPYING according to the LICENSE update [ci skip] Apr 17, 2017
CTestConfig.cmake fixes #372 and #374 Jun 9, 2014
INSTALL.md up INSTALL [ci-skip] Jun 29, 2014
LICENSE.adoc use adoc file [ci skip] Apr 15, 2017
README.adoc add news and docs [ci skip] Apr 17, 2017
configure various changes related to #862 and #863 Apr 8, 2017
feelpp.sublime-project add sublime text project file Aug 5, 2014

README.adoc

Feel++: Finite Element Embedded Library in C++

Feel++ is a C++ library for arbitrary order Galerkin methods (e.g. finite element method and spectral element methods) continuous or discontinuous in 1D 2D and 3D.

Feel++ Documentation

Gitter Discussion Forum

We encourage you to ask questions and discuss any aspects of the project on the Feel++ Gitter forum. New contributors are always welcome!

Join the chat at

Continuous Integration

Feel++ maintains various branches. At the core, the development model is greatly inspired by existing models out there. The central repo holds two main branches with an infinite lifetime: master and develop

master

Main branch where the source code of HEAD always reflects a production-ready state.

develop

Main branch where the source code of HEAD always reflects a state with the latest delivered development changes for the next release. Some would call this the “integration branch”. This is where any automatic nightly builds are built from.

feature/*

Feature branches (or sometimes called topic branches) are used to develop new features for the upcoming or a distant future release. When starting development of a feature, the target release in which this feature will be incorporated may well be unknown at that point. The essence of a feature branch is that it exists as long as the feature is in development, but will eventually be merged back into develop (to definitely add the new feature to the upcoming release) or discarded (in case of a disappointing experiment).

Platform & Compiler master develop feature/hdg-sc feature/quad

Travis/Ubuntu & Clang >= 3.5

Build Status

Build Status

Build Status

Build Status

Buildkite Ubuntu 16.10

Build Status

Build Status

Build Status

Build Status

What is Feel++?

Feel++ is a C++ library for arbitrary order Galerkin methods (e.g. finite and spectral element methods) continuous or discontinuous in 1D 2D and 3D. The objectives of this framework is quite ambitious; ambitions which could be express in various ways such as :

  • the creation of a versatile mathematical kernel solving easily problems using different techniques thus allowing testing and comparing methods, e.g. cG versus dG,

  • the creation of a small and manageable library which shall nevertheless encompass a wide range of numerical methods and techniques,

  • build mathematical software that follows closely the mathematical abstractions associated with partial differential equations (PDE),

  • the creation of a library entirely in C++ allowing to create complex and typically multi-physics applications such as fluid-structure interaction or mass transport in haemodynamic.

Some basic installation procedure are available in the INSTALL file, the detailled process is available here.

Releases

Here are the latest releases of Feel++

Features

  • 1D 2D and 3D (including high order) geometries and also lower topological dimension 1D(curve) in 2D and 3D or 2D(surface) in 3D

  • continuous and discontinuous arbitrary order Galerkin Methods in 1D, 2D and 3D including finite and spectral element methods

  • domain specific embedded language in C++ for variational formulations

  • interfaced with PETSc for linear and non-linear solvers

  • seamless parallel computations using PETSc

  • interfaced with SLEPc for large-scale sparse standard and generalized eigenvalue solvers

  • supports Gmsh for mesh generation

  • supports Gmsh for post-processing (including on high order geometries)

  • supports Paraview and CEI/Ensight for post-processing and the following file formats: ensight gold, gmsh, xdmf.

Contributing

In the spirit of free software, everyone is encouraged to help improve this project. If you discover errors or omissions in the source code, documentation, or website content, please don’t hesitate to submit an issue or open a pull request with a fix. New contributors are always welcome!

Here are some ways you can contribute:

  • by using develop versions

  • by reporting bugs

  • by suggesting new features

  • by writing or editing documentation

  • by writing specifications

  • by writing code — No patch is too small.

    • fix typos

    • add comments

    • write examples!

    • write tests!

  • by refactoring code

  • by fixing issues

  • by reviewing Pull Requests

The Contributing guide provides information on how to create, style, and submit issues, feature requests, code, and documentation to the Feel++ Project.

Getting Help

The Feel++ project is developed to help you easily do (i) modelisation simulation and optimisation and (ii) high performance computing. But we can’t do it without your feedback! We encourage you to ask questions and discuss any aspects of the project on the discussion list, on Twitter or in the chat room.

Twitter

#feelpp hashtag or @feelpp mention

Chat (Gitter)

Gitter

Further information and documentation about Feel++ can be found on the project’s website.

Home | News | Docs

The Feel++ organization on GitHub hosts the project’s source code, issue tracker, and sub-projects.

Source repository (git)

https://github.com/feelpp/feelpp

Issue tracker

https://github.com/feelpp/feelpp/issues

Feel++ organization on GitHub

https://github.com/feelpp

Copyright © 2011-2017 Feel++ Consortium. Free use of this software is granted under the terms of the GPL License.

See the LICENSE file for details.

Authors

Feel++ is led by Christophe Prud’homme and has received contributions from many other individuals. The project was initiated in 2006 by Christophe Prud’homme and based initially on lifeV and completely re-written since then.

Examples

Laplacian in 2D using P3 Lagrange basis functions

Here is a full example to solve -\Delta u = f \mbox{ in } \Omega,\quad u=g \mbox{ on } \partial \Omega

#include <feel/feel.hpp>

int main(int argc, char**argv )
{
    using namespace Feel;
    Environment env( _argc=argc, _argv=argv,
                     _desc=feel_options(),
                     _about=about(_name="qs_laplacian",
                                  _author="Feel++ Consortium",
                                  _email="feelpp-devel@feelpp.org"));

    auto mesh = unitSquare();
    auto Vh = Pch<1>( mesh );
    auto u = Vh->element();
    auto v = Vh->element();

    auto l = form1( _test=Vh );
    l = integrate(_range=elements(mesh),
                  _expr=id(v));

    auto a = form2( _trial=Vh, _test=Vh );
    a = integrate(_range=elements(mesh),
                  _expr=gradt(u)*trans(grad(v)) );
    a+=on(_range=boundaryfaces(mesh), _rhs=l, _element=u,
          _expr=constant(0.) );
    a.solve(_rhs=l,_solution=u);

    auto e = exporter( _mesh=mesh, _name="qs_laplacian" );
    e->add( "u", u );
    e->save();
    return 0;
}

Bratu equation in 2D

Here is a full non-linear example - the Bratu equation - to solve -\Delta u + e^u = 0 \mbox{ in } \Omega,\quad u=0 \mbox{ on } \partial \Omega.

#include <feel/feel.hpp>

inline
Feel::po::options_description
makeOptions()
{
    Feel::po::options_description bratuoptions( "Bratu problem options" );
    bratuoptions.add_options()
    ( "lambda", Feel::po::value<double>()->default_value( 1 ),
                "exp() coefficient value for the Bratu problem" )
    ( "penalbc", Feel::po::value<double>()->default_value( 30 ),
                 "penalisation parameter for the weak boundary conditions" )
    ( "hsize", Feel::po::value<double>()->default_value( 0.1 ),
               "first h value to start convergence" )
    ( "export-matlab", "export matrix and vectors in matlab" )
    ;
    return bratuoptions.add( Feel::feel_options() );
}

/**
 * Bratu Problem
 *
 * solve \f$ -\Delta u + \lambda \exp(u) = 0, \quad u_\Gamma = 0\f$ on \f$\Omega\f$
 */
int
main( int argc, char** argv )
{

    using namespace Feel;
    Environment env( _argc=argc, _argv=argv,
                     _desc=makeOptions(),
                     _about=about(_name="bratu",
                                  _author="Christophe Prud'homme",
                                  _email="christophe.prudhomme@feelpp.org"));
    auto mesh = unitSquare();
    auto Vh = Pch<3>( mesh );
    auto u = Vh->element();
    auto v = Vh->element();
    double penalbc = option(_name="penalbc").as<double>();
    double lambda = option(_name="lambda").as<double>();

    auto Jacobian = [=](const vector_ptrtype& X, sparse_matrix_ptrtype& J)
        {
            auto a = form2( _test=Vh, _trial=Vh, _matrix=J );
            a = integrate( elements( mesh ), gradt( u )*trans( grad( v ) ) );
            a += integrate( elements( mesh ), lambda*( exp( idv( u ) ) )*idt( u )*id( v ) );
            a += integrate( boundaryfaces( mesh ),
               ( - trans( id( v ) )*( gradt( u )*N() ) - trans( idt( u ) )*( grad( v )*N()  + penalbc*trans( idt( u ) )*id( v )/hFace() ) );
        };
    auto Residual = [=](const vector_ptrtype& X, vector_ptrtype& R)
        {
            auto u = Vh->element();
            u = *X;
            auto r = form1( _test=Vh, _vector=R );
            r = integrate( elements( mesh ), gradv( u )*trans( grad( v ) ) );
            r +=  integrate( elements( mesh ),  lambda*exp( idv( u ) )*id( v ) );
            r +=  integrate( boundaryfaces( mesh ),
               ( - trans( id( v ) )*( gradv( u )*N() ) - trans( idv( u ) )*( grad( v )*N() ) + penalbc*trans( idv( u ) )*id( v )/hFace() ) );
        };
    u.zero();
    backend()->nlSolver()->residual = Residual;
    backend()->nlSolver()->jacobian = Jacobian;
    backend()->nlSolve( _solution=u );

    auto e = exporter( _mesh=mesh );
    e->add( "u", u );
    e->save();
}