SpM Handbook 1.2.4
Loading...
Searching...
No Matches
spm_driver.jl
1#!/usr/bin/env julia
2#=
3 @file spm_driver.jl
4
5 @brief SpM example to generate a sparse matrix from the spm drivers
6
7 @copyright 2019-2024 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
8 Univ. Bordeaux. All rights reserved.
9
10 @version 1.2.4
11 @author Mathieu Faverge
12 @author Selmane Lebdaoui
13 @author Tony Delarue
14 @date 2024-05-29
15
16 @ingroup examples_julia
17 @code
18
19=#
20using Pkg
21Pkg.activate("spm")
22Pkg.instantiate()
23using CBinding
24using spm
25
26if spm.spm_mpi_enabled
27 using MPI
28 MPI.Init()
29end
30
31#
32# Two solutions to select the outpu file to pass to output functions
33#
34my_stdout = Libc.FILE( Libc.RawFD(1), "w" ) # Select the stdout or stderr through 1 or 2
35#my_stdout = Ptr{Cvoid}(0) # Give a null pointer to use the default
36
37A = spm.spmatrix_t(zero)
38Aptr = Ptr{spm.spmatrix_t}(pointer_from_objref(A))
39
40# Example from a HB file
41#spm.spmReadDriver( spm.SpmDriverHB, "__SPM_DIR__/tests/matrix/orsirr.rua", Aptr )
42
43# Example from the Laplacian generator driver
44spm.spmReadDriver( spm.SpmDriverLaplacian, "10:10:10:2.:1.", Aptr )
45
46spm.spmPrintInfo( Aptr, my_stdout )
47
48# Scale A for low-rank: A / ||A||_f
49norm = spm.spmNorm( spm.SpmFrobeniusNorm, Aptr )
50spm.spmScal( 1. / norm, Aptr )
51
52# Generate b and x0 vectors such that A * x0 = b
53
54nrhs = 10
55n = A.nexp
56X0 = zeros( Cdouble, (n, nrhs) )
57B = zeros( Cdouble, (n, nrhs) )
58X = zeros( Cdouble, (n, nrhs) )
59
60spm.spmGenRHS( spm.SpmRhsRndX, nrhs, Aptr, X, n, B, n )
61
62# Copy x0 into x for backup
63X0 = copy(X)
64
65# Check that A * x = b
66eps = 1.e-15 # Set to 1e-7 for single precision
67spm.spmCheckAxb( eps, nrhs, Aptr, X0, n, B, n, X, n )
68
69spm.spmExit( Aptr )
70
71#=
72 @endcode
73=#