// This is a simple standalone example. See README.txt #include #include #include #include #include "magma_v2.h" #include "magmasparse.h" // ------------------------------------------------------------ // This is an example how magma can be integrated into another software. int main( int argc, char** argv ) { if(argc<4){printf("please enter mat (mtx) rhs (dat/txt) and Nrow\n");exit(0);} int m=atoi(argv[3]); int i, n=1,counter=0; double *sol; // Initialize MAGMA and create some LA structures. magma_init(); magma_dopts opts; magma_queue_t queue; magma_queue_create( 0, &queue ); magma_d_matrix A={Magma_CSR}, dA={Magma_CSR}; magma_d_matrix b={Magma_CSR}, db={Magma_CSR}; magma_d_matrix x={Magma_CSR}, dx={Magma_CSR}; sol = (double*) calloc(m, sizeof(double)); magma_d_csr_mtx(&A,argv[1], queue ); magma_d_vread(&b,m, argv[2],queue ); magma_dvset( m, 1, sol, &x, queue ); // Choose a solver, preconditioner, etc. - see documentation for options. opts.solver_par.solver = Magma_PGMRES; opts.solver_par.restart = 30; opts.solver_par.maxiter = 1000; opts.solver_par.rtol = 1e-10; opts.solver_par.atol = 1e-16; opts.precond_par.solver = Magma_ILU; opts.precond_par.levels = 0; opts.precond_par.trisolver = Magma_CUSOLVE;// Magma_ISAI;// // //// // Initialize the solver. magma_dsolverinfo_init( &opts.solver_par, &opts.precond_par, queue ); // Copy the system to the device (optional, only necessary if using the GPU) magma_dmtransfer( A, &dA, Magma_CPU, Magma_DEV, queue ); magma_dmtransfer( b, &db, Magma_CPU, Magma_DEV, queue ); magma_dmtransfer( x, &dx, Magma_CPU, Magma_DEV, queue ); //first system magma_dcumilusetup(dA,&opts.precond_par, queue); printf("ILU1\n"); system("nvidia-smi|grep memoryCheck"); magma_dfgmres(dA,db,&dx,&opts.solver_par, &opts.precond_par,queue); printf("FGMRES1\n"); system("nvidia-smi|grep memoryCheck"); //second system magma_dcumilusetup(dA,&opts.precond_par, queue); printf("ILU2\n"); system("nvidia-smi|grep memoryCheck"); magma_dfgmres(dA,db,&dx,&opts.solver_par, &opts.precond_par,queue); printf("FGMRES2\n"); system("nvidia-smi|grep memoryCheck"); // Then copy the solution back to the host... magma_dmtransfer( dx, &x, Magma_DEV,Magma_CPU, queue ); // and back to the application code magma_dvget( x, &m, &n, &sol, queue ); // Free the allocated memory... magma_dmfree( &dx, queue ); magma_dmfree( &db, queue ); magma_dmfree( &dA, queue ); // and finalize MAGMA. magma_queue_destroy( queue ); magma_finalize(); printf("3\n"); system("nvidia-smi|grep memoryCheck"); // From here on, the application code may continue with the solution in sol... for (i = 0; i < 10; ++i) { printf("%e\n", sol[i]); } return 0; }