Intro
We consider the Laplacian problem on a 2D square mesh of size \(N_X \times N_Y\). Unknowns are related to the mesh nodes indexed by \((i,j)\). We use a 5-Points stencil to discretize the problem.
First of all Alien provides tools to initialize the library, the MPI parallel environment to parametrize IO. The following code illustrates how to initialize :
- the ParallelMng object
- the TraceMng object
Environment::initialize(argc, argv);
auto parallel_mng = Environment::parallelMng();
auto trace_mng = Environment::traceMng();
auto comm_size = Environment::parallelMng()->commSize();
auto comm_rank = Environment::parallelMng()->commRank();
trace_mng->info() << "NB PROC = " << comm_size;
trace_mng->info() << "RANK = " << comm_rank;
Arccore::StringBuilder filename("tutorial.log");
Arccore::ReferenceCounter<Arccore::ITraceStream> ofile;
if (comm_size > 1) {
filename += comm_rank;
ofile = Arccore::ITraceStream::createFileStream(filename.toString());
trace_mng->setRedirectStream(ofile.get());
}
trace_mng->finishInitialize();
Alien::setTraceMng(trace_mng);
Alien::setVerbosityLevel(Alien::Verbosity::Debug);
Space
The Space concept enable to modelize the mathematical algebraic real space \(R^N\) of dimension \(N\)
To build this concept several tools are provided:
- the IndexeManager package provides helper tools to manage Integer IndexSets
- the Distribution package provides helper tools to manage the partition of IndexSets between MPI processes
int Nx = 10;
int Ny = 10;
int local_ny = Ny / comm_size;
int r = Ny % comm_size;
std::vector<int> y_offset(comm_size + 1);
y_offset[0] = 0;
for (int ip = 0; ip < r; ++ip)
y_offset[ip + 1] = y_offset[ip] + local_ny + 1;
for (int ip = r; ip < comm_size; ++ip)
y_offset[ip + 1] = y_offset[ip] + local_ny;
auto node_uid = [&](int i, int j) { return j * Nx + i; };
Alien::UniqueArray<UID> uid;
Alien::UniqueArray<LID> lid;
int first_j = y_offset[comm_rank];
int last_j = y_offset[comm_rank + 1];
int index = 0;
for (int j = first_j; j < last_j; ++j) {
for (int i = 0; i < Nx; ++i) {
uid.add(node_uid(i, j));
lid.add(index);
++index;
}
}
auto indexSetU = index_manager.buildScalarIndexSet("U", lid, family, 0);
index_manager.prepare();
auto global_size = index_manager.globalSize();
auto local_size = index_manager.localSize();
trace_mng->info() << "GLOBAL SIZE : " << global_size;
trace_mng->info() << "LOCAL SIZE : " << local_size;
auto mdist =
trace_mng->info() << "MATRIX DISTRIBUTION INFO";
trace_mng->info() << "GLOBAL ROW SIZE : " << mdist.globalRowSize();
trace_mng->info() << "LOCAL ROW SIZE : " << mdist.localRowSize();
trace_mng->info() << "GLOBAL COL SIZE : " << mdist.globalColSize();
trace_mng->info() << "LOCAL COL SIZE : " << mdist.localColSize();
trace_mng->info() << "VECTOR DISTRIBUTION INFO";
trace_mng->info() << "GLOBAL SIZE : " << vdist.globalSize();
trace_mng->info() << "LOCAL SIZE : " << vdist.localSize();
Computes a matrix distribution.
Implementation of an algebraic space.
Computes a vector distribution.
Matrix
The Matrix concept represents a set of \(N_i\) linear equations (rows) \((y_i)\) of \(N_j\) unknowns \((x_j)\) (columns). This represents a linear application \(S_X \mapsto S_Y\) with \(x \in S_X\), \(y \in S_Y\) and \(x \mapsto y=A*x\). Usually the dimension of \(S_X\) and \(S_Y\) are equal, \(N_i=N_j\). In that case the matrix is square.
alien_info([&] { cout() << "DIRECT ONE STEP FILLING PHASE";}) ;
auto tag = Alien::DirectMatrixOptions::eResetValues;
{
builder.reserve(5);
builder.allocate();
for (int j = first_j; j < last_j; ++j) {
for (int i = 0; i < Nx; ++i) {
auto n_uid = node_uid(i, j);
auto n_lid = uid2lid[n_uid];
auto irow = allUIndex[n_lid];
builder(irow, irow) = 4;
if (j > 0) {
auto off_uid = node_uid(i, j - 1);
auto off_lid = uid2lid[off_uid];
auto jcol = allUIndex[off_lid];
if (jcol != -1)
builder(irow, jcol) = -1;
}
if (i > 0) {
auto off_uid = node_uid(i - 1, j);
auto off_lid = uid2lid[off_uid];
auto jcol = allUIndex[off_lid];
if (jcol != -1)
builder(irow, jcol) = -1;
}
if (i < Nx - 1) {
auto off_uid = node_uid(i + 1, j);
auto off_lid = uid2lid[off_uid];
auto jcol = allUIndex[off_lid];
if (jcol != -1)
builder(irow, jcol) = -1;
}
if (j < Ny - 1) {
auto off_uid = node_uid(i, j + 1);
auto off_lid = uid2lid[off_uid];
auto jcol = allUIndex[off_lid];
if (jcol != -1)
builder(irow, jcol) = -1;
}
}
}
}
Vector
The Vector concept represents the set of unknowns \(x=(x_i)\) element of a linear space \(S_X\).
{
for (int j = first_j; j < last_j; ++j) {
for (int i = 0; i < Nx; ++i) {
auto n_uid = node_uid(i, j);
auto n_lid = uid2lid[n_uid];
auto irow = allUIndex[n_lid];
writer[irow] = 1. / (1. + i + j);
}
}
}
{
for (int i = 0; i < reader.size(); ++i) {
trace_mng->info() << "B[" << i << "]=" << reader[i];
}
}
Linear Systems resolution
A linear system is represented by a matrix \(A\), and two vectors \(B\) and \(X\) where \(B\) is the system right hand side and \(X\) the solution.
Solving the linear system consists in finding the solution X such that \(A*X=B\) applying a linear solver algorithm.
auto solver = createSolver() ;
solver->init() ;
solver->solve(matrixA, vectorB, vectorX);
{
alien_info()([&]{ cout()<<"SOLVER HAS SUCCEEDED";}) ;
SimpleCSRLinearAlgebra alg;
alg.mult(matrixA, vectorX, vectorR);
alg.axpy(-1., vectorB, vectorR);
Real res = alg.norm2(vectorR);
alien_info([&] cout() << "RES : " << res;}) ;
}
else
alien_info()([&]{ cout()<<"SOLVER FAILED";}) ;
solver->getSolverStat().print(Universe().traceMng(), status, "Linear Solver : ");
}
solver->end();
Structure to store a solver status.
bool succeeded
Whether or not the solver succeeded.