Mercurial > fem-fenics-eugenio
view doc/doc.tex @ 194:c29ac833819f
Updated version of the documentation.
author | gedeone-octave <marcovass89@hotmail.it> |
---|---|
date | Sun, 17 Nov 2013 23:52:24 +0000 |
parents | 01b5ec639095 |
children | 6b82ac887c99 |
line wrap: on
line source
\documentclass[a4paper,10pt]{book} \usepackage[utf8x]{inputenc} \usepackage{graphicx} \usepackage{amsmath} \usepackage{amssymb} \usepackage{caption} \usepackage{subcaption} \usepackage{cclicenses} \usepackage{url} \usepackage{listings} \usepackage{hyperref} \usepackage{listings} \usepackage{color} \usepackage[english]{babel} \usepackage[T1]{fontenc} \definecolor{dkgreen}{rgb}{0,0.6,0} \definecolor{gray}{rgb}{0.5,0.5,0.5} \definecolor{mauve}{rgb}{0.58,0,0.82} \lstset{frame=tb, language=Octave, aboveskip=3mm, belowskip=3mm, showstringspaces=false, columns=flexible, basicstyle={\small\ttfamily}, numbers=left, numberstyle=\tiny\color{gray}, keywordstyle=\color{blue}, commentstyle=\color{dkgreen}, stringstyle=\color{mauve}, breaklines=true, breakatwhitespace=true tabsize=3 } \author{Marco Vassallo} \title{ \textbf{Fem-fenics} \\ \bigskip \textit{General purpose Finite Element library \\ for GNU-Octave}\\ \bigskip \textsc{ Work in progress \\(help and remarks are welcome)}} \begin{document} \maketitle \tableofcontents \chapter{Introduction} Fem-Fenics is an open-source package for the resolution of partial differential equations with Octave. The project has been developed during the Google Summer of Code 2013 with the help and the sustain of the GNU-Octave community under the supervision of prof. De Falco.\\ The report is structured as follows: \begin{itemize} \item in chapter \ref{intr} we provide a simple reference guide for beginners \item in chapter \ref{impl} is presented a detailed explanation of the relevant parts of the program. In this way, the interested reader can see what there is ``behind'' and expecially anyone interested in it can learn quickly how it is possible to extend the code and contribute to the project. \item in chapter \ref{exem} more examples are provided. For a lot of them, we present the octave script alongside with the code for Fenics (in C++ and/or Python) in order to provide the user with a quick reference guide. \end{itemize} If you think that going inside the report could be boring, it is available a wiki at \begin{center} \url{http://wiki.octave.org/Fem-fenics} \end{center} while if you want to see how the project has grown during the time you can give a look at \begin{center} \url{http://gedeone-gsoc.blogspot.com/} \end{center} Finally, the API is available at the following address \begin{center} \url{http://octave.sourceforge.net/fem-fenics/overview.html} \end{center} \chapter{Introduction to Fem-fenics}\label{intr} \section{Installation} Fem-fenics is an external package for Octave, which means that it can be installed only once that Octave has been successfully installed on the PC. Furthermore, as Fem-fenics is based on Fenics, it is also needed a running version of the latter. They can be easily installed following the guidelines provided on the official Octave \cite{instoctave} and Fenics \cite{instfenics} websites. Once that Octave and Fenics are correctly installed, to install Fem-fenics launch Octave (which now is provided with a new amazing GUI) and type \begin{verbatim} >> pkg install fem-fenics -forge \end{verbatim} That's all! For any problem during the installation don't hesitate to contact us. To be sure that everything is working fine, load the fem-fenics pkg and run one of the examples provided within the package: \begin{verbatim} >> pkg load fem-fenics >> femfenics_examples() \end{verbatim} For a description of the examples, look at chapter \ref{exem}. \\ \\ \textbf{NOTE} For completing the installation process successfully, the form compiler FFC and the header file dolfin.h should also be available on the machine. They are managed automatically by Fenics if it is installed as a binary package or with Dorsal. If it has been done manually, please be sure that they are available before starting the installation of Fem-fenics. \section{General layout and first example} A generic problem has to be solved in two steps: \begin{enumerate} \item a \textbf{.ufl file} where the abstract problem is described: this file has to be written in Unified Form Language (UFL), which is a domain specific language for defining discrete variational forms and functionals in a notation close to pen-and-paper formulation. UFL is easy to learn, and the User manual provides explanations and examples \cite{ufl}. \item a script file \textbf{.m} where the abstract problem is imported and a specific problem is implemented and solved: this is the script file where the fem-fenics functions described in the following chapters are used. \end{enumerate} We provide immediately a simple example in order to familiarize the user with the code. \paragraph{The Poisson equation} In this example, we show how it is possible to solve the Poisson equation with mixed Boundary Conditions. If we indicate with $\Omega$ the domain and with $\Gamma = \Gamma_{N} \cup \Gamma_{D}$ the boundaries, the problem can be expressed as \begin{align*} \Delta u &= f \qquad \text{on } \Omega \\ u &= 0 \qquad \text{on } \Gamma_{D} \\ \nabla u \cdot n &= g \qquad \text{on } \Gamma_{N} \end{align*} where $f, \, g$ are data which represent the source for and the flux of the scalar variable $u$. A possible variational formulation of the problem is: \\ find $u \in H_{0, \Gamma_{D}}^{1} :$ \begin{align*} a(u, v) &= L(v) \qquad \forall v \in H_{0, \Gamma_{D}}^{1} \\ a(u, v) &= \int_{\Omega} \nabla u \cdot \nabla v \\ L(v) &= \int_{\Omega} f v + \int_{\Gamma_{N}} g v \\ \end{align*} The abstract problem can thus be written in the \verb|Poisson.ufl| file immediately. The only thing that has to be specified at this stage is the space of Finite Elements used for the discretization of $H_{0, \Gamma_{D}}^{1}$. In this case, we choose the space of continuous lagrangian polynomial of degree one \verb|FiniteElement("Lagrange", triangle, 1)|, but many more possibilities are available. \subparagraph{Poisson.ufl} \begin{lstlisting} element = FiniteElement("Lagrange", triangle, 1) u = TrialFunction(element) v = TestFunction(element) f = Coefficient(element) g = Coefficient(element) a = inner(grad(u), grad(v))*dx L = f*v*dx + g*v*ds \end{lstlisting} It is always a good idea to check if the ufl code is correctly written before importing it into Octave: \begin{verbatim} >> ffc -l dolfin Poisson.ufl \end{verbatim} shouldn't produce any error. We can now implement and solve a specific instance of the Poisson problem in Octave. The parameters are setted as follow \begin{itemize} \item $\Omega = [0, 1]\times[0, 1]$ \item $\Gamma_{D} = {(0, y) \cup(1, y)} \ \subset \partial\Omega$ \item $\Gamma_{N} = {(x, 0) \cup(x, 1)} \ \subset \partial\Omega$ \item $f = 10 \exp \dfrac{(x-0.5)^{2} + (y-0.5)^{2}}{0.02}$ \item $g = \sin(5x)$ \end{itemize} As a first thing we need to load into Octave the pkgs previously installed \begin{verbatim} pkg load fem-fenics msh \end{verbatim} The ufl file can thus be imported inside Octave. For every specific element defined inside the ufl file there is a specific function which stores it for later use \begin{itemize} \item \verb$ufl_import_FunctionSpace ('Poisson')$ is a function which looks for the finite element space defined inside the file called Laplace.ufl; if everything is ok, it generates a function which we will use later \item \verb$ufl_import_BilinearForm ('Poisson')$ is a function which looks for the rhs of the equation, i.e. for the bilinear form defined inside Laplace.ufl; \item \verb$ufl_import_LinearForm ('Poisson')$ is a function which looks for the linear form. \end{itemize} In some cases one could be interested in using these functions separately but if, as in our example, all the three elements are defined in the same ufl file (and only in this case), the \verb$import_ufl_Problem ('Poisson')$ can be used; it generates at once all the three functions described above \begin{verbatim} ufl_import_Problem ('Poisson'); \end{verbatim} To set the concrete elements which define the problem, the first things to do is to create a mesh. It can be managed easily using the msh pkg. For a uniform structured mesh \begin{verbatim} x = y = linspace (0, 1, 33); msho = msh2m_structured_mesh (x, y, 1, 1:4); \end{verbatim} Once that the mesh is available, we can thus initialize the fem-fenics mesh using the function \verb$Mesh ()$: \begin{verbatim} mesh = Mesh (msho); \end{verbatim} If instead of an Octave mesh we had a mesh stored in a different format, the program dolfin-convert can try to convert it to the dolfin xml format. Furthermore, \verb$Mesh ()$ accepts as arguments also a string with the name of the dolfin xml file which contains the mesh. For example, for a mesh saved in the gmsh format, do as follows: \begin{verbatim} Shell: dolfin-convert msh.gmsh msh.xml \end{verbatim} and then inside the Octave script: \begin{verbatim} mshd = fem_init_mesh ('msh.xml'); \end{verbatim} To initialize the functional space, we have to specify as argument only the fem-fenics mesh, because the finite element type and the polynomial degree have yet been specified in the ufl file: \begin{verbatim} V = FunctionSpace('Poisson', mesh); \end{verbatim} Essential BC can now be applied using \verb$DirichletBC ()$; this function receives as argument the functional space, a function handle which specifies the value to set, and the label of the sides where the BC applies. In this case, homogenous boundary conditions hold on the left and right side of the square \begin{verbatim} bc = DirichletBC(V, @(x, y) 0.0, [2; 4]); \end{verbatim} The last thing to do before solving the problem, is to set the coefficients specified in the ufl file. It is important that they are called in the same way as in the ufl file: the source term 'f' and the normal flux 'g'. To set them, the function \verb$Expression ()$ can be used passing as argument a string, which specifies the name of the coefficient, and a function handle with the value required: \begin{verbatim} f = Expression ('f', @(x,y) 10*exp(-((x - 0.5)^2 + (y - 0.5)^2) / 0.02)); g = Expression ('g', @(x,y) sin (5.0 * x)); \end{verbatim} Another possibility for dealing with the coefficients defined in the ufl file would have been to use the function \verb$Constant ()$ or \verb$Function ()$. The coefficients can thus be used together with the FunctionSpace to set the Bilinear and the Linear form \begin{verbatim} a = BilinearForm ('Poisson', V, V); L = LinearForm ('Poisson', V, f, g); \end{verbatim} The discretized representation of our operator is obtained using the functions \verb$assemble ()$ or \verb$assemble_system ()$, which also allow to specify the BC(s) to apply. Whenever possible, it is better to use the \verb$assemble_system ()$ function because it keeps the symmetry of the matrix while setting the entries for the BC: \begin{verbatim} [A, b] = assemble_system (a, L, bc); \end{verbatim} Here A is a sparse matrix and b is a column vector. All the functionalities available within Octave can now be exploited to solve the linear system. The easisest possibility is the backslash command: \begin{verbatim} u = A \ b; \end{verbatim} Once that the solution has been obtained, the vector is converted into a fem-fenics function and plotted \verb$plot ()$ or saved \verb$save ()$ in the vtu format. \begin{verbatim} u = Function ('u', V, sol); save (u, 'poisson') plot (u); \end{verbatim} The complete code for the Poisson problem is reported below, while in figure \ref{Poissonfig} the output is presented. \begin{figure} \begin{center} \includegraphics[height=7 cm,keepaspectratio=true]{./Fem-fenics_poisson.png} \caption{The result for the Poisson equation} \label{Poissonfig} \end{center} \end{figure} \subparagraph{Poisson.m} \begin{lstlisting} #load the pkg and import the ufl problem pkg load fem-fenics msh import_ufl_Problem ('Poisson') # Create the mesh and define function space x = y = linspace (0, 1, 33); mesh = Mesh(msh2m_structured_mesh (x, y, 1, 1:4)); V = FunctionSpace('Poisson', mesh); # Define boundary condition and source term bc = DirichletBC(V, @(x, y) 0.0, [2;4]); f = Expression ('f', @(x,y) 10*exp(-((x - 0.5)^2 + (y - 0.5)^2) / 0.02)); g = Expression ('g', @(x,y) sin (5.0 * x)); #Create the Bilinear and the Linear form a = BilinearForm ('Poisson', V, V); L = LinearForm ('Poisson', V, f, g); #Extract the matrix and compute the solution [A, b] = assemble_system (a, L, bc); sol = A \ b; u = Function ('u', V, sol); # Save solution in VTK format and plot it save (u, 'poisson') plot (u); \end{lstlisting} \chapter{Implementation}\label{impl} Fem-fenics aims to fill a gap in Octave: even if there are packages for the creation of mesh \cite{msh}, for the postprocessing of data \cite{fpl} and for the resolution of some specific pde \cite{secs1d} \cite{bim}, no general purpose finite element library is available. The goal of the project is thus to provide a package which can be used to solve user defined problems and which is able to exploit the functionality provided with Octave. \begin{figure} \begin{center} \includegraphics[height=10 cm,keepaspectratio=true]{./code_layout.png} \caption{General layout of the package} \label{Codelayout} \end{center} \end{figure} Instead of writing a library from scratch, an interface to one of the finite element library which are already available has been created. Among the many libraries taken into account, the one whch was best suited for our purposes seemed to be the FEniCS project. It ``is a collection of free, open source, software components with the common goal to enable automated solution of pde.'' In particular, Dolfin is the C++/Python interface of FEniCS, providing a consistent Problem Solving Environment for ODE and PDE. The idea has been to create wrappers in Octave for Dolfin, in a similar way to what it has been done for Python. This is a very natural choice, because Octave is mainly written in script language and in C++. It is in fact possible to implement an Octave interpreter function in C++ through the native oct-file interface or, coversely, to use Octave's Matrix/Array Classes in a C++ application. The works can be summirized as follows (fig. \ref{Codelayout}): the elements already available in Octave for the resolution of PDE (Mesh and Linear Algebra) have been exploited, and wrappers to FEniCS functions added. To allow exchanges between this programs, the necessary functions for converting an Octave mesh/matrix into a FEniCS one and viceversa have been written. Two main ideas have guided us throughout the realization of the pkg: \begin{itemize} \item keep the syntax as close as possible to the original one in Fenics \item make the interface as simple as possible. \end{itemize} \section{General layout of a class} Seven new classes are implemented for dealing with FEniCS objects and for using them inside Octave: \begin{itemize} \item \textbf{boundarycondition} stores and builds a dolfin::DirichletBC \item \textbf{coefficient} stores an expression object which is used for the evaluation of user defined values \item \textbf{expression} is needed for internal use only as explained below \item \textbf{form} stores a general dolfin::Form and can be used either for a dolfn::BilinearForm as well as for a dolfin::LinearForm \item \textbf{function} for the dolfin::Function objects \item \textbf{functionspace} stores the user defined FunctionSpace \item \textbf{mesh} converts a PDE-tool like mesh structure in a dolfin::Mesh \end{itemize} The general layout of a class is quite simple and its main purpose is the storage of the corresponding FEniCS objects. They are managed throughout \verb$boost::shared_ptr< >$ to the corresponding FEniCS type. All the classes also implement a constructor which takes as argument the corresponding dolfin type, and a default constructor which is necessary to register a type in the Octave interpreter. All these classes derive from \verb$octave_base_value$. For example, the form class implementation follows, while classes which differs from the general layout are presented below. \begin{lstlisting} #ifndef _FORM_OCTAVE_ #define _FORM_OCTAVE_ #include <memory> #include <vector> #include <dolfin.h> #include <octave/oct.h> class form : public octave_base_value { public: form () : octave_base_value () {} form (const dolfin::Form _frm) : octave_base_value (), frm (new dolfin::Form (_frm)) {} form (boost::shared_ptr <const dolfin::Form> _frm) : octave_base_value (), frm (_frm) {} void print (std::ostream& os, bool pr_as_read_syntax = false) const { os << "Form " << ": is a form of rank " << frm->rank () << " with " << frm->num_coefficients () << " coefficients" << std::endl; } ~form(void) {} bool is_defined (void) const { return true; } const dolfin::Form & get_form (void) const { return (*frm); } const boost::shared_ptr <const dolfin::Form> & get_pform (void) const { return frm; } private: boost::shared_ptr <const dolfin::Form> frm; DECLARE_OCTAVE_ALLOCATOR; DECLARE_OV_TYPEID_FUNCTIONS_AND_DATA; }; static bool form_type_loaded = false; DEFINE_OCTAVE_ALLOCATOR (form); DEFINE_OV_TYPEID_FUNCTIONS_AND_DATA (form, "form", "form"); #endif \end{lstlisting} \subsection{mesh} In addition to usual methods, the mesh class implemens functionalities which allow to deal with meshes as they are currently available with the msh pkg, i.e. in the (p, e, t) format. Dolfin instead use a xml format for storing informations relative to a mesh. We have therefore added a constructor \begin{verbatim} mesh (Array<double>& p, Array<octave_idx_type>& e, Array<octave_idx_type>& t); \end{verbatim} which accept as input a mesh in (p, e, t) format and convert it into a xml one. Before exploring the code in more details,the main differences between the two storing formats are presented using the very simple but rather instructive example of a unit square mesh with just two elements, fig. \ref{mesh}. \paragraph{pet} \begin{figure} \begin{center} \includegraphics[height=5 cm,keepaspectratio=true]{./mesh_1.png} \caption{The (very) simple mesh for our example} \label{mesh} \end{center} \end{figure} A mesh is represented using the three matrices p, e, t, and, using msh, we can easily obtain the mesh for our example typing \begin{verbatim} mesh = msh2m_structured_mesh ([0 1], [0 1], 1, [11 12 12 13]) \end{verbatim} The matrix p stores information about the coordinates of the vertices \begin{verbatim} >> mesh.p \end{verbatim} $ \begin{array}{rrrrl} \; 0 & 0 & 1 & 1 & \quad \text{x-coordinates} \\ \; 0 & 1 & 0 & 1 & \quad \text{y-coordinates} \\ \end{array} $ \newline \newline Thus the vertex in the $n^{th}$ column will be labelled as the vertex number $n$, and so on.. The matrix t stores information about the connectivity \begin{verbatim} >> mesh.t \end{verbatim} $ \begin{array}{rrl} \; 1 & 1 & \quad \text{number of the first vertex of the element} \\ \; 3 & 4 & \quad \text{number of the second vertex of the element} \\ \; 4 & 2 & \quad \text{number of the third vertex of the element} \\ \; 0 & 0 & \\ \end{array} $ \newline \newline The first element is thus the one obtained connecting vertices 1-3-4 and so on.. The matrix e stores information related to every side edge, like the number of the vertices of the boundary elements, and the number of the geometrical border containing the edge. (This last information could be very useful later in this project, when we will have to find a way to apply boundary conditions to our problem) \begin{verbatim} >> mesh.e \end{verbatim} $ \begin{array}{rrrrl} \; 1 & 3 & 2 & 1 & \quad \text{first vertex of the side edge} \\ \; 3 & 4 & 4 & 2 & \quad \text{second vertex of the side edge} \\ \; 0 & 0 & 0 & 0 & \text{} \\ \; 0 & 0 & 0 & 0 & \text{} \\ \; 11 & 12 & 12 & 13 & \quad \text{label of the geometrical border containing the edge} \\ \; 0 & 0 & 0 & 0 & \text{}\\ \; 1 & 1 & 1 & 1 & \text{} \\ \end{array} $ \newline \newline The side edge between vertex 1-3 is labelled 11, between 3-4 12... \paragraph{dolfin xml} A mesh is an object of the dolfin::Mesh class which stores information only about the coordinates of the vertices (like p) and the information about the connectivity (like t). A mesh can thus be manipulated using the functions and the methods of the class. The information about boundaries is not directly stored in the mesh. The mesh used in th example is stored as \begin{verbatim} <?xml version="1.0"?> <dolfin xmlns:dolfin="http://fenicsproject.org"> <mesh celltype="triangle" dim="2"> <vertices size="4"> <vertex index="0" x="0.000e+00" y="0.000e+00" /> <vertex index="1" x="0.000e+00" y="1.000e+00" /> <vertex index="2" x="1.000e+00" y="0.000e+00" /> <vertex index="3" x="1.000e+00" y="1.000e+00" /> </vertices> <cells size="2"> <triangle index="0" v0="0" v1="2" v2="3" /> <triangle index="1" v0="0" v1="1" v2="3" /> </cells> </mesh> </dolfin> \end{verbatim} \subsection{Shared pointer} This is done for two main reasons: in our program we have to refer in several places to resources which are built dynamically (e.g. a FunctionSpace contain a reference to the Mesh) and we want that it is destroyed only when the last references is destroyed. they are widely used inside DOLFIN and it is thus easier to deal with it if we use the same types. \iffalse \paragraph{functionspace} A dolfin::FunctionSpace is defined by specifying a mesh and the type of the finite element which we want to use. The mesh is handled by our class, while the FE are specified inside the .ufl file. Possible choices are: Argyris ARG Arnold–Winther AW Brezzi–Douglas–Marini BDM Crouzeix–Raviart CR Discontinuous Lagrange DG Hermite HER Lagrange CG Mardal–Tai–Winther MTW Morley MOR Nédélec 1st kind H (curl) N1curl Nédélec 2nd kind H (curl) N2curl Raviart–Thomas RT Once that the .ufl file has been compiled, in the output file it is defined a namespace where we can find the corresponding FunctionSpace and which only needs a mesh to be initialized. \paragraph{boundary condition} This class is used for dealing with essential BC. We need three elements for the definition of an object of type dolfin::DirichletBC : a dolfin::FunctionSpace where the BC is defined, which is handled by our class functionspace described above a dolfin::Function or a dolfin::Expression which represents the value that we want to assign. This is done deriving a new class from dolfin::Expression and overloading the eval function. a dolfin::subdomain which specifies where we want to apply the BC. \paragraph{coefficient} This class is used to store a std::string where we save the name of the coefficient, and an object of type expression, which is a class derived from dolfin::Expression. The class expression would be explained in more detail in the next post. \fi \section{General layout of a function} \iffalse The functions which we describe below are used to set and get the values for the variable that are commonly used in a FEM problem. The main functions which are available right now are: \verb$fem_init_mesh$: this function can take as input either a string which specify the name of the file to read or a mesh produced with the msh pkg \verb$fem_bc$: take as input the functional space, a function handle which specifies the value that we want to set on the boundary, and the label of the boundaries where we want to apply the bc coefficient: take as input the name of one of the coefficients defined in the .ufl file and a function handle with the value that we want to set The following functions work fine but at the moment are directly compiled including the output file obtained from the .ufl file (which should be avoided as discussed below) : \verb$fem_fs$: take as input the mesh where we want to define the functional space \verb$fem_rhs$: take as input the function space where the bilinear form is defined, a list of coefficients and a list of BCs that we want to apply. It returns a sparse matrix which contains the discretization of the bilinear operator \verb$fem_lhs$: take the same input as the function \verb$fem_rhs$. If you have essential BC in the problem, you need to pass them also to this function in such a way that the DOF corresponding to essential nodes are set to the right value The general layout of a function is very simple and is composed of 4 steps which we describe using an example: \begin{lstlisting} DEFUN_DLD (fem_fs, args, , "initialize a fs from a mesh") { // 1 read data const mesh & msho = static_cast<const mesh&> (args(0).get_rep ()); // 2 convert the data from octave to dolfin const dolfin::Mesh & mshd = msho.get_msh (); // 3 build the new object using dolfin boost::shared_ptr <const dolfin::FunctionSpace> g (new Laplace::FunctionSpace (mshd)); // 4 convert the new object from dolfin to Octave and return it octave_value retval = new functionspace(g); return retval; } \end{lstlisting} All the functions presented above follow this general structure, and thus here we present in detail only functions which present some differences. \paragraph{fem bc and fem coeff} These two functions take as input a function handle which cannot be directly evaluated by a dolfin function to set, respectively, the value on the boundary or the value of the coefficient. We have thus derived from dolfin::Expression a class "expression" which has as private member an octave function handle and which overloads the function eval(). In this way, an object of the class expression can be initialized throughout a function handle and can be used inside dolfin because "it is" a dolfin::Expression. \begin{lstlisting} class expression : public dolfin::Expression { public: expression (octave_fcn_handle & _f) : dolfin::Expression (), f (new octave_fcn_handle (_f)) {} void eval (dolfin::Array<double>& values, const dolfin::Array<double>& x) const { octave_value_list b; b.resize (x.size ()); for (int i = 0; i < x.size (); ++i) b(i) = x[i]; octave_value_list res = feval (f->function_value (), b); values[0] = res(0).double_value (); } private: octave_fcn_handle * f; }; \end{lstlisting} \paragraph{fem lhs and fem rhs} For these functions, we need to set the possible coefficient of the (bi)linear form and to apply the BC. For the coefficient: \begin{lstlisting} // read the coefficient const coefficient & cf = static_cast <const coefficient&> (args (i).get_rep ()); // get the correct position of the coefficient in the declaration of the problem std::size_t n = a.coefficient_number (cf.get_str ()); // extract the information.. we use again the expression class descripted above const boost::shared_ptr<const expression> & pexp = cf.get_expr (); //set the coefficient for the (bi)linear form a.set_coefficient (n, pexp); \end{lstlisting} For the BC At the moment we are simply using the method apply() of the class DirichletBC, but it does not preserve the symmetry of the matrix. We have thus planned to insert later a function which instead of the method apply() calls the function \verb$assemble_system()$ , which preserves the symmetry but which builds together the lhs and the rhs. The sparse matrix is then built from the dolfin::Matrix following the general guidelines of Octave. \paragraph{other function} SubSpace allows to extract a subspace from a vectorial one. For example, if our space is P2 x P0 we can extract the one or the other and then apply BC only where it is necessary. \verb$fem_eval$ takes as input a Function and a coordinate and returns a vector representing the value of the function at this point. for dealing with form of rank 0, i.e. with functional, we have now added the functions \verb$fem_create_functional$ to create it from a .ufl file. We have thus extended the function assemble which returns the corresponding double value. \verb$plot_2d$ and \verb$plot_3d$: these functions allow us to plot a function specifying a mesh and the value of the function at every node of the mesh. This is something which could be useful also outside of fem-fenics. \section{Implementation Details} The relevant implementation details which the user should know are: all the objects are managed using \verb$boost::shared_ptr <>$. It means that the same resource can be shared by more objects and useless copies should be avoided. For example, if we have two different functional spaces in the same problem, like with Navier-Stokes for the velocity and the pressure, the mesh is shared between them and no one has its own copy. The BC are imposed directly to the mesh setting to zero all the off diagonal elements in the corresponding line. This means that we could loose the symmetry of the matrix, if any. The simmetry is lost because we are using the method apply() of the class dolfin::DirichletBC. We have thus planned to insert later a function which instead of the method apply() calls the function \verb$assemble_system()$ , which preserves the symmetry of the system but which builds together the lhs and the rhs. The coefficient of the variational problem can be specified using either a fem-coefficient or a fem-function. They are different objects which behave in different ways: a fem-coefficient object overloads the eval() method of the dolfin::Expression class and it is evaluated at run time using the octave function feval (). A fem-function instead doesn't need to be evaluated because it is assembled copying element-by-element the values contained in the input vector. The relevant implementation details which the user should know are: all the objects are managed using \verb$boost::shared_ptr <>$. It means that the same resource can be shared by more objects and useless copies should be avoided. For example, if we have two different functional spaces in the same problem, like with Navier-Stokes f or the velocity and the pressure, the mesh is shared between them and no one has its own copy. The essential BC are imposed directly to the matrix with the command assemble(), which sets to zero all the off diagonal elements in the corresponding line, sets to 1 the diagonal element and sets to the exact value the rhs. This means that we could loose the symmetry of the matrix, if any. To avoid this problem and preserve the symmetry of the system it is available the \verb$assemble_system()$ command which builds at once the lhs and the rhs. The coefficient of the variational problem can be specified using either an Expression(), a Constant() or a Function(). They are different objects which behave in different ways: an Expession or a Constant object overloads the eval() method of the dolfin::Expression class and it is evaluated at run time using the octave function feval (). A Function instead doesn't need to be evaluated because it is assembled copying element-by-element the values contained in the input vector. We have split the construction of the form into two steps: We set all the coefficients of the form using the function which we create on the fly. They will be named \verb$ProblemName_BilinearForm$ or \verb$ProblemName_LinearForm$. Then we apply specific BC to the form using the assemble() function and we get back the matrix. If we are assembling the whole system and we want to keep the symmetry of the matrix (if any), we can instead use the command \verb$assemble_system$ (). Finally, if we are solving a non-linear problem and we need to apply essential BC, we should provide to the function also the vector with the tentative solution in order to modify the entries corresponding to the boundary values. This will be illustrated below in the HyperElasticity example. \fi \subsection{Mesh generation and conversion} \subsection{Sparse Matrices} \subsection{Polymorphism} \subsection{Code release} \subsection{Code on the fly} \iffalse For the creation on the fly of the code from the header file, Juan Pablo provided me a python code which makes a great job. I have spent some time adapting it to my problem, but when I finally got a working code, we realized that it was probably enough to use the functions available inside Octave because the problem was rather simple. The pyton code is however available here, while the final solution adopted can be found there. For the problem related with \verb$fem_init_env ()$, we are still working on it. \fi \subsection{Autoconf} In this section we want to discuss how we can write a config.ac and a Makefile.in files which: \begin{itemize} \item check if a program is available and stop if it is not \item check if a header file is available and issue a warning if not, but go ahead with the compilation \end{itemize} To reach this goal, we need two components: \paragraph{configure.ac} Is a file which checks whether the program/header is available or not and sets consequently the values of some variables. \begin{lstlisting} # Checks if the program mkoctfile is available and sets the variable HAVE_MKOCTFILE consequently AC_CHECK_PROG([HAVE_MKOCTFILE], [mkoctfile], [yes], [no]) # if mkoctfile is not available, it issues an error and stops the compilation if [test $HAVE_MKOCTFILE = "no"]; then AC_MSG_ERROR([mkoctfile required to install $PACKAGE_NAME]) fi #Checks if the header dolfin.h is available; if it is available, the value of the ac_dolfin_cpp_flags is substituted with -DHAVE_DOLFIN_H, otherwise it is left empty and a warning message is printed AC_CHECK_HEADER([dolfin.h], [AC_SUBST(ac_dolfin_cpp_flags,-DHAVE_DOLFIN_H) AC_SUBST(ac_dolfin_ld_flags,-ldolfin)], [AC_MSG_WARN([dolfin headers could not be found, some functionalities will be disabled, don't worry your package will still be working, though.])] ). # It generates the Makefile, using the template described below AC_CONFIG_FILES([Makefile]) \end{lstlisting} \paragraph{Makefile.ac} This file is a template for the Makefile, which will be automatically generated when the configure.ac file is executed. The values of the variable \verb$ac_dolfin_cpp_flags$ and \verb$ac_dolfin_ld_flags$ are substituted with the results obtained above: \begin{lstlisting} CPPFLAGS += @ac_dolfin_cpp_flags@ LDFLAGS += @ac_dolfin_ld_flags@ \end{lstlisting} In this way, if dolfin.h is available, CPPFLAGS contains also the flag \verb$-DHAVE_DOLFIN_H.$ \paragraph {program.cc} Our .cc program, should thus include the header dolfin.h only if \verb$-DHAVE_DOLFIN_H$ is defined at compilation time. For example \begin{lstlisting} #ifdef HAVE_DOLFIN_H #include <dolfin.h> #endif int main () { #ifndef HAVE_DOLFIN_H error("program: the program was built without support for dolfin"); #else /* Body of your function */ #endif return 0; } \end{lstlisting} \paragraph {Warning} If in the Makefile.in you write something like \begin{lstlisting} HAVE_DOLFIN_H = @HAVE_DOLFIN_H@ ifdef HAVE_DOLFIN_H CPPFLAGS += -DHAVE_DOLFIN_H LIBS += -ldolfin endif \end{lstlisting} it doesn't work because the variable \verb$HAVE_DOLFIN_H$ seems to be always defined, even if the header is not available. \chapter{More Advanced Examples}\label{exem} \iffalse With the following examples, we can see directly in action the new features and understand how they work. Navier-Stokes: we learn how to deal with a vector-field problem and how we can save the solution using the \verb$fem_save$ () function. We also use the fem pkg to generate a mesh using gmesh. Mixed-Poisson: we solve the Poisson problem presented in the previous posts using a mixed formulation, and we see how we can extract a scalar field from a vector one. HyperElasticity: we exploit the fsolve () command to solve a non-linear problem. In particular, we see how to use the assemble() function to apply BC also in this situation. Advection-Diffusion: we solve a time dependent problem using the lsode () command and save the solution using the pkg flp. For each problem, we refer the reader to the complete desciption on FEniCS or bim web-page, while here we highlight only the implementation detail relevant for our pkg. \fi \section{Navier-Stokes equation with Chorin-Temam projection algorithm} \paragraph{TentativeVelocity.ufl} \begin{lstlisting} # Copyright (C) 2010 Anders Logg # Define function spaces (P2-P1) V = VectorElement("CG", triangle, 2) Q = FiniteElement("CG", triangle, 1) # Define trial and test functions u = TrialFunction(V) v = TestFunction(V) # Define coefficients k = Constant(triangle) u0 = Coefficient(V) f = Coefficient(V) nu = 0.01 # Define bilinear and linear forms eq = (1/k)*inner(u - u0, v)*dx + inner(grad(u0)*u0, v)*dx + \ nu*inner(grad(u), grad(v))*dx - inner(f, v)*dx a = lhs(eq) L = rhs(eq) \end{lstlisting} \paragraph{PressureUpdate.ufl} \begin{lstlisting} # Copyright (C) 2010 Anders Logg # Define function spaces (P2-P1) V = VectorElement("CG", triangle, 2) Q = FiniteElement("CG", triangle, 1) # Define trial and test functions p = TrialFunction(Q) q = TestFunction(Q) # Define coefficients k = Constant(triangle) u1 = Coefficient(V) # Define bilinear and linear forms a = inner(grad(p), grad(q))*dx L = -(1/k)*div(u1)*q*dx \end{lstlisting} \paragraph{VelocityUpdate.ufl} \begin{lstlisting} # Copyright (C) 2010 Anders Logg # Define function spaces (P2-P1) V = VectorElement("CG", triangle, 2) Q = FiniteElement("CG", triangle, 1) # Define trial and test functions u = TrialFunction(V) v = TestFunction(V) # Define coefficients k = Constant(triangle) u1 = Coefficient(V) p1 = Coefficient(Q) # Define bilinear and linear forms a = inner(u, v)*dx L = inner(u1, v)*dx - k*inner(grad(p1), v)*dx \end{lstlisting} \paragraph{NS.m} \begin{lstlisting} pkg load fem-fenics msh import_ufl_Problem ("TentativeVelocity"); import_ufl_Problem ("VelocityUpdate"); import_ufl_Problem ("PressureUpdate"); # We can either load the mesh from the file as in Dolfin but # we can also use the msh pkg to generate the L-shape domain L-shape-domain; mesh = Mesh (msho); # Define function spaces (P2-P1). V = FunctionSpace ('VelocityUpdate', mesh); Q = FunctionSpace ('PressureUpdate', mesh); # Set parameter values and define coefficients dt = 0.01; T = 3.; k = Constant ('k', dt); f = Constant ('f', [0; 0]); u0 = Expression ('u0', @(x,y) [0; 0]); # Define boundary conditions noslip = DirichletBC (V, @(x,y) [0; 0], [3, 4]); outflow = DirichletBC (Q, @(x,y) 0, 2); # Assemble matrices a1 = BilinearForm ('TentativeVelocity', V, V, k); a2 = BilinearForm ('PressureUpdate', Q, Q); a3 = BilinearForm ('VelocityUpdate', V, V); A1 = assemble (a1, noslip); A3 = assemble (a3, noslip); # Time-stepping t = dt; i = 0; while t < T # Update pressure boundary condition inflow = DirichletBC (Q, @(x,y) sin(3.0*t), 1); # Compute tentative velocity step L1 = LinearForm ('TentativeVelocity', V, k, u0, f); b1 = assemble (L1, noslip); utmp = A1 \ b1; u1 = Function ('u1', V, utmp); # Pressure correction L2 = LinearForm ('PressureUpdate', Q, u1, k); [A2, b2] = assemble_system (a2, L2, inflow, outflow); ptmp = A2 \ b2; p1 = Function ('p1', Q, ptmp); # Velocity correction L3 = LinearForm ('VelocityUpdate', V, k, u1, p1); b3 = assemble (L3, noslip); ut = A3 \ b3; u1 = Function ('u0', V, ut); # Save to file save (p1, sprintf ("p_%3.3d", ++i)); save (u1, sprintf ("u_%3.3d", i)); # Move to next time step u0 = u1; t += dt end \end{lstlisting} \paragraph{L-shape-domain.m} \begin{lstlisting} name = [tmpnam ".geo"]; fid = fopen (name, "w"); fputs (fid,"Point (1) = {0, 0, 0, 0.1};\n"); fputs (fid,"Point (2) = {1, 0, 0, 0.1};\n"); fputs (fid,"Point (3) = {1, 0.5, 0, 0.1};\n"); fputs (fid,"Point (4) = {0.5, 0.5, 0, 0.1};\n"); fputs (fid,"Point (5) = {0.5, 1, 0, 0.1};\n"); fputs (fid,"Point (6) = {0, 1, 0,0.1};\n"); fputs (fid,"Line (1) = {5, 6};\n"); fputs (fid,"Line (2) = {2, 3};\n"); fputs (fid,"Line(3) = {6,1,2};\n"); fputs (fid,"Line(4) = {5,4,3};\n"); fputs (fid,"Line Loop(7) = {3,2,-4,1};\n"); fputs (fid,"Plane Surface(8) = {7};\n"); fclose (fid); msho = msh2m_gmsh (canonicalize_file_name (name)(1:end-4),... "scale", 1,"clscale", .2); unlink (canonicalize_file_name (name)); \end{lstlisting} \section{A penalization method to take into account obstacles in incompressible viscous flows} \newpage \bibliographystyle{unsrt} \bibliography{doc} \end{document}