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 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215
  
     | 
    
      <!DOCTYPE html PUBLIC "-//W3C//DTD HTML 4.01 Transitional//EN" "http://www.w3.org/TR/html4/loose.dtd">
<html>
<!-- Copyright (C) 1996, 1997, 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009, 2010, 2011, 2012, 2013, 2014, 2015, 2016 The GSL Team.
Permission is granted to copy, distribute and/or modify this document
under the terms of the GNU Free Documentation License, Version 1.3 or
any later version published by the Free Software Foundation; with the
Invariant Sections being "GNU General Public License" and "Free Software
Needs Free Documentation", the Front-Cover text being "A GNU Manual",
and with the Back-Cover Text being (a) (see below). A copy of the
license is included in the section entitled "GNU Free Documentation
License".
(a) The Back-Cover Text is: "You have the freedom to copy and modify this
GNU Manual." -->
<!-- Created by GNU Texinfo 5.1, http://www.gnu.org/software/texinfo/ -->
<head>
<title>GNU Scientific Library – Reference Manual: Sparse Linear Algebra Examples</title>
<meta name="description" content="GNU Scientific Library – Reference Manual: Sparse Linear Algebra Examples">
<meta name="keywords" content="GNU Scientific Library – Reference Manual: Sparse Linear Algebra Examples">
<meta name="resource-type" content="document">
<meta name="distribution" content="global">
<meta name="Generator" content="makeinfo">
<meta http-equiv="Content-Type" content="text/html; charset=utf-8">
<link href="index.html#Top" rel="start" title="Top">
<link href="Function-Index.html#Function-Index" rel="index" title="Function Index">
<link href="Sparse-Linear-Algebra.html#Sparse-Linear-Algebra" rel="up" title="Sparse Linear Algebra">
<link href="Sparse-Linear-Algebra-References-and-Further-Reading.html#Sparse-Linear-Algebra-References-and-Further-Reading" rel="next" title="Sparse Linear Algebra References and Further Reading">
<link href="Iterating-the-Sparse-Linear-System.html#Iterating-the-Sparse-Linear-System" rel="previous" title="Iterating the Sparse Linear System">
<style type="text/css">
<!--
a.summary-letter {text-decoration: none}
blockquote.smallquotation {font-size: smaller}
div.display {margin-left: 3.2em}
div.example {margin-left: 3.2em}
div.indentedblock {margin-left: 3.2em}
div.lisp {margin-left: 3.2em}
div.smalldisplay {margin-left: 3.2em}
div.smallexample {margin-left: 3.2em}
div.smallindentedblock {margin-left: 3.2em; font-size: smaller}
div.smalllisp {margin-left: 3.2em}
kbd {font-style:oblique}
pre.display {font-family: inherit}
pre.format {font-family: inherit}
pre.menu-comment {font-family: serif}
pre.menu-preformatted {font-family: serif}
pre.smalldisplay {font-family: inherit; font-size: smaller}
pre.smallexample {font-size: smaller}
pre.smallformat {font-family: inherit; font-size: smaller}
pre.smalllisp {font-size: smaller}
span.nocodebreak {white-space:nowrap}
span.nolinebreak {white-space:nowrap}
span.roman {font-family:serif; font-weight:normal}
span.sansserif {font-family:sans-serif; font-weight:normal}
ul.no-bullet {list-style: none}
-->
</style>
</head>
<body lang="en" bgcolor="#FFFFFF" text="#000000" link="#0000FF" vlink="#800080" alink="#FF0000">
<a name="Sparse-Linear-Algebra-Examples"></a>
<div class="header">
<p>
Next: <a href="Sparse-Linear-Algebra-References-and-Further-Reading.html#Sparse-Linear-Algebra-References-and-Further-Reading" accesskey="n" rel="next">Sparse Linear Algebra References and Further Reading</a>, Previous: <a href="Sparse-Iterative-Solvers.html#Sparse-Iterative-Solvers" accesskey="p" rel="previous">Sparse Iterative Solvers</a>, Up: <a href="Sparse-Linear-Algebra.html#Sparse-Linear-Algebra" accesskey="u" rel="up">Sparse Linear Algebra</a>   [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
<hr>
<a name="Examples-32"></a>
<h3 class="section">43.3 Examples</h3>
<a name="index-sparse-linear-algebra_002c-examples"></a>
<p>This example program demonstrates the sparse linear algebra routines on
the solution of a simple 1D Poisson equation on <em>[0,1]</em>:
</p><div class="example">
<pre class="example">u''(x) = f(x) = -\pi^2 \sin(\pi x)
</pre></div>
<p>with boundary conditions <em>u(0) = u(1) = 0</em>. The analytic solution of
this simple problem is <em>u(x) = \sin{\pi x}</em>. We will solve this
problem by finite differencing the left hand side to give
</p><div class="example">
<pre class="example">1/h^2 ( u_(i+1) - 2 u_i + u_(i-1) ) = f_i
</pre></div>
<p>Defining a grid of <em>N</em> points, <em>h = 1/(N-1)</em>. In the finite
difference equation above, <em>u_0 = u_{N-1} = 0</em> are known from
the boundary conditions, so we will only put the equations for
<em>i = 1, ..., N-2</em> into the matrix system. The resulting
<em>(N-2) \times (N-2)</em> matrix equation is
An example program which constructs and solves this system is given below.
The system is solved using the iterative GMRES solver. Here is
the output of the program:
</p>
<div class="example">
<pre class="example">iter 0 residual = 4.297275996844e-11
Converged
</pre></div>
<p>showing that the method converged in a single iteration.
The calculated solution is shown in the following plot.
</p>
<div class="example">
<pre class="verbatim">#include <stdio.h>
#include <stdlib.h>
#include <math.h>
#include <gsl/gsl_math.h>
#include <gsl/gsl_vector.h>
#include <gsl/gsl_spmatrix.h>
#include <gsl/gsl_splinalg.h>
int
main()
{
  const size_t N = 100;                       /* number of grid points */
  const size_t n = N - 2;                     /* subtract 2 to exclude boundaries */
  const double h = 1.0 / (N - 1.0);           /* grid spacing */
  gsl_spmatrix *A = gsl_spmatrix_alloc(n ,n); /* triplet format */
  gsl_spmatrix *C;                            /* compressed format */
  gsl_vector *f = gsl_vector_alloc(n);        /* right hand side vector */
  gsl_vector *u = gsl_vector_alloc(n);        /* solution vector */
  size_t i;
  /* construct the sparse matrix for the finite difference equation */
  /* construct first row */
  gsl_spmatrix_set(A, 0, 0, -2.0);
  gsl_spmatrix_set(A, 0, 1, 1.0);
  /* construct rows [1:n-2] */
  for (i = 1; i < n - 1; ++i)
    {
      gsl_spmatrix_set(A, i, i + 1, 1.0);
      gsl_spmatrix_set(A, i, i, -2.0);
      gsl_spmatrix_set(A, i, i - 1, 1.0);
    }
  /* construct last row */
  gsl_spmatrix_set(A, n - 1, n - 1, -2.0);
  gsl_spmatrix_set(A, n - 1, n - 2, 1.0);
  /* scale by h^2 */
  gsl_spmatrix_scale(A, 1.0 / (h * h));
  /* construct right hand side vector */
  for (i = 0; i < n; ++i)
    {
      double xi = (i + 1) * h;
      double fi = -M_PI * M_PI * sin(M_PI * xi);
      gsl_vector_set(f, i, fi);
    }
  /* convert to compressed column format */
  C = gsl_spmatrix_ccs(A);
  /* now solve the system with the GMRES iterative solver */
  {
    const double tol = 1.0e-6;  /* solution relative tolerance */
    const size_t max_iter = 10; /* maximum iterations */
    const gsl_splinalg_itersolve_type *T = gsl_splinalg_itersolve_gmres;
    gsl_splinalg_itersolve *work =
      gsl_splinalg_itersolve_alloc(T, n, 0);
    size_t iter = 0;
    double residual;
    int status;
    /* initial guess u = 0 */
    gsl_vector_set_zero(u);
    /* solve the system A u = f */
    do
      {
        status = gsl_splinalg_itersolve_iterate(C, f, tol, u, work);
        /* print out residual norm ||A*u - f|| */
        residual = gsl_splinalg_itersolve_normr(work);
        fprintf(stderr, "iter %zu residual = %.12e\n", iter, residual);
        if (status == GSL_SUCCESS)
          fprintf(stderr, "Converged\n");
      }
    while (status == GSL_CONTINUE && ++iter < max_iter);
    /* output solution */
    for (i = 0; i < n; ++i)
      {
        double xi = (i + 1) * h;
        double u_exact = sin(M_PI * xi);
        double u_gsl = gsl_vector_get(u, i);
        printf("%f %.12e %.12e\n", xi, u_gsl, u_exact);
      }
    gsl_splinalg_itersolve_free(work);
  }
  gsl_spmatrix_free(A);
  gsl_spmatrix_free(C);
  gsl_vector_free(f);
  gsl_vector_free(u);
  return 0;
} /* main() */
</pre></div>
<hr>
<div class="header">
<p>
Next: <a href="Sparse-Linear-Algebra-References-and-Further-Reading.html#Sparse-Linear-Algebra-References-and-Further-Reading" accesskey="n" rel="next">Sparse Linear Algebra References and Further Reading</a>, Previous: <a href="Sparse-Iterative-Solvers.html#Sparse-Iterative-Solvers" accesskey="p" rel="previous">Sparse Iterative Solvers</a>, Up: <a href="Sparse-Linear-Algebra.html#Sparse-Linear-Algebra" accesskey="u" rel="up">Sparse Linear Algebra</a>   [<a href="Function-Index.html#Function-Index" title="Index" rel="index">Index</a>]</p>
</div>
</body>
</html>
 
     |