Appendix C
Dynamical link

Now, it’s possible to add built-in functionnalites in FreeFem++ under the three environnents Linux, Windows and MacOS X 10.3 or newer. It is agood idea to, first try the example load.edp in directory example++-load.

But unfortunately, you need to install a c++ compiler (generally g++/gcc compiler) to compile your function.

Windows
Install the cygwin environnent or the mingw
MacOs
Install the developer tools xcode on the apple DVD
Linux/Unix
Install the correct compiler (gcc for instance)

Now, assume , you are in a shell window (a cygwin window under Windows) in the directory example++-load. Remark that in the sub directory include they are all the FreeFem++ include file to make the link with FreeFem++.

Note C.1 If you try to load dynamically a file with command load "xxx"

The compilation of your module: the script ff-c++ compiles and makes the link with FreeFem++, but be careful, the script has no way to known if you try to compile for a pure Windows environment or for a cygwin environment so to build the load module under cygwin you must add the -cygwin parameter.

C.1 A first example myfunction.cpp

The following defines a new function call myfunction with no parameter, but using the x,y current value.


#include  <iostream>
#include  <cfloat>
using namespace std;
#include "error.hpp"
#include "AFunction.hpp"
#include "rgraph.hpp"
#include "RNM.hpp"
#include "fem.hpp"
#include "FESpace.hpp"
#include "MeshPoint.hpp"
using namespace Fem2D;
double myfunction(Stack stack)
{
   // to get FreeFem++ data
  MeshPoint &mp= ⋆MeshPointStack(stack); // the struct to get x,y, normal , value
  double x= mp.P.x;  // get the current x value
  double y= mp.P.y;  // get the current y value
   // cout << "x = " << x << " y=" << y << endl;
  return sin(x)⋆cos(y);
}

Now the Problem is to build the link with FreeFem++, to do that we need two classes, one to call the function myfunction

All FreeFem++ evaluable expression must be a struct/class C++ which derive from E_F0. By default this expression does not depend of the mesh position, but if they derive from E_F0mps the expression depends of the mesh position, and for more details see [12].


// A class build the link with FreeFem++
// generaly this class are already in AFunction.hpp
// but unfortunatly, I have no simple function with no parameter
// in FreeFem++ depending of the mesh,
template<class R>
class  OneOperator0s : public OneOperator {
   // the class to defined a evaluated a new function
   // It must devive from E_F0 if it is mesh independent
   // or from E_F0mps if it is mesh dependent
  class E_F0_F :public  E_F0mps { public:
    typedef  R (⋆func)(Stack stack) ;
    func f;  // the pointeur to the fnction myfunction
    E_F0_F(func ff)  : f(ff) {}
     // the operator evaluation in FreeFem++
    AnyType operator()(Stack stack)  const {return SetAny<R>( f(stack)) ;}
  };
  typedef  R (⋆func)(Stack ) ;
  func  f;
public:
   // the function which build the FreeFem++ byte code
  E_F0 ⋆ code(const basicAC_F0 & ) const { return  new E_F0_F(f);}
   // the constructor to say ff is a function without parameter
   // and returning a R
  OneOperator0s(func  ff): OneOperator(map_type[typeid(R).name()]),f(ff){}
};

To finish we must add this new function in FreeFem++ table , to do that include :


class Init { public:  Init();  };
Init init;
Init::Init(){
  Global.Add("myfunction","(",new OneOperator0s<double>(myfunction));
}

It will be called automatically at load module time.

To compile and link, do


Brochet% ff-c++ myfunction.cpp
g++ -c -g -Iinclude myfunction.cpp
g++ -bundle -undefined dynamic_lookup -g myfunction.o -o ./myfunction.dylib

To, try the simple example under Linux or MacOS, do


Brochet% FreeFem++-nw load.edp
-- FreeFem++ v 1.4800028 (date Tue Oct  4 11:56:46 CEST 2005)
 file : load.edp
 Load: lg_fem lg_mesh eigenvalue  UMFPACK
    1 :  // Example of dynamic function load
    2 :  // --------------------------------
    3 :  // Id : freefem + +doc.tex,v1.862009051307 : 12 : 04hechtExp
    4 :
    5 :  load "myfunction"
     lood: myfunction
