FreeFem++ version (2d and 3d)
Université Pierre et Marie Curie Laboratoire Jacques-Louis Lions
Financé par l'ANR / FF2A3, Réf. : ANR-07-CIS7-002-01

FreeFem++ team : Olivier Pironneau, Frédéric Hecht, Antoine Le Hyaric, Jacques Morice ,

Please send questions, comments and bug reports to Frédéric Hecht or to the FreeFem++ mailing list.

 
Introduction

 
FreeFem++ is an implementation of a language dedicated to the finite element method. It enables you to solve Partial Differential Equations (PDE) easily.

Problems involving PDE (2d, 3d) from several branches of physics such as fluid-structure interactions require interpolations of data on several meshes and their manipulation within one program. FreeFem++ includes a fast 2^d-tree-based interpolation algorithm and a language for the manipulation of data on multiple meshes (as a follow up of bamg).

FreeFem++ is written in C++ and the FreeFem++ language is a C++ idiom. It runs on any Unix-like OS (with g++ version 3 or higher, X11R6 or OpenGL with GLUT) Linux, FreeBSD, Solaris 10, Microsoft Windows ( 2000, NT, XP, Vista,7 ) and MacOS X (native version using OpenGL). FreeFem++ replaces the older freefem and freefem+.

 
Quick Start

 

The FreeFem++-cs (IDE: an integrated development environment) by Antoine Le Hyaric is a excellent tool to develop FreeFem++ scripts.

 

Unix
Linux, FreeBSD, ...
Windows
98, NT, 2000, XP
MacOS X 10.6, 10.7
Intel
Documentation

 

Download this sources archive,
decompress, compile;
Download this EXE
and run it
Download
pkg
install the package
Full documentation :
freefem++doc.pdf

 

, install and run to install mpi with mpi

 

  • To compile FreeFem++, go into the FreeFem++ directory and type :
    ./configure --enable-download
    make
    make install (as root).

  • For all other platforms and versions, please refer to the full list of downloads.

     

     
    News

    1/10/2005 Read or participate in general discussions about FreeFEM++ with the FreeFEM++ mailing list

    23/9/2005 FreeFem++ and OpenFEM presentation day at IHP (Paris) : details, slides and examples

    16/10/2008 Ateliers de Simualtion Numérique, (Setif,Algérie)  : details, par Olivier Pantz

    14-15/09/2009 The first tutorial and Workshop on FreeFem++ , held September 14th and 15th, 2009, in Paris at IHP, amphi Hermite ,11 Rue Pierre et Marie Curie. All Presentation, examples, and data are here

    1-2/09/2010 The second tutorial and Workshop on FreeFem++ , held September 1st and 2nd, 2010, in Paris at Université Pierre et Marie Curie, Barre 16-15, 3ieme, 4 place Jussieu, Paris All Presentation, examples, and data are here

    5-7/12/2011 The THIRD tutorial and Workshop on FreeFem++ , held december 5th,6th and 7th, 2010, in Paris at Université Pierre et Marie Curie, Barre 16-15, 3ieme, 4 place Jussieu, Paris All Presentation, examples, and data are here

    Some Freefem++ presentation (with useful information):

    Examples in 2d and 3d.
    • A small movie (340Kb) : Cool air (green) comes from the lower left and mix with hot air (magenta), the right boundary is free. This is Navier-Stokes-Boussinesq integrated with P1-bubblle P1 mixte finite element.
    • A very small example 2d of how to solve the Poisson equation on a L shape :
      border aaa(t=0,1){x=t;y=0;};
      border bbb(t=0,0.5){x=1;y=t;};
      border ccc(t=0,0.5){x=1-t;y=0.5;};
      border ddd(t=0.5,1){x=0.5;y=t;};
      border eee(t=0.5,1){x=1-t;y=1;};
      border fff(t=0,1){x=0;y=1-t;};
      mesh Th = buildmesh (aaa(6) + bbb(4) + ccc(4) +ddd(4) + eee(4) + fff(6));
      fespace Vh(Th,P1);   //  to change P1 in P2 to make P2 finite element.
      Vh u=0,v;
      func f= 1;
      func g= 0;
      int i=0;
      real error=0.1, coef= 0.1^(1./5.);
      problem Probem1(u,v,solver=CG,eps=-1.0e-6) =
          int2d(Th)(  dx(u)*dx(v) + dy(u)*dy(v)) 
        + int2d(Th) ( v*f ) 
        + on(aaa,bbb,ccc,ddd,eee,fff,u=g)  ;
        
      for (i=0;i< 10;i++)
      {   
        real d = clock();
        Probem1; //  solve the problem 
        plot(u,Th,wait=1);
        Th=adaptmesh(Th,u,inquire=1,err=error);
        error = error * coef;
      } ;
      
      

      Solution on adapted mesh and associated mesh.

    A 3d example
    • A very small example of how to solve the Stokes equation 3d on cube shape :
      load "msh3" load "medit"  // dynamics load tools for 3d.
      int nn=8;
      mesh Th2=square(nn,nn);
      fespace Vh2(Th2,P2);  Vh2 ux,uz,p2;
      int[int] rup=[0,2],  rdown=[0,1], rmid=[1,1,2,1,3,1,4,1];
      real zmin=0,zmax=1;
      mesh3 Th=buildlayers(Th2,nn,
        zbound=[zmin,zmax],  reffacemid=rmid, 
        reffaceup = rup,     reffacelow = rdown);
        
      medit("c10x10x10",Th);  // see the 3d mesh with medit software
      Fe macro Grad(u) [dx(u),dy(u),dz(u)] // EOM
       macro div(u1,u2,u3) (dx(u1)+dy(u2)+dz(u3))  //EOM
      
      VVh [u1,u2,u3,p];
      VVh [v1,v2,v3,q];
        
      solve vStokes([u1,u2,u3,p],[v1,v2,v3,q]) = 
        int3d(Th,qforder=3)( Grad(u1)'*Grad(v1) +  Grad(u2)'*Grad(v2) +  Grad(u3)'*Grad(v3)
                        - div(u1,u2,u3)*q - div(v1,v2,v3)*p + 1e-10*q*p ) 
        + on(2,u1=1.,u2=0,u3=0) + on(1,u1=0,u2=0,u3=0) ;
       plot(p,wait=1, nbiso=5);  // a 3d plot of iso  pressure. in progress... march 2009
       //  to see the 10 cut plan in 2d 
      for(int i=1;i<10;i++)
      {
       real yy=i/10.; // compute yy.
        // do 3d -> 2d interpolation.
       ux= u1(x,yy,y); uz= u3(x,yy,y);  p2= p(x,yy,y);
       plot([ux,uz],p2,cmm=" cut y = "+yy,wait= 1);
      }
      

      Solution on cup plan y=0.5 and mesh 10x10x10 and associated mesh.

    Ongoing Work
    • FreeVol: Finite Vol thecnics in freefem for hyperbolic PDEs
    • 3D implementation: new solver, new mesh tools, new kind of finite element
    • Stabilize and test all // linear solver interface

     
    Full List of Downloads

    The current version of Freefem++ is . You can get the latest source from an anonymous Mercurial SCM copy with the following unix shell commands :

    hg clone  http://www.freefem.org/ff++/ff++
    

    Remark in the new Mac tested on MacOS 10.6 and 10.7 version (download), the MPI tools are include except the openmpi library, so to download see this page .

    To recompile freefem++ with mpi interface freefem++ mpi under Windows See this page .

    Self-contained archives for all other systems :

  • Setup file for generic hardware and Windows 95,98,NT,2000, XP FreeFem++-3.18.exe 17861.5 Kb Sun Jan 15 16:10:51 2012
  • All the sources (tar+gzip) freefem++-3.18.tar.gz 23702.2 Kb Sun Jan 15 14:42:57 2012
  • -MacOS 10.[67] MacOS 10.6 or better, Cocoa, MPI, OpenGL (pkg) FreeFem++-3.17-MacOS_10.[67].pkg 41308.4 Kb Fri Dec 2 14:50:53 2011
  • -MacOS 10.6 MacOS 10.6 or better, Cocoa, MPI, OpenGL (pkg) FreeFem++-3.18-MacOS_10.6.pkg 41300.4 Kb Mon Jan 16 15:33:57 2012
  • -snow-leopard MacOsX (tar+ gzip file) FreeFem++v3.13-3-snow-leopard_MacOsX.tgz 16090.8 Kb Thu Jun 30 12:08:29 2011
  • -Universal MacOsX (tar+ gzip file) FreeFem++v3.13-Universal_MacOsX.tgz 20848.2 Kb Fri Jun 10 09:40:16 2011
  • -lion MacOsX (tar+ gzip file) FreeFem++v3.14-2-lion_MacOsX.tgz 32631.4 Kb Thu Sep 29 15:40:55 2011
  • INNOVATION INNOVATION 40.7 Kb Mon Jan 16 11:08:04 2012
  • HISTORY HISTORY 330.5 Kb Thu Jan 20 13:34:07 2011
  • HISTORY BEFORE 2005 HISTORY_BEFORE_2005 22 Kb Wed Dec 1 15:47:23 2010
  • New manual (pdf file) freefem++doc.pdf 48479.6 Kb Mon Jan 16 13:16:52 2012
  • This directory contains all the different versions of FreeFem.

     
    Page written by

    Frédéric Hecht et Antoine Le Hyaric

    Last modified : Mon Jan 16 15:50:02 2012