We created interfaces in 4 different programming languages, you can click for know how use it
SMG2S install command will generate a shared library libsmg2s.so into ${INSTALL_DIRECTORY}/lib
. It can be used to profit the C wrapper of SMG2S. The compile command for C codes using SMG2S :
Include header file:
#include <interface/C/c_wrapper.h>
create Nilpotency object :
struct NilpotencyInt *n;
n = newNilpotencyInt();
NilpType1(n, 2, 10);
After that, you need to create the parallel Sparse Matrix Object Mt like this :
struct parMatrixSparseDoubleInt *m;
m = newParMatrixSparseDoubleInt();
Generate by SMG2S :
smg2s(m, 10, n, 3 ," ",MPI_COMM_WORLD);
GRelease Nilpotency Object and parMatrixSparse Object :
ReleaseNilpotencyInt(&n);
ReleaseParMatrixSparseDoubleInt(&m);
Generate the shared library and install the python module of smg2s
cd ./interface/python;
mpicxx -fpic -c smg2s_wrap.cxx -I/apps/python/3/include/python3.5m -std=c++0x
mpicxx -shared smg2s_wrap.o -o _smg2s.so
python setup.py install
Before the utilisation, make sure that mpi4py installed.
A little example of usge :
from mpi4py import MPI
import smg2s
import sys
size = MPI.COMM_WORLD.Get_size()
rank = MPI.COMM_WORLD.Get_rank()
name = MPI.Get_processor_name()
sys.stdout.write(
"Hello, World! I am process %d of %d on %s.\n"
% (rank, size, name))
if rank == 0:
print ('INFO ]> Starting ...')
print("INFO ]> The MPI Comm World Size is %d" %size)
#bandwidth for the lower band of initial matrix
lbandwidth = 3
#create the nilpotent matrix
nilp = smg2s.NilpotencyInt()
#setup the nilpotent matrix: 2 = continous 1 nb, 10 = matrix size
nilp.NilpType1(2,10)
if rank == 0:
print("Nilptency matrix continuous one nb = %d" %nilp.nbOne)
Mt = smg2s.parMatrixSparseDoubleInt()
#Generate Mt by SMG2S
#vector.txt is the file that stores the given spectral distribution in local filesystem.
Mt=smg2s.smg2sDoubleInt(10,nilp,lbandwidth,"vector.txt", MPI.COMM_WORLD)
Include header file
#include <interface/PETSc/petsc_interface.h>
Create parMatrixSparse type matrix :
parMatrixSparse<std::complex<double>,int> *Mt;
Restore this matrix into CSR format :
Mt->Loc_ConvertToCSR();
Create PETSc MAT type :
MatCreate(PETSC_COMM_WORLD,&A);
Convert to PETSc MAT format :
A = ConvertToPETSCMat(Mt);
Include header file
#include <interface/Trilinos/trilinos_interface.hpp>
Create parMatrixSparse type matrix :
parMatrixSparse<std::complex<double>,int> *Mt;
Create Trilinos/Teptra MAT type :
RCP<CrsMatrix<ScalarType> > K ;
Convert to Trilinos MAT format :
K = ConvertToTrilinosMat(Mt);