load: dlopen(./myfunction) = 0xb01cc0
    6 :  mesh Th=square(5,5);
    7 :  fespace Vh(Th,P1);
    8 :  Vh uh=myfunction();  // warning do not forget ()
    9 :  cout << uh[].min << " " << uh[].max << endl;
   10 :  sizestack + 1024 =1240  ( 216 )
 -- square mesh : nb vertices  =36 ,  nb triangles = 50 ,  nb boundary edges 20
   Nb of edges on Mortars  = 0
   Nb of edges on Boundary = 20, neb = 20
 Nb Of Nodes = 36
 Nb of DF = 36
0 0.841471
times: compile 0.05s, execution -3.46945e-18s
 CodeAlloc : nb ptr  1394,  size :71524
Bien: On a fini Normalement

Under Windows, launch FreeFem++ with the mouse on the example.

C.2 Example Discrete Fast Fourier Transform

This will add FFT to FreeFem++,taken from http://www.fftw.org/. To download and install under download/include just go in download/fftw and trymake.

The 1D dfft (fast discret fourier transform) for a simple array f of size n is defined by the following formula

             ∑n-1
df ft(f,ε)k =    fieε2πikj∕n
              j=0

The 2D DFT for an array of size N = n ×m is

                  m∑-1 n∑-1             ′
df ft(f,m, ε)k+nl =       fi+njeε2πi(kj∕n+lj∕m)
                  j′=0 j=0
Remark: the value n is given by size(f)∕m, and the numbering is row-major order.

So the classical discrete DFT is ^f = dfft(f,-1)√ --
  n and the reverse dFT f = dfft(^f ,1)√ --
  n

Remark: the 2D Laplace operator is

           √ --m∑-1 n∑-1            ′
f (x,y) = 1∕  N        f^i+njeε2πi(xj+yj)
               j′=0 j=0
and we have
fk+nl = f(k∕n,l∕m)

So

             2   2    2
^Δfkl = -((2π)((~k) + (~l)))^fkl
where k~ = k if k n∕2 else ~k = k -n and l~= l if l m∕2 else ~l= l -m.

And to get a real function we need all symetric modes around to zero, so n and m must be odd.

To compile and make a new library

%  ff-c++  dfft.cpp ../download/install/lib/libfftw3.a -I../download/install/include  
export MACOSX_DEPLOYMENT_TARGET=10.3  
g++ -c -Iinclude -I../download/install/include dfft.cpp  
g++ -bundle -undefined dynamic_lookup dfft.o -o ./dfft.dylib ../download/install/lib/libfftw3.a

To test ,


-- FreeFem++ v 1.4800028 (date Mon Oct 10 16:53:28 EEST 2005)
 file : dfft.edp
 Load: lg_fem cadna lg_mesh eigenvalue  UMFPACK
    1 :  // Example of dynamic function load
    2 :  // --------------------------------
    3 :  // Id : freefem + +doc.tex,v1.862009051307 : 12 : 04hechtExp
    4 :  // Discret Fast Fourier Transform
    5 :  // -------------------------------
    6 :  load "dfft" lood: init dfft
