The \akantu source code can be requested using the form accessible at the URL
\url{http://lsms.epfl.ch/akantu}. There, you will be asked to accept the LGPL
license terms.
\section{Compiling \akantu}
\akantu is a \code{cmake} project, so to configure it, you can either
follow the usual way:
\begin{command}
> cd akantu
> mkdir build
> cd build
> ccmake ..
[ Set the options that you need ]
> make
> make install
\end{command}
\noindent Or, use the \code{Makefile} we added for your convenience to
handle the \code{cmake} configuration
\begin{command}
> cd akantu
> make config
> make
> make install
\end{command}
\noindent All the \akantu options are documented in Appendix
\ref{app:package-dependencies}.
\section{Writing a \texttt{main} Function\label{sect:common:main}}
\label{sec:writing_main}
First of all, \akantu needs to be initialized. The memory management
included in the core library handles the correct allocation and
de-allocation of vectors, structures and/or objects. Moreover, in
parallel computations, the initialization procedure performs the
communication setup. This is achieved by a pair of functions
(\code{initialize} and \code{finalize}) that are used as follows:
\begin{cpp}
#include "aka_common.hh"
#include "..."
using namespace akantu;
int main(int argc, char *argv[]) {
initialize("input_file.dat", argc, argv);
// your code
...
finalize();
}
\end{cpp}
The \code{initialize} function takes the text inpute file and the program
parameters which can be parsed by \akantu in due form (see \ref{sect:parser}).
Obviously it is necessary to include all files needed in main. In this manual
all provided code implies the usage of \code{akantu} as namespace.
\section{Creating and Loading a Mesh\label{sect:common:mesh}}
In its current state, \akantu supports three types of meshes: Gmsh~\cite{gmsh},
Abaqus~\cite{abaqus} and Diana~\cite{diana}. Once a \code{Mesh} object is
created with a given spatial dimension, it can be filled by reading a mesh input file.
The method \code{read} of the class \code{Mesh} infers the mesh type from the file extension. If a non-standard file extension is used, the mesh type has to be specified.
\begin{cpp}
UInt spatial_dimension = 2;
Mesh mesh(spatial_dimension);
// Reading Gmsh files
mesh.read("my_gmsh_mesh.msh");
mesh.read("my_gmsh_mesh", _miot_gmsh);
// Reading Abaqus files
mesh.read("my_abaqus_mesh.inp");
mesh.read("my_abaqus_mesh", _miot_abaqus);
// Reading Diana files
mesh.read("my_diana_mesh.dat");
mesh.read("my_diana_mesh", _miot_diana);
\end{cpp}
The Gmsh reader adds the geometrical and physical tags as mesh data. The physical
values are stored as a \code{UInt} data called \code{tag\_0}, if a string
name is provided it is stored as a \code{std::string} data named
\code{physical\_names}. The geometrical tag is stored as a \code{UInt} data named
\code{tag\_1}.
The Abaqus reader stores the \code{ELSET} in ElementGroups and the \code{NSET}
in NodeGroups. The material assignment can be retrieved from the
\code{std::string} mesh data named \code{abaqus\_material}.
% \akantu supports meshes generated with Gmsh~\cite{gmsh}, a free
% software available at \url{http://geuz.org/gmsh/} where a detailed
% documentation can be found. Consequently, this manual will not provide
% Gmsh usage directions. Gmsh outputs meshes in \code{.msh} format that can be read
% by \akantu. In order to import a mesh, it is necessary to create
% a \code{Mesh} object through the following function calls:
% \begin{cpp}
% UInt spatial_dimension = 2;
% Mesh mesh(spatial_dimension);
% \end{cpp}
% The only parameter that has to be specified by the user is the spatial
% dimension of the problem. Now it is possible to read a \code{.msh} file with
% a \code{MeshIOMSH} object that takes care of loading a mesh to memory.
% This step is carried out by:
% \begin{cpp}
% mesh.read("square.msh");
% \end{cpp}
% where the \code{MeshIOMSH} object is first created before being
% used to read the \code{.msh} file. The mesh file name as well as the \code{Mesh}
% object must be specified by the user.
% The \code{MeshIOMSH} object can also write mesh files. This feature
% is useful to save a mesh that has been modified during a
% simulation. The \code{write} method takes care of it:
% \begin{cpp}
% mesh_io.write("square_modified.msh", mesh);
% \end{cpp}
% which works exactly like the \code{read} method.
% \akantu supports also meshes generated by
% DIANA (\url{http://tnodiana.com}), but only in reading mode. A similar
% procedure applies where the only
% difference is that the \code{MeshIODiana} object should be used
% instead of the \code{MeshIOMSH} one. Additional mesh readers can be
% introduced into \akantu by coding new \code{MeshIO} classes.
\section{Using \texttt{Arrays}}
Data in \akantu can be stored in data containers implemented by
the \code{Array} class. In its most basic usage, the \code{Array} class
implemented in \akantu is similar to the \code{vector} class of
the Standard Template Library (STL) for C++. A simple \code{Array}
containing a sequence of \code{nb\_element} values (of a given type) can be generated with:
\begin{cpp}
Array<type> example_array(nb_element);
\end{cpp}
where \code{type} usually is \code{Real}, \code{Int}, \code{UInt} or
\code{bool}. Each value is associated to an index, so that data can be
accessed by typing:
\begin{cpp}
auto & val = example_array(index)
\end{cpp}
\code{Arrays} can also contain tuples of values for each index. In
that case, the number of components per tuple must be specified at the
\code{Array} creation. For example, if we want to create an
\code{Array} to store the coordinates (sequences of three values) of
In this case the $x$ position of the eighth node number will be given by
\code{position(7, 0)} (in C++, numbering starts at 0 and not
1). If the number of components for the sequences is not specified, the
default value of 1 is used.
Here is a list of some basic operations that can be performed on \code{Array}:
\begin{itemize}
\item \code{resize(size)} change the size of the \code{Array}.
\item \code{clear()} set all entries of the \code{Array} to zero.
\item \code{set(t)} set all entries of the \code{Array} to \code{t}.
\item \code{copy(const Array<T> $\&$ other)} copy another \code{Array} into the current one. The two \code{Array} should have the same number of components.
\item \code{push$\_$back(tuple)} append a tuple with the correct number of components at the end of the \code{Array}.
\item \code{erase(i) erase the value at the i-th position.}
\item \code{find(value)} search \code{value} in the current \code{Array}. Return position index of the first occurence or $-1$ if not found.
\item \code{storage()} Return the address of the allocated memory of the \code{Array}.
\end{itemize}
\subsection{\texttt{Arrays} iterators}
It is very common in \akantu to loop over arrays to perform a specific
treatment. This ranges from geometric calculation on nodal quantities
to tensor algebra (in constitutive laws for example).
The \code{Array} object has the possibility to request iterators
in order to make the writing of loops easier and enhance readability.
For instance, a loop over the nodal coordinates can be performed like:
\begin{cpp}
//accessing the nodal coordinates Array (spatial_dimension components)
const auto & nodes = mesh.getNodes();
//creating the iterators
auto it = nodes.begin(spatial_dimension);
auto end = nodes.end(spatial_dimension);
for (; it != end; ++it){
const auto & coords = (*it);
//do what you need
....
}
\end{cpp}
In that example, each \code{coords} is a \code{Vector<Real>} containing
geometrical array of size \code{spatial\_dimension} and the iteration is
conveniently performed by the \code{Array} iterator.
With the switch to \code{c++14} this can be also written as
\begin{cpp}
//accessing the nodal coordinates Array (spatial_dimension components)
const auto & nodes = mesh.getNodes();
for (const auto & coords : make_view(nodes, spatial_dimension) {
//do what you need
....
}
\end{cpp}
The \code{Array} object is intensively used to store second order
tensor values. In that case, it should be specified that the returned
object type is a matrix when constructing the iterator. This is done
when calling the \code{begin} function. For instance, assuming that we
have a \code{Array} storing stresses, we can loop over the stored
tensors by:
\begin{cpp}
//creating the iterators
auto it = stresses.begin(spatial_dimension, spatial_dimension);
auto end = stresses.end(spatial_dimension, spatial_dimension);
for (; it != end; ++it){
Matrix<Real> & stress = (*it);
//do what you need
....
}
\end{cpp}
In that last example, the \code{Matrix} objects are
A group of nodes is stored in \code{NodeGroup} objects.
They are quite simple objects which store the indexes
of the selected nodes in a \code{Array<UInt>}.
Nodes are selected by adding them when calling
\code{NodeGroup::add}. For instance you can select nodes
having a positive $X$ coordinate with the following code:
\begin{cpp}
const auto & nodes = mesh.getNodes();
auto & group = mesh.createNodeGroup("XpositiveNode");
for (auto && data : enumerate(make_view(nodes, spatial_dimension))){
auto node = std::get<0>(data);
const auto & position = std::get<1>(data);
if (position(0) > 0) group.add(node);
}
\end{cpp}
\subsection{The \texttt{ElementGroup} object}
A group of elements is stored in \code{ElementGroup} objects.
Since a group can contain elements of various types
the \code{ElementGroup} object stores indexes in
a \code{ElementTypeMapArray<UInt>} object.
Then elements can be added to the group by calling \code{addElement}.
For instance, selecting the elements for which the barycenter of the nodes
has a positive $X$ coordinate can be made with:
\begin{cpp}
auto & group = mesh.createElementGroup("XpositiveElement");
Vector<Real> barycenter(spatial_dimension);
for(auto type : mesh.elementTypes()){
UInt nb_element = mesh.getNbElement(type);
for(UInt e = 0; e < nb_element; ++e) {
Element element{type, e, _not_ghost};
mesh.getBarycenter(element, barycenter);
if (barycenter(_x) > 0.) group.add(element);
}
}
\end{cpp}
\section{Compiling your simulation}
The easiest way to compile your simulation is to create a \code{cmake} project by putting all your code in some directory of your choosing. Then, make sure that you have \code{cmake} installed and create a \code{CMakeLists.txt} file. An example of a minimal \code{CMakeLists.txt} file would look like this:
\begin{cmake}
project(my_simu)
cmake_minimum_required(VERSION 3.0.0)
find_package(Akantu REQUIRED)
add_akantu_simulation(my_simu my_simu.cc)
\end{cmake}
%
Then create a directory called \code{build} and inside it execute \code{ccmake ..} which opens a configuration screen. If you installed \akantu in a standard directory such as \code{/usr/local} (using \code{make install}), you should be able to compile by hitting the \code{c} key, setting \code{CMAKE\textunderscore{}BUILD\textunderscore{}TYPE} to \code{Release}, hitting the \code{c} key again a few times and then finishing the configuration with the \code{g} key. After that, you can type \code{make} to build your simulation. If you change your simulation code later on, you only need to type \code{make} again.
If you get an error that \code{FindAkantu.cmake} was not found, you have to set the \code{Akantu\textunderscore{}DIR} variable, which will appear after dismissing the error message. If you built \akantu without running \code{make install}, the variable should be set to the \code{build} subdirectory of the \akantu source code. If you installed it in \code{\$PREFIX}, set the variable to \code{\$PREFIX/share/cmake/Akantu} instead.