File: example_sparse.cpp

package info (click to toggle)
magma-rocm 2.9.0%2Bds-3
  • links: PTS, VCS
  • area: main
  • in suites: sid
  • size: 83,540 kB
  • sloc: cpp: 709,115; fortran: 121,916; ansic: 32,343; python: 25,603; f90: 15,208; makefile: 945; xml: 253; csh: 232; sh: 203; perl: 104
file content (112 lines) | stat: -rw-r--r-- 3,704 bytes parent folder | download | duplicates (6)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
// This is a simple standalone example. See README.txt

#include <stdlib.h>
#include <stdio.h>
#include <string.h>
#include <math.h>

#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 )
{
    // The software does e.g. discretization of a PDE,
    // ends up with a sparse linear system in CSR format and a RHS.
    // Let's assume this system is a diagonal system of size m.
    
    int i, m=200, n=1;
    int *col, *row;
    double *val, *rhs, *sol;
    
    row = (int*) calloc(m+1, sizeof(int));
    col = (int*) calloc(m,   sizeof(int));
    
    val = (double*) calloc(m, sizeof(double));
    rhs = (double*) calloc(m, sizeof(double));
    sol = (double*) calloc(m, sizeof(double));
    
    for (i = 0; i < m; ++i) {
        col[i] = i;
        row[i] = i;
        val[i] = 55.0;
        rhs[i] = 3.0;
        sol[i] = 0.0;
    }
    row[m] = m;
    
    // Initialize MAGMA and create some LA structures.
    magma_init();
    magma_dopts opts;
    magma_queue_t queue;
    magma_queue_create( 0, &queue );
    
    memset(&opts, 0, sizeof(magma_dopts));
    
    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};
    
    // Pass the system to MAGMA.
    magma_dcsrset( m, m, row, col, val, &A, queue );
    magma_dvset( m, 1, rhs, &b, queue );
    
    // Choose a solver, preconditioner, etc. - see documentation for options.
    opts.solver_par.solver     = Magma_PIDRMERGE;
    opts.solver_par.restart    = 8;
    opts.solver_par.maxiter    = 1000;
    opts.solver_par.rtol       = 1e-10;
    opts.solver_par.maxiter    = 1000;
    opts.precond_par.solver    = Magma_ILU;
    opts.precond_par.levels    = 0;
    opts.precond_par.trisolver = Magma_CUSOLVE;
    
    // Initialize the solver.
    magma_dsolverinfo_init( &opts.solver_par, &opts.precond_par, queue );
    
    // Copy the system to the device
    magma_dmtransfer( A, &dA, Magma_CPU, Magma_DEV, queue );
    magma_dmtransfer( b, &db, Magma_CPU, Magma_DEV, queue );
    // initialize an initial guess for the iteration vector
    magma_dvinit( &dx, Magma_DEV, b.num_rows, b.num_cols, 0.0, queue );

    // Generate the preconditioner.
    magma_d_precondsetup( dA, db, &opts.solver_par, &opts.precond_par, queue );

    // In case we only wanted to generate a preconditioner, we are done.
    // The preconditioner in the opts.precond_par structure - in this case an ILU.
    // The lower ILU(0) factor is in opts.precond_par.L and
    // the upper ILU(0) factor is in opts.precond_par.U (in this case on the device).
    
    // If we want to solve the problem, run:
    magma_d_solver( dA, db, &dx, &opts, queue );
    
    // Then copy the solution back to the host...
    magma_dmtransfer( dx, &x, Magma_DEV, Magma_CPU, queue );
    // and back to the application code
    magma_dvcopy( x, &m, &n, sol, queue );
    
    // Free the allocated memory...
    magma_dmfree( &dx, queue );
    magma_dmfree( &db, queue ); 
    magma_dmfree( &dA, queue );
    magma_dmfree( &b, queue );  // won't do anything as MAGMA does not own the data. 
    magma_dmfree( &A, queue );  // won't do anything as MAGMA does not own the data. 
    // and finalize MAGMA.
    magma_queue_destroy( queue );
    magma_finalize();
    
    // From here on, the application code may continue with the solution in sol...
    for (i = 0; i < 20; ++i) {
        printf("%.4f\n", sol[i]);
    }
    free(val);
    free(col);
    free(row);
    free(sol);
    free(rhs);
    
    return 0;
}