load: dlopen(dfft.dylib) = 0x2b0c700
    7 :
    8 : int nx=32,ny=16,N=nx⋆ny;
    9 :  // warning the fourier space is not exactly the unite square
                       // due to periodic condition
   10 : mesh Th=square(nx-1,ny-1,[(nx-1)⋆x/nx,(ny-1)⋆y/ny]);
   11 :  // warring the numbering is of the vertices (x,y) is
   12 :  // given by i = x∕nx + nx * y∕ny
   13 :
   14 : fespace Vh(Th,P1);
   15 :
   16 : func f1 = cos(2⋆x⋆2⋆pi)⋆cos(3⋆y⋆2⋆pi);
   17 : Vh<complex> u=f1,v;
   18 : Vh w=f1;
   19 :
   20 :
   21 : Vh  ur,ui;
   22 :  // in dfft the matrix n,m is in row-major order ann array n,m is
   23 :  // store j + m⋆ i ( the transpose of the square numbering )
   24 :  v[]=dfft(u[],ny,-1);
   25 :  u[]=dfft(v[],ny,+1);
   26 :  u[] /= complex(N);
   27 :  v = f1-u;
   28 : cout << " diff = "<< v[].max << " " <<  v[].min << endl;
   29 : assert( norm(v[].max) < 1e-10 &&  norm(v[].min) < 1e-10) ;
   30 :   // ------- a more hard example ----
   31 :   // Lapacien en FFT
   32 :   // -Δu = f with biperiodic condition
   33 : func f = cos(3⋆2⋆pi⋆x)⋆cos(2⋆2⋆pi⋆y);  //
   34 : func ue =  +(1./(square(2⋆pi)⋆13.))⋆cos(3⋆2⋆pi⋆x)⋆cos(2⋆2⋆pi⋆y);   //
   35 : Vh<complex> ff = f;
   36 : Vh<complex> fhat;
   37 : fhat[] = dfft(ff[],ny,-1);
   38 :
   39 : Vh<complex> wij;
   40 :  // warning in fact we take mode between -nx/2, nx/2 and -ny/2,ny/2
   41 :  // thank to the operator ?:
   42 : wij = square(2.⋆pi)⋆(square(( x<0.5?x⋆nx:(x-1)⋆nx))
            + square((y<0.5?y⋆ny:(y-1)⋆ny)));
   43 : wij[][0] = 1e-5;  // to remove div / 0
   44 : fhat[] = fhat[]./ wij[];   //
   45 : u[]=dfft(fhat[],ny,1);
   46 : u[] /= complex(N);
   47 : ur = real(u);  // the solution
   48 : w = real(ue);  // the exact solution
   49 : plot(w,ur,value=1 ,cmm=" ue   ", wait=1);
   50 : w[] -= ur[];  // array sub
   51 : real err= abs(w[].max)+abs(w[].min) ;
   52 : cout << " err = " << err << endl;
   53 : assert( err  < 1e-6);
   54 :  sizestack + 1024 =3544  ( 2520 )
----------CheckPtr:-----init execution ------ NbUndelPtr  2815  Alloc: 111320  NbPtr 6368
 -- square mesh : nb vertices  =512 ,  nb triangles = 930 ,  nb boundary edges 92
   Nb of edges on Mortars  = 0
   Nb of edges on Boundary = 92, neb = 92
 Nb Of Nodes = 512
 Nb of DF = 512
0x2d383d8 -1 16 512 n: 16 m:32
 dfft 0x402bc08 = 0x4028208 n = 16 32 sign = -1
 --- --- ---0x2d3ae08 1 16 512 n: 16 m:32
 dfft 0x4028208 = 0x402bc08 n = 16 32 sign = 1
 --- --- --- diff = (8.88178e-16,3.5651e-16) (-6.66134e-16,-3.38216e-16)
0x2d3cfb8 -1 16 512 n: 16 m:32
 dfft 0x402de08 = 0x402bc08 n = 16 32 sign = -1
 --- --- ---0x2d37ff8 1 16 512 n: 16 m:32
 dfft 0x4028208 = 0x402de08 n = 16 32 sign = 1
 --- --- --- err = 3.6104e-12
times: compile 0.13s, execution 2.05s
----------CheckPtr:-----end execution -- ------ NbUndelPtr  2815  Alloc: 111320  NbPtr 26950
 CodeAlloc : nb ptr  1693,  size :76084
Bien: On a fini Normalement
                CheckPtr:Nb of undelete pointer is 2748 last 114
                CheckPtr:Max Memory used    228.531 kbytes  Memory undelete 105020

C.3 Load Module for Dervieux’ P0-P1 Finite Volume Method

the associed edp file is examples++-load/convect_dervieux.edp


// Implementation of P1-P0 FVM-FEM
// ---------------------------------------------------------------------
// Id : freefem + +doc.tex,v1.862009051307 : 12 : 04hechtExp
// compile and link with ff-c++ mat_dervieux.cpp (i.e. the file name without .cpp)
#include  <iostream>
#include  <cfloat>
#include  <cmath>
using namespace std;
#include "error.hpp"
#include "AFunction.hpp"
#include "rgraph.hpp"
#include "RNM.hpp"
// remove problem of include
#undef  HAVE_LIBUMFPACK
#undef HAVE_CADNA
#include "MatriceCreuse_tpl.hpp"
#include "MeshPoint.hpp"
#include "lgfem.hpp"
#include "lgsolver.hpp"
#include "problem.hpp"
class MatrixUpWind0 :  public E_F0mps { public:
  typedef Matrice_Creuse<R> ⋆ Result;
  Expression emat,expTh,expc,expu1,expu2;
  MatrixUpWind0(const basicAC_F0 & args)
  {
    args.SetNameParam();
    emat =args[0];  // the matrix expression
    expTh= to<pmesh>(args[1]);   // a the expression to get the mesh
    expc = CastTo<double>(args[2]);  // the expression to get c (must be a double)
     // a array expression [ a, b]
    const E_Array ⋆ a= dynamic_cast<const E_Array⋆>((Expression) args[3]);
    if (a->size() != 2) CompileError("syntax:  MatrixUpWind0(Th,rhi,[u1,u2])");
    int err =0;
    expu1= CastTo<double>((⋆a)[0]);  // fist exp of the array (must be a double)
    expu2= CastTo<double>((⋆a)[1]);  // second exp of the array (must be a double)
  }
  ~MatrixUpWind0()
  {
  }
  static ArrayOfaType  typeargs()
  { return  ArrayOfaType(atype<Matrice_Creuse<R>⋆>(),
    atype<pmesh>(),atype<double>(),atype<E_Array>());}
  static  E_F0 ⋆ f(const basicAC_F0 & args){ return new MatrixUpWind0(args);}
  AnyType operator()(Stack s) const ;
};
int   fvmP1P0(double q[3][2], double u[2],double c[3], double a[3][3], double where[3] )
{                                // computes matrix a on a triangle for the Dervieux FVM
  for(int i=0;i<3;i++) for(int j=0;j<3;j++) a[i][j]=0;
    for(int i=0;i<3;i++){
        int ip = (i+1)%3, ipp =(ip+1)%3;
        double unL =-((q[ip][1]+q[i][1]-2⋆q[ipp][1])⋆u[0]
                -(q[ip][0]+q[i][0]-2⋆q[ipp][0])⋆u[1])/6;
        if(unL>0) { a[i][i] += unL; a[ip][i]-=unL;}
            else{ a[i][ip] += unL; a[ip][ip]-=unL;}
        if(where[i]&&where[ip]){         // this is a boundary edge
            unL=((q[ip][1]-q[i][1])⋆u[0] -(q[ip][0]-q[i][0])⋆u[1])/2;
            if(unL>0) { a[i][i]+=unL; a[ip][ip]+=unL;}
        }
    }
  return 1;
}
// the evaluation routine
AnyType MatrixUpWind0::operator()(Stack stack) const
{
  Matrice_Creuse<R> ⋆ sparse_mat =GetAny<Matrice_Creuse<R>⋆ >((⋆emat)(stack));
  MatriceMorse<R> ⋆ amorse =0;
  MeshPoint ⋆mp(MeshPointStack(stack)) , mps=⋆mp;
  Mesh ⋆ pTh = GetAny<pmesh>((⋆expTh)(stack));
  ffassert(pTh);
  Mesh & Th (⋆pTh);
  {
    map< pair<int,int>, R> Aij;
    KN<double> cc(Th.nv);
    double infini=DBL_MAX;
    cc=infini;
    for (int it=0;it<Th.nt;it++)
      for (int iv=0;iv<3;iv++)
    {
      int i=Th(it,iv);
      if ( cc[i]==infini) {  // if nuset the set
        mp->setP(&Th,it,iv);
        cc[i]=GetAny<double>((⋆expc)(stack));
      }
    }
    for (int k=0;k<Th.nt;k++)
      {
    const Triangle & K(Th[k]);
    const Vertex & A(K[0]), &B(K[1]),&C(K[2]);
    R2 Pt(1./3.,1./3.);
    R u[2];
    MeshPointStack(stack)->set(Th,K(Pt),Pt,K,K.lab);
    u[0] = GetAny< R>( (⋆expu1)(stack) ) ;
    u[1] = GetAny< R>( (⋆expu2)(stack) ) ;
    int ii[3] ={  Th(A), Th(B),Th(C)};
    double q[3][2]= { { A.x,A.y} ,{B.x,B.y},{C.x,C.y} } ;  // coordinates of 3 vertices (input)
    double c[3]={cc[ii[0]],cc[ii[1]],cc[ii[2]]};
    double a[3][3], where[3]={A.lab,B.lab,C.lab};
    if (fvmP1P0(q,u,c,a,where) )
      {
        for (int i=0;i<3;i++)
          for (int j=0;j<3;j++)
        if (fabs(a[i][j]) >= 1e-30)
          { Aij[make_pair(ii[i],ii[j])]+=a[i][j];
            }
      }
      }
    amorse=  new MatriceMorse<R>(Th.nv,Th.nv,Aij,false);
  }
  sparse_mat->pUh=0;
  sparse_mat->pVh=0;
  sparse_mat->A.master(amorse);
  sparse_mat->typemat=(amorse->n == amorse->m) ? TypeSolveMat(TypeSolveMat::GMRES) : TypeSolveMat(TypeSolveMat::NONESQUARE);  // none square matrice (morse)
  ⋆mp=mps;
  if(verbosity>3) { cout << "  End Build MatrixUpWind : " << endl;}
  return sparse_mat;
}
class Init { public:
  Init();
};
 Init init;
 Init::Init()
   {
     cout << " lood: init Mat Chacon " << endl;
     Global.Add("MatUpWind0","(", new OneOperatorCode<MatrixUpWind0 >( ));
   }

C.4 Add a new finite element

First read the section 12 of the appendix, we add two new finite elements examples in the directory examples++-load.

Bernardi Raugel The Bernardi-Raugel finite element is make to solve Navier Stokes equation in u,p formulation, and the velocity space PKbr is minimal to prove the inf-sup condition with piecewise contante pressure by triangle.

the finite element space Vh is

           1   2                   br
Vh = {u ∈ H (Ω );  ∀K  ∈ Th,u|K ∈ P K}
where
Pbr = span {λK ek}i=1,2,3,k=1,2 ∪ {λK λK nK }i=1,2,3
  K          i               i  i+1 i+2
with notation 4 = 1,5 = 2 and where λiK are the barycentric coordonnate of the triangle K, (ek)k=1,2 the canonical basis of ℝ2 and nkK the outer normal of triangle K opposite to vertex k.


// The P2BR finite element : the Bernadi Raugel Finite Element
// F. Hecht, decembre 2005
// -------------
// See Bernardi, C., Raugel, G.: Analysis of some finite elements for the Stokes problem. Math. Comp. 44, 71-79 (1985).
// It is a 2d coupled FE
// the Polynomial space is P12 + 3 normals bubbles edges function (P2)
// the degre of freedom is 6 values at of the 2 componantes at the 3 vertices
// and the 3 flux on the 3 edges
// So 9 degrees of freedom and N= 2.
// ----------------------- related files:
// to check and validate : testFE.edp
// to get a real example : NSP2BRP0.edp
// ------------------------------------------------------------
// -----------------------
#include "error.hpp"
#include "AFunction.hpp"
#include "rgraph.hpp"
using namespace std;
#include "RNM.hpp"
#include "fem.hpp"
#include "FESpace.hpp"
#include "AddNewFE.h"
namespace  Fem2D {
  class TypeOfFE_P2BRLagrange : public  TypeOfFE { public:
    static int Data[];
    TypeOfFE_P2BRLagrange(): TypeOfFE(6+3+0,
      2,
      Data,
      4,
      1,
      6+3⋆(2+2),  // nb coef to build interpolation
      9,  // np point to build interpolation
      0)
    {
 ....   // to long see the source
     }
    void FB(const bool ⋆ whatd, const Mesh & Th,const Triangle & K,const R2 &P, RNMK_ & val) const;
    void TypeOfFE_P2BRLagrange::Pi_h_alpha(const baseFElement & K,KN_<double> & v) const;
  } ;
   // on what nu df on node node of df
 int TypeOfFE_P2BRLagrange::Data[]={
   0,0, 1,1, 2,2,  3,4,5,
   0,1, 0,1, 0,1,  0,0,0,
   0,0, 1,1, 2,2,  3,4,5,
   0,0, 0,0, 0,0,  0,0,0,
   0,1, 2,3, 4,5,  6,7,8,
   0,0
};
void TypeOfFE_P2BRLagrange::Pi_h_alpha(const baseFElement & K,KN_<double> & v) const
  {
    const Triangle & T(K.T);
    int k=0;
     // coef pour les 3 sommets fois le 2 composantes
    for (int i=0;i<6;i++)
      v[k++]=1;
     // integration sur les aretes
    for (int i=0;i<3;i++)
      {
        R2 N(T.Edge(i).perp());
N  ⋆= T.EdgeOrientation(i)⋆0.5 ;
        v[k++]= N.x;
        v[k++]= N.y;
        v[k++]= N.x;
        v[k++]= N.y;
      }
  }
  void TypeOfFE_P2BRLagrange::FB(const bool ⋆ whatd,const Mesh & ,const Triangle & K,const R2 & P,RNMK_ & val) const
  {
 ....   // to long see the source
  }
// ---- cooking to add the finite elemet to freefem table --------
// a static variable to def the finite element
  static TypeOfFE_P2BRLagrange P2LagrangeP2BR;
   // now adding FE in FreeFEm++ table
  static AddNewFE P2BR("P2BR",&P2LagrangeP2BR);
// --- end cooking
}  // end FEM2d namespace
 

A way to check the finite element


load "BernadiRaugel"
// a macro the compute numerical derivative
macro DD(f,hx,hy) ( (f(x1+hx,y1+hy)-f(x1-hx,y1-hy))/(2⋆(hx+hy)))  //
mesh Th=square(1,1,[10⋆(x+y/3),10⋆(y-x/3)]);
real x1=0.7,y1=0.9, h=1e-7;
int it1=Th(x1,y1).nuTriangle;
fespace Vh(Th,P2BR);
Vh [a1,a2],[b1,b2],[c1,c2];
for (int i=0;i<Vh.ndofK;++i)
cout << i << " " << Vh(0,i) << endl;
for (int i=0;i<Vh.ndofK;++i)
{
  a1[]=0;
  int j=Vh(it1,i);
  a1[][j]=1;  // a bascis functions
  plot([a1,a2], wait=1);
  [b1,b2]=[a1,a2];  // do the interpolation
  c1[] = a1[] - b1[];
  cout << " ---------" << i << " " << c1[].max << " " << c1[].min << endl;
  cout << " a = " << a1[] <<endl;
  cout << " b = " << b1[] <<endl;
  assert(c1[].max < 1e-9 && c1[].min > -1e-9); // check if the interpolation is correct
  // check the derivative and numerical derivative
  cout << " dx(a1)(x1,y1) = " << dx(a1)(x1,y1) << " == " << DD(a1,h,0) << endl;
  assert( abs(dx(a1)(x1,y1)-DD(a1,h,0) ) < 1e-5);
  assert( abs(dx(a2)(x1,y1)-DD(a2,h,0) ) < 1e-5);
  assert( abs(dy(a1)(x1,y1)-DD(a1,0,h) ) < 1e-5);
  assert( abs(dy(a2)(x1,y1)-DD(a2,0,h) ) < 1e-5);
}
 

A real example using this finite element, just a small modification of the NSP2P1.edp examples, just the begenning is change to


load "BernadiRaugel"
real s0=clock();
mesh Th=square(10,10);
fespace Vh2(Th,P2BR);
fespace Vh(Th,P0);
Vh2 [u1,u2],[up1,up2];
Vh2 [v1,v2];

And the plot instruction is also change because the pressure is constant, and we cannot plot isovalues.

Morley See the example.

C.5 Add a new sparse solver

I will show the sketch of the code, see the full code from SuperLU.cpp or NewSolve.cpp.

First the include files:


#include  <iostream>
using namespace std;
#include "rgraph.hpp"
#include "error.hpp"
#include "AFunction.hpp"
// #include "lex.hpp"
#include "MatriceCreuse_tpl.hpp"
#include "slu_ddefs.h"
#include "slu_zdefs.h"

A small template driver to unified the double and Complex version.


template <class R> struct SuperLUDriver
{
};
template <> struct SuperLUDriver<double>
{
  ....  double version
};
template <> struct SuperLUDriver<Complex>
{
....  Complex version
};

To get Matrix value, we have just to remark that the Morse Matrice the storage, is the SLU_NR format is the compressed row storage, this is the transpose of the compressed column storage.

So if AA is a MatriceMorse you have with SuperLU notatiopn.


     n=AA.n;
     m=AA.m;
     nnz=AA.nbcoef;
     a=AA.a;
     asub=AA.cl;
     xa=AA.lg;
     options.Trans = TRANS;
     Dtype_t R_SLU = SuperLUDriver<R>::R_SLU_T();
     Create_CompCol_Matrix(&A, m, n, nnz, a, asub, xa, SLU_NC, R_SLU, SLU_GE);

To get vector infomation, to solver the linear solver x = A-1b


   void Solver(const MatriceMorse<R> &AA,KN_<R> &x,const KN_<R> &b) const
 {
 ....
    Create_Dense_Matrix(&B, m, 1, b, m, SLU_DN, R_SLU, SLU_GE);
    Create_Dense_Matrix(&X, m, 1, x, m, SLU_DN, R_SLU, SLU_GE);
 ....
 }

The two BuildSolverSuperLU function, to change the default sparse solver variable

DefSparseSolver<double>::solver


MatriceMorse<double>::VirtualSolver ⋆
BuildSolverSuperLU(const MatriceMorse<double> ⋆A,int strategy,double tgv, double eps,
                   double  tol_pivot,double tol_pivot_sym,
            int NbSpace,int itmax ,const  void ⋆ precon, ?void ⋆ stack)
{
    if(verbosity>9)
       cout << " BuildSolverSuperLU<double>" << endl;
    return new SolveSuperLU<double>(⋆A,strategy,tgv,eps,tol_pivot,tol_pivot_sym);
}
MatriceMorse<Complex>::VirtualSolver ⋆
BuildSolverSuperLU(const MatriceMorse<Complex> ⋆A,int strategy,double tgv, double eps,
           double tol_pivot,double tol_pivot_sym,
    int NbSpace,int itmax ,const  void ⋆ precon, void ⋆ stack)
{
    if(verbosity>9)
    cout << " BuildSolverSuperLU<Complex>" << endl;
    return new SolveSuperLU<Complex>(⋆A,strategy,tgv,eps,tol_pivot,tol_pivot_sym);
}

The link to FreeFem++


class Init { public:
    Init();
};
  To set the  2 default sparse solver double and complex:
DefSparseSolver<double>::SparseMatSolver SparseMatSolver_R ; ;
DefSparseSolver<Complex>::SparseMatSolver SparseMatSolver_C;

To save the default solver type


TypeSolveMat::TSolveMat  TypeSolveMatdefaultvalue=TypeSolveMat::defaultvalue;

To reset to the default solver, call this function:


bool SetDefault()
{
    if(verbosity>1)
  cout << " SetDefault sparse to default" << endl;
    DefSparseSolver<double>::solver =SparseMatSolver_R;
    DefSparseSolver<Complex>::solver =SparseMatSolver_C;
    TypeSolveMat::defaultvalue =TypeSolveMat::SparseSolver;
}

To set the default solver to superLU, call this function:


bool SetSuperLU()
{
    if(verbosity>1)
  cout << " SetDefault sparse solver to SuperLU" << endl;
    DefSparseSolver<double>::solver  =BuildSolverSuperLU;
    DefSparseSolver<Complex>::solver =BuildSolverSuperLU;
    TypeSolveMat::defaultvalue =TypeSolveMatdefaultvalue;
}

To add new function/name defaultsolver,defaulttoSuperLUin freefem++, and set the default solver to the new solver., just do:


Init init;
Init::Init()
{
  SparseMatSolver_R= DefSparseSolver<double>::solver;
  SparseMatSolver_C= DefSparseSolver<Complex>::solver;
  if(verbosity>1)
    cout << "\n Add: SuperLU,  defaultsolver defaultsolverSuperLU" << endl;
  TypeSolveMat::defaultvalue=TypeSolveMat::SparseSolver;
  DefSparseSolver<double>::solver =BuildSolverSuperLU;
  DefSparseSolver<Complex>::solver =BuildSolverSuperLU;
   // test if the name "defaultsolver" exist in freefem++
  if(! Global.Find("defaultsolver").NotNull() )
    Global.Add("defaultsolver","(",new OneOperator0<bool>(SetDefault));
  Global.Add("defaulttoSuperLU","(",new OneOperator0<bool>(SetSuperLU));
}

To compile superlu.cpp, just do:

  1. download the SuperLu 3.0 package and do
    curl   http://crd.lbl.gov/~xiaoye/SuperLU/superlu_3.0.tar.gz  -o superlu_3.0.tar.gz  
    tar xvfz superlu_3.0.tar.gz  
    go SuperLU_3.0 directory  
    $EDITOR  make.inc  
    make

  2. In directoy include do to have a correct version of SuperLu header due to mistake in case of inclusion of double and Complex version in the same file.
    tar xvfz ../SuperLU_3.0-include-ff.tar.gz

    I will give a correct one to compile with freefm++.

    To compile the freefem++ load file of SuperLu with freefem do: some find like :

    ff-c++ SuperLU.cpp  -L$HOME/work/LinearSolver/SuperLU_3.0/ -lsuperlu_3.0

    And to test the simple example:

    A example:


    load "SuperLU"
    verbosity=2;
    for(int i=0;i<3;++i)
    {
    // if i == 0 then SuperLu solver
    // i == 1 then GMRES solver
    // i == 2 then Default solver
      {
        matrix A =
          [[ 0, 1, 0, 10],
           [ 0,  0,  2, 0],
           [ 0, 0, 0,  3],
           [ 4,0 , 0, 0]];
        real[int] xx = [ 4,1,2,3], x(4), b(4);
        b = A⋆xx;
        cout << b << " " << xx << endl;
        set(A,solver=sparsesolver);
        x = A^-1⋆b;
        cout << x << endl;
      }
      {
        matrix<complex> A =
          [[ 0, 1i, 0, 10],
           [ 0 ,  0,  2i, 0],
           [ 0, 0, 0,  3i],
           [ 4i,0 , 0, 0]];
        complex[int] xx = [ 4i,1i,2i,3i], x(4), b(4);
        b = A⋆xx;
        cout << b << " " << xx << endl;
        set(A,solver=sparsesolver);
        x = A^-1⋆b;
        cout << x << endl;
      }
      if(i==0)defaulttoGMRES();
      if(i==1)defaultsolver();
    }

    To Test do for exemple:

    FreeFem++ SuperLu.edp