<?xml-stylesheet type="text/xsl" href="../styles/pmathml.xsl"?>
<!-- saved from url=(0014)about:internet -->
<html xmlns="http://www.w3.org/1999/xhtml" xmlns:dsi="http://www.w3.org/1999/xlink" xmlns:m="http://www.w3.org/1998/Math/MathML" xml:space="preserve"><head><meta http-equiv="Content-Type" content="text/html; charset=US-ASCII"/><title>F08HEF (DSBTRD) : NAG Library, Mark 23</title><link rel="stylesheet" href="../styles/libdoc.css" type="text/css"/><script type="text/javascript">
   function showLevel(_levelId){
    var thisLevel = document.getElementById(_levelId);
    var thisplus = document.getElementById( _levelId.concat('plus'));
    var thisminus = document.getElementById( _levelId.concat('minus'));
    if(thisLevel.style.display != "block"){
     thisLevel.style.display = "block";
     thisplus.style.display = "none";
     thisminus.style.display = "inline";
     }
    else{
     thisLevel.style.display = "none";
     thisminus.style.display = "none";
     thisplus.style.display = "inline";
     }
    }
  </script></head><body><hr/><div><a class="rout" href="../../pdf/F08/f08hef.pdf">F08HEF (DSBTRD) (PDF version)</a></div><div><a class="chap" href="f08conts.xml">F08 Chapter Contents</a></div><div><a class="chapint" href="f08intro.xml">F08 Chapter Introduction</a></div>
<div><a class="htmltoc" href="../FRONTMATTER/manconts.xml">NAG Library Manual</a></div><hr/><h1 class="libdoc">NAG Library Routine Document<br/><br/>F08HEF (DSBTRD)</h1><div class="paramtext"><div class="header"><b>Note:</b>&#160; before using this routine, please read the Users' Note for your implementation to check the interpretation of <span class="bitalic">bold italicised</span> terms and other implementation-dependent details.</div></div> 
<div class="htmltoc">
<h2 class="htmltoc"><span class="htmltochead" onclick="showLevel('htmltoc');"><span class="htmltocplus" id="htmltocplus">+</span><span class="htmltocminus" id="htmltocminus">&#8722;</span></span>&#160;Contents</h2>
<div class="htmltocitem" id="htmltoc">
<div class="htmltoc">
<span class="htmltocplus">&#160;&#160;&#160;</span>
<a class="htmltoc" href="#purpose">1&#160;&#160;<b>Purpose</b></a>
</div><div class="htmltoc">
<span class="htmltocplus">&#160;&#160;&#160;</span>
<a class="htmltoc" href="#specification">2&#160;&#160;<b>Specification</b></a>
</div><div class="htmltoc">
<span class="htmltocplus">&#160;&#160;&#160;</span>
<a class="htmltoc" href="#description">3&#160;&#160;<b>Description</b></a>
</div><div class="htmltoc">
<span class="htmltocplus">&#160;&#160;&#160;</span>
<a class="htmltoc" href="#references">4&#160;&#160;<b>References</b></a>
</div><div class="htmltoc">
<span class="htmltocplus">&#160;&#160;&#160;</span>
<a class="htmltoc" href="#parameters">5&#160;&#160;<b>Parameters</b></a>
</div><div class="htmltoc">
<span class="htmltocplus">&#160;&#160;&#160;</span>
<a class="htmltoc" href="#errors">6&#160;&#160;<b>Error Indicators and Warnings</b></a>
</div><div class="htmltoc">
<span class="htmltocplus">&#160;&#160;&#160;</span>
<a class="htmltoc" href="#accuracy">7&#160;&#160;<b>Accuracy</b></a>
</div><div class="htmltoc">
<span class="htmltocplus">&#160;&#160;&#160;</span>
<a class="htmltoc" href="#fcomments">8&#160;&#160;<b>Further Comments</b></a>
</div><div class="htmltoc">
<span class="htmltoc" onclick="showLevel('tocexample');"><span class="htmltocplus" id="tocexampleplus">+</span><span class="htmltocminus" id="tocexampleminus">&#8722;</span></span>
<a class="htmltoc" href="#example">9&#160;&#160;<b>Example</b></a>
<div class="htmltocitem" id="tocexample">
<div class="htmltoc">
<span class="htmltocplus">&#160;&#160;&#160;</span>
<a class="htmltoc" href="#examtext">9.1&#160;&#160;<b>Program Text</b></a>
</div><div class="htmltoc">
<span class="htmltocplus">&#160;&#160;&#160;</span>
<a class="htmltoc" href="#examdata">9.2&#160;&#160;<b>Program Data</b></a>
</div><div class="htmltoc">
<span class="htmltocplus">&#160;&#160;&#160;</span>
<a class="htmltoc" href="#examresults">9.3&#160;&#160;<b>Program Results</b></a>
</div>
</div>
</div>
</div>
</div><h2 class="standard"><a class="sec" name="purpose" id="purpose"/>1&#160;&#160;Purpose</h2>
<div class="paramtext">F08HEF (DSBTRD) reduces a real symmetric band matrix to tridiagonal form.</div><h2 class="standard"><a class="sec" name="specification" id="specification"/>2&#160;&#160;Specification</h2><table class="fspec"><tr><td class="tdfspec1">
<div class="left-tablediv"><table class="fspec1"><tbody>
<tr>
<td class="tdfspec1" valign="top" align="left">SUBROUTINE&#160;F08HEF&#160;(</td>
<td class="tdfspec2" valign="top" align="left"><a class="arg" href="#VECT">VECT</a>, <a class="arg" href="#UPLO">UPLO</a>, <a class="arg" href="#N">N</a>, <a class="arg" href="#KD">KD</a>, <a class="arg" href="#AB">AB</a>, <a class="arg" href="#LDAB">LDAB</a>, <a class="arg" href="#D">D</a>, <a class="arg" href="#E">E</a>, <a class="arg" href="#Q">Q</a>, <a class="arg" href="#LDQ">LDQ</a>, <a class="arg" href="#WORK">WORK</a>, <a class="arg" href="#INFO">INFO</a>)</td>
</tr>
</tbody>
</table></div>
<div class="left-tablediv"><table class="fspec3"><tbody>
<tr>
<td class="tdfspec1" valign="top" align="left">INTEGER&#160;</td>
<td class="tdfspec2" valign="top" align="left">N, KD, LDAB, LDQ, INFO</td>
</tr>
<tr>
<td class="tdfspec1" valign="top" align="left">REAL&#160;(KIND=nag_wp)&#160;</td>
<td class="tdfspec2" valign="top" align="left">AB(LDAB,*), D(N), E(N-1), Q(LDQ,*), WORK(N)</td>
</tr><tr>
<td class="tdfspec1" valign="top" align="left">CHARACTER(1)&#160;</td>
<td class="tdfspec2" valign="top" align="left">VECT, UPLO</td></tr></tbody>
</table></div>
</td></tr></table>
<div class="paramtext">The routine may be called by its 
    LAPACK
    name <span class="bitalic">dsbtrd</span>.</div><h2 class="standard"><a class="sec" name="description" id="description"/>3&#160;&#160;Description</h2>
<div class="paramtext">F08HEF (DSBTRD) reduces a symmetric band matrix <m:math><m:mi>A</m:mi></m:math>&#160;to symmetric tridiagonal form <m:math><m:mi>T</m:mi></m:math>&#160;by an orthogonal similarity transformation: 

<div class="formula"><table class="formula"><tr><td class="formula"><m:math display="block">
 <m:mi>T</m:mi>
 <m:mo>=</m:mo>
 <m:msup><m:mi>Q</m:mi><m:mi mathvariant="normal">T</m:mi></m:msup>
 <m:mi>A</m:mi>
 <m:mi>Q</m:mi>
 <m:mtext>.</m:mtext>
</m:math></td><td class="formula2"/></tr></table></div></div><div class="paramtext">The orthogonal matrix <m:math><m:mi>Q</m:mi></m:math>&#160;is determined as a product of Givens rotation matrices, and may be formed explicitly by the routine if required.</div><div class="paramtext">The routine uses a vectorizable form of the reduction, due to <a class="ref" href="#ref449">Kaufman (1984)</a>.</div><h2 class="standard"><a class="sec" name="references" id="references"/>4&#160;&#160;References</h2><div class="paramtext"><a name="ref449" id="ref449"/>Kaufman L (1984)  Banded eigenvalue solvers on vector machines <i>ACM Trans. Math. Software</i> <b>10</b> 73&#8211;86 </div>
<div class="paramtext"><a name="ref111" id="ref111"/>Parlett B N (1998)  <i>The Symmetric Eigenvalue Problem</i> SIAM, Philadelphia </div><h2 class="standard"><a class="sec" name="parameters" id="parameters"/>5&#160;&#160;Parameters</h2>
<dl><dt class="paramhead"><a name="VECT" id="VECT"/>1: &#160;&#160;&#8194; VECT &#8211; CHARACTER(1)<span class="pclass">Input</span></dt><dd>
<div class="paramtext"><i>On entry</i>: indicates whether <m:math><m:mi>Q</m:mi></m:math>&#160;is to be returned.

<dl>
<dt class="paramval"><m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#VECT"><m:mi mathcolor="#EE0000" mathvariant="bold">VECT</m:mi></m:maction><m:mo>=</m:mo><m:mtext>'V'</m:mtext></m:math></dt>
<dd><m:math><m:mi>Q</m:mi></m:math>&#160;is returned.</dd>
<dt class="paramval"><m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#VECT"><m:mi mathcolor="#EE0000" mathvariant="bold">VECT</m:mi></m:maction><m:mo>=</m:mo><m:mtext>'U'</m:mtext></m:math></dt>
<dd><m:math><m:mi>Q</m:mi></m:math>&#160;is updated (and the array <a class="arg" href="#Q">Q</a> must contain a matrix on entry).</dd>
<dt class="paramval"><m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#VECT"><m:mi mathcolor="#EE0000" mathvariant="bold">VECT</m:mi></m:maction><m:mo>=</m:mo><m:mtext>'N'</m:mtext></m:math></dt>
<dd><m:math><m:mi>Q</m:mi></m:math>&#160;is not required.</dd></dl>
</div><div class="paramtext"><i>Constraint</i>:
  <m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#VECT"><m:mi mathcolor="#EE0000" mathvariant="bold">VECT</m:mi></m:maction><m:mo>=</m:mo><m:mtext>'V'</m:mtext></m:math>, <m:math><m:mtext>'U'</m:mtext></m:math>&#160;or <m:math><m:mtext>'N'</m:mtext></m:math>.
</div>
</dd><dt class="paramhead"><a name="UPLO" id="UPLO"/>2: &#160;&#160;&#8194; UPLO &#8211; CHARACTER(1)<span class="pclass">Input</span></dt><dd>
<div class="paramtext"><i>On entry</i>: indicates whether the upper or lower triangular part of <m:math><m:mi>A</m:mi></m:math>&#160;is stored.

<dl>
<dt class="paramval"><m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#UPLO"><m:mi mathcolor="#EE0000" mathvariant="bold">UPLO</m:mi></m:maction><m:mo>=</m:mo><m:mtext>'U'</m:mtext></m:math></dt>
<dd>The upper triangular part of <m:math><m:mi>A</m:mi></m:math>&#160;is stored.</dd>
<dt class="paramval"><m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#UPLO"><m:mi mathcolor="#EE0000" mathvariant="bold">UPLO</m:mi></m:maction><m:mo>=</m:mo><m:mtext>'L'</m:mtext></m:math></dt>
<dd>The lower triangular part of <m:math><m:mi>A</m:mi></m:math>&#160;is stored.</dd></dl>
</div><div class="paramtext"><i>Constraint</i>:
  <m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#UPLO"><m:mi mathcolor="#EE0000" mathvariant="bold">UPLO</m:mi></m:maction><m:mo>=</m:mo><m:mtext>'U'</m:mtext></m:math>&#160;or <m:math><m:mtext>'L'</m:mtext></m:math>.
</div>
</dd><dt class="paramhead"><a name="N" id="N"/>3: &#160;&#160;&#8194; N &#8211; INTEGER<span class="pclass">Input</span></dt><dd>
<div class="paramtext"><i>On entry</i>: 
<m:math><m:mi>n</m:mi></m:math>, the order of the matrix <m:math><m:mi>A</m:mi></m:math>.</div><div class="paramtext"><i>Constraint</i>:
  <m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#N"><m:mi mathcolor="#EE0000" mathvariant="bold">N</m:mi></m:maction><m:mo>&#8805;</m:mo><m:mn>0</m:mn></m:math>.
</div>
</dd><dt class="paramhead"><a name="KD" id="KD"/>4: &#160;&#160;&#8194; KD &#8211; INTEGER<span class="pclass">Input</span></dt><dd><div class="paramtext"><i>On entry</i>: if <m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#UPLO"><m:mi mathcolor="#EE0000" mathvariant="bold">UPLO</m:mi></m:maction><m:mo>=</m:mo><m:mtext>'U'</m:mtext></m:math>, the number of superdiagonals, <m:math><m:msub><m:mi>k</m:mi><m:mi>d</m:mi></m:msub></m:math>, of the matrix <m:math><m:mi>A</m:mi></m:math>.
<div class="paramtext">If <m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#UPLO"><m:mi mathcolor="#EE0000" mathvariant="bold">UPLO</m:mi></m:maction><m:mo>=</m:mo><m:mtext>'L'</m:mtext></m:math>, the number of subdiagonals, <m:math><m:msub><m:mi>k</m:mi><m:mi>d</m:mi></m:msub></m:math>, of the matrix <m:math><m:mi>A</m:mi></m:math>.</div></div><div class="paramtext"><i>Constraint</i>:
  <m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#KD"><m:mi mathcolor="#EE0000" mathvariant="bold">KD</m:mi></m:maction><m:mo>&#8805;</m:mo><m:mn>0</m:mn></m:math>.
</div>
</dd><dt class="paramhead"><a name="AB" id="AB"/>5: &#160;&#160;&#8194; AB(<a class="arg" href="#LDAB">LDAB</a>,<m:math><m:mo>*</m:mo></m:math>) &#8211; REAL&#160;(KIND=nag_wp)&#160;array<span class="pclass">Input/Output</span></dt><dd>
<div class="paramtext"><b>Note:</b> the second dimension of the array <a class="arg" href="#AB">AB</a>
must be at least
<m:math><m:mrow><m:mi>max</m:mi><m:mspace width="0.125em"/><m:mfenced separators=""><m:mn>1</m:mn><m:mo>,</m:mo><m:maction actiontype="link" dsi:type="simple" dsi:href="#N"><m:mi mathcolor="#EE0000" mathvariant="bold">N</m:mi></m:maction></m:mfenced></m:mrow></m:math>.</div>
<div class="paramtext"><i>On entry</i>: the upper or lower triangle of the <m:math><m:mi>n</m:mi></m:math>&#160;by <m:math><m:mi>n</m:mi></m:math>&#160;symmetric
band matrix <m:math><m:mi>A</m:mi></m:math>.

<div class="paramtext">The matrix is stored in rows <m:math><m:mn>1</m:mn></m:math>&#160;to <m:math><m:msub><m:mi>k</m:mi><m:mi>d</m:mi></m:msub><m:mo>+</m:mo><m:mn>1</m:mn></m:math>, more precisely,<ul class="listind"><li class="listind">if <m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#UPLO"><m:mi mathcolor="#EE0000" mathvariant="bold">UPLO</m:mi></m:maction><m:mo>=</m:mo><m:mtext>'U'</m:mtext></m:math>, the elements of the upper triangle of <m:math><m:mi>A</m:mi></m:math>&#160;within the band must be stored with element <m:math><m:msub><m:mi>A</m:mi><m:mrow><m:mi>i</m:mi><m:mi>j</m:mi></m:mrow></m:msub></m:math>&#160;in <m:math><m:mrow><m:maction actiontype="link" dsi:type="simple" dsi:href="#AB"><m:mi mathcolor="#EE0000" mathvariant="bold">AB</m:mi></m:maction><m:mfenced separators="," open="(" close=")"><m:mrow><m:msub><m:mi>k</m:mi><m:mi>d</m:mi></m:msub><m:mo>+</m:mo><m:mn>1</m:mn><m:mo>+</m:mo><m:mi>i</m:mi><m:mo>-</m:mo><m:mi>j</m:mi></m:mrow><m:mi>j</m:mi></m:mfenced></m:mrow><m:mtext>&#8203; for &#8203;</m:mtext><m:mrow><m:mi>max</m:mi><m:mspace width="0.125em"/><m:mfenced separators=""><m:mn>1</m:mn><m:mo>,</m:mo><m:mrow><m:mi>j</m:mi><m:mo>-</m:mo><m:msub><m:mi>k</m:mi><m:mi>d</m:mi></m:msub></m:mrow></m:mfenced></m:mrow><m:mo>&#8804;</m:mo><m:mi>i</m:mi><m:mo>&#8804;</m:mo><m:mi>j</m:mi></m:math>;</li><li class="listind">if <m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#UPLO"><m:mi mathcolor="#EE0000" mathvariant="bold">UPLO</m:mi></m:maction><m:mo>=</m:mo><m:mtext>'L'</m:mtext></m:math>, the elements of the lower triangle of <m:math><m:mi>A</m:mi></m:math>&#160;within the band must be stored with element <m:math><m:msub><m:mi>A</m:mi><m:mrow><m:mi>i</m:mi><m:mi>j</m:mi></m:mrow></m:msub></m:math>&#160;in <m:math><m:mrow><m:maction actiontype="link" dsi:type="simple" dsi:href="#AB"><m:mi mathcolor="#EE0000" mathvariant="bold">AB</m:mi></m:maction><m:mfenced separators="," open="(" close=")"><m:mrow><m:mn>1</m:mn><m:mo>+</m:mo><m:mi>i</m:mi><m:mo>-</m:mo><m:mi>j</m:mi></m:mrow><m:mi>j</m:mi></m:mfenced></m:mrow><m:mtext>&#8203; for &#8203;</m:mtext><m:mi>j</m:mi><m:mo>&#8804;</m:mo><m:mi>i</m:mi><m:mo>&#8804;</m:mo><m:mrow><m:mi>min</m:mi><m:mspace width="0.125em"/><m:mfenced separators=""><m:mi>n</m:mi><m:mo>,</m:mo><m:mrow><m:mi>j</m:mi><m:mo>+</m:mo><m:msub><m:mi>k</m:mi><m:mi>d</m:mi></m:msub></m:mrow></m:mfenced></m:mrow><m:mtext>.</m:mtext></m:math></li></ul></div>
</div>
<div class="paramtext"><i>On exit</i>: 
<a class="arg" href="#AB">AB</a> is overwritten by values generated during the reduction to tridiagonal form.
<div class="paramtext">The first superdiagonal and the diagonal of the tridiagonal matrix <m:math><m:mi>T</m:mi></m:math>&#160;are returned in <a class="arg" href="#AB">AB</a> using the same storage format as described above.</div>
</div></dd><dt class="paramhead"><a name="LDAB" id="LDAB"/>6: &#160;&#160;&#8194; LDAB &#8211; INTEGER<span class="pclass">Input</span></dt><dd>
<div class="paramtext"><i>On entry</i>: the first dimension of the array <a class="arg" href="#AB">AB</a> as declared in the (sub)program from which F08HEF (DSBTRD) is called.</div><div class="paramtext"><i>Constraint</i>:
  <m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#LDAB"><m:mi mathcolor="#EE0000" mathvariant="bold">LDAB</m:mi></m:maction><m:mo>&#8805;</m:mo><m:mrow><m:mi>max</m:mi><m:mspace width="0.125em"/><m:mfenced separators=""><m:mn>1</m:mn><m:mo>,</m:mo><m:mrow><m:maction actiontype="link" dsi:type="simple" dsi:href="#KD"><m:mi mathcolor="#EE0000" mathvariant="bold">KD</m:mi></m:maction><m:mo>+</m:mo><m:mn>1</m:mn></m:mrow></m:mfenced></m:mrow></m:math>.
</div>
</dd><dt class="paramhead"><a name="D" id="D"/>7: &#160;&#160;&#8194; D(<a class="arg" href="#N">N</a>) &#8211; REAL&#160;(KIND=nag_wp)&#160;array<span class="pclass">Output</span></dt><dd><div class="paramtext"><i>On exit</i>: the diagonal elements of the tridiagonal matrix <m:math><m:mi>T</m:mi></m:math>.</div>
</dd><dt class="paramhead"><a name="E" id="E"/>8: &#160;&#160;&#8194; E(<m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#N"><m:mi mathcolor="#EE0000" mathvariant="bold">N</m:mi></m:maction><m:mo>-</m:mo><m:mn>1</m:mn></m:math>) &#8211; REAL&#160;(KIND=nag_wp)&#160;array<span class="pclass">Output</span></dt><dd><div class="paramtext"><i>On exit</i>: the off-diagonal elements of the tridiagonal matrix <m:math><m:mi>T</m:mi></m:math>.</div>
</dd><dt class="paramhead"><a name="Q" id="Q"/>9: &#160;&#160;&#8194; Q(<a class="arg" href="#LDQ">LDQ</a>,<m:math><m:mo>*</m:mo></m:math>) &#8211; REAL&#160;(KIND=nag_wp)&#160;array<span class="pclass">Input/Output</span></dt><dd>
<div class="paramtext"><b>Note:</b> the second dimension of the array <a class="arg" href="#Q">Q</a>
must be at least
<m:math><m:mrow><m:mi>max</m:mi><m:mspace width="0.125em"/><m:mfenced separators=""><m:mn>1</m:mn><m:mo>,</m:mo><m:maction actiontype="link" dsi:type="simple" dsi:href="#N"><m:mi mathcolor="#EE0000" mathvariant="bold">N</m:mi></m:maction></m:mfenced></m:mrow></m:math>&#160;if <m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#VECT"><m:mi mathcolor="#EE0000" mathvariant="bold">VECT</m:mi></m:maction><m:mo>=</m:mo><m:mtext>'V'</m:mtext></m:math>&#160;or <m:math><m:mtext>'U'</m:mtext></m:math>&#160;and at least <m:math><m:mn>1</m:mn></m:math>&#160;if <m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#VECT"><m:mi mathcolor="#EE0000" mathvariant="bold">VECT</m:mi></m:maction><m:mo>=</m:mo><m:mtext>'N'</m:mtext></m:math>.</div>
<div class="paramtext"><i>On entry</i>: if <m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#VECT"><m:mi mathcolor="#EE0000" mathvariant="bold">VECT</m:mi></m:maction><m:mo>=</m:mo><m:mtext>'U'</m:mtext></m:math>, <a class="arg" href="#Q">Q</a> must contain the matrix formed in a previous stage of the reduction (for example, the reduction of a banded symmetric-definite generalized eigenproblem); otherwise <a class="arg" href="#Q">Q</a> need not be set.</div>
<div class="paramtext"><i>On exit</i>: if <m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#VECT"><m:mi mathcolor="#EE0000" mathvariant="bold">VECT</m:mi></m:maction><m:mo>=</m:mo><m:mtext>'V'</m:mtext></m:math>&#160;or <m:math><m:mtext>'U'</m:mtext></m:math>, the <m:math><m:mi>n</m:mi></m:math>&#160;by <m:math><m:mi>n</m:mi></m:math>&#160;matrix <m:math><m:mi>Q</m:mi></m:math>.
<div class="paramtext">If <m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#VECT"><m:mi mathcolor="#EE0000" mathvariant="bold">VECT</m:mi></m:maction><m:mo>=</m:mo><m:mtext>'N'</m:mtext></m:math>, <a class="arg" href="#Q">Q</a> is not referenced.</div>
</div>
</dd><dt class="paramhead"><a name="LDQ" id="LDQ"/>10: &#8194; LDQ &#8211; INTEGER<span class="pclass">Input</span></dt><dd>
<div class="paramtext"><i>On entry</i>: the first dimension of the array <a class="arg" href="#Q">Q</a> as declared in the (sub)program from which F08HEF (DSBTRD) is called.</div><div class="paramtext"><i>Constraints</i>:
   <div class="paramtext"/><ul class="listcons">
<li class="listcons">if <m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#VECT"><m:mi mathcolor="#EE0000" mathvariant="bold">VECT</m:mi></m:maction><m:mo>=</m:mo><m:mtext>'V'</m:mtext></m:math>&#160;or <m:math><m:mtext>'U'</m:mtext></m:math>, <m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#LDQ"><m:mi mathcolor="#EE0000" mathvariant="bold">LDQ</m:mi></m:maction><m:mo>&#8805;</m:mo><m:mrow><m:mi>max</m:mi><m:mspace width="0.125em"/><m:mfenced separators=""><m:mn>1</m:mn><m:mo>,</m:mo><m:maction actiontype="link" dsi:type="simple" dsi:href="#N"><m:mi mathcolor="#EE0000" mathvariant="bold">N</m:mi></m:maction></m:mfenced></m:mrow></m:math>;</li>
<li class="listcons">if <m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#VECT"><m:mi mathcolor="#EE0000" mathvariant="bold">VECT</m:mi></m:maction><m:mo>=</m:mo><m:mtext>'N'</m:mtext></m:math>, <m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#LDQ"><m:mi mathcolor="#EE0000" mathvariant="bold">LDQ</m:mi></m:maction><m:mo>&#8805;</m:mo><m:mn>1</m:mn></m:math>.</li>
</ul></div>
</dd><dt class="paramhead"><a name="WORK" id="WORK"/>11: &#8194; WORK(<a class="arg" href="#N">N</a>) &#8211; REAL&#160;(KIND=nag_wp)&#160;array<span class="pclass">Workspace</span></dt><dt class="paramhead"><a name="INFO" id="INFO"/>12: &#8194; INFO &#8211; INTEGER<span class="pclass">Output</span></dt><dd><div class="paramtext"><i>On exit</i>: <m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#INFO"><m:mi mathcolor="#EE0000" mathvariant="bold">INFO</m:mi></m:maction><m:mo>=</m:mo><m:mn>0</m:mn></m:math>&#160;unless the routine detects an error (see <a class="sec" href="#errors">Section 6</a>).</div></dd></dl><h2 class="standard"><a class="sec" name="errors" id="errors"/>6&#160;&#160;Error Indicators and Warnings</h2>
<div class="paramtext">Errors or warnings detected by the routine:</div>
<dl class="ifail">
<dt class="errorhead"><a name="INlt0" id="INlt0"/><m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#INFO"><m:mi mathcolor="#EE0000" mathvariant="bold">INFO</m:mi></m:maction><m:mo>&lt;</m:mo><m:mn>0</m:mn></m:math></dt>
<dd><div class="paramtext">If <m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#INFO"><m:mi mathcolor="#EE0000" mathvariant="bold">INFO</m:mi></m:maction><m:mo>=</m:mo><m:mo>-</m:mo><m:mi>i</m:mi></m:math>, argument <m:math><m:mi>i</m:mi></m:math>&#160;had an illegal value. An explanatory message is output, and execution of the program is terminated.</div></dd>
</dl><h2 class="standard"><a class="sec" name="accuracy" id="accuracy"/>7&#160;&#160;Accuracy</h2>
<div class="paramtext">The computed tridiagonal matrix <m:math><m:mi>T</m:mi></m:math>&#160;is exactly similar to a nearby matrix <m:math><m:mfenced separators=""><m:mi>A</m:mi><m:mo>+</m:mo><m:mi>E</m:mi></m:mfenced></m:math>, where

<div class="formula"><table class="formula"><tr><td class="formula"><m:math display="block">
 <m:msub><m:mfenced open="&#8214;" close="&#8214;" separators=""><m:mi>E</m:mi></m:mfenced><m:mn>2</m:mn></m:msub><m:mo>&#8804;</m:mo>
 <m:mi>c</m:mi>
 <m:mfenced separators=""><m:mi>n</m:mi></m:mfenced>
 <m:mi>&#949;</m:mi>
 <m:msub><m:mfenced open="&#8214;" close="&#8214;" separators=""><m:mi>A</m:mi></m:mfenced><m:mn>2</m:mn></m:msub>
 <m:mtext>,</m:mtext>
</m:math></td><td class="formula2"/></tr></table></div><m:math><m:mi>c</m:mi><m:mfenced separators=""><m:mi>n</m:mi></m:mfenced></m:math>&#160;is a modestly increasing function of <m:math><m:mi>n</m:mi></m:math>, and <m:math><m:mi>&#949;</m:mi></m:math>&#160;is the <span class="bitalic">machine precision</span>.</div><div class="paramtext">The elements of <m:math><m:mi>T</m:mi></m:math>&#160;themselves may be sensitive to small perturbations in <m:math><m:mi>A</m:mi></m:math>&#160;or to rounding errors in the computation, but this does not affect the stability of the eigenvalues and eigenvectors.</div><div class="paramtext">The computed matrix <m:math><m:mi>Q</m:mi></m:math>&#160;differs from an exactly orthogonal matrix by a matrix <m:math><m:mi>E</m:mi></m:math>&#160;such that

<div class="formula"><table class="formula"><tr><td class="formula"><m:math display="block">
 <m:msub><m:mfenced open="&#8214;" close="&#8214;" separators=""><m:mi>E</m:mi></m:mfenced><m:mn>2</m:mn></m:msub>
 <m:mo>=</m:mo>
 <m:mrow><m:mi mathvariant="italic">O</m:mi><m:mfenced separators=""><m:mi>&#949;</m:mi></m:mfenced></m:mrow>
 <m:mtext>,</m:mtext>
</m:math></td><td class="formula2"/></tr></table></div>

where <m:math><m:mi>&#949;</m:mi></m:math>&#160;is the <span class="bitalic">machine precision</span>.</div><h2 class="standard"><a class="sec" name="fcomments" id="fcomments"/>8&#160;&#160;Further Comments</h2>
<div class="paramtext">The total number of floating point operations is approximately <m:math><m:mn>6</m:mn><m:msup><m:mi>n</m:mi><m:mn>2</m:mn></m:msup><m:mi>k</m:mi></m:math>&#160;if <m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#VECT"><m:mi mathcolor="#EE0000" mathvariant="bold">VECT</m:mi></m:maction><m:mo>=</m:mo><m:mtext>'N'</m:mtext></m:math>&#160;with <m:math><m:mn>3</m:mn><m:msup><m:mi>n</m:mi><m:mn>3</m:mn></m:msup><m:mfenced separators=""><m:mi>k</m:mi><m:mo>-</m:mo><m:mn>1</m:mn></m:mfenced><m:mo>/</m:mo><m:mi>k</m:mi></m:math>&#160;additional operations if <m:math><m:maction actiontype="link" dsi:type="simple" dsi:href="#VECT"><m:mi mathcolor="#EE0000" mathvariant="bold">VECT</m:mi></m:maction><m:mo>=</m:mo><m:mtext>'V'</m:mtext></m:math>.</div><div class="paramtext">The complex analogue of this routine is <a class="rout" href="../F08/f08hsf.xml">F08HSF (ZHBTRD)</a>.</div><h2 class="standard"><a class="sec" name="example" id="example"/>9&#160;&#160;Example</h2>
<div class="paramtext">This example computes all the eigenvalues and eigenvectors of the matrix <m:math><m:mi>A</m:mi></m:math>, where

<div class="formula"><table class="formula"><tr><td class="formula"><m:math display="block">
 <m:mi>A</m:mi>
 <m:mo>=</m:mo>
 <m:mfenced><m:mtable columnalign="right">
<m:mtr>
   <m:mtd><m:mn>4.99</m:mn></m:mtd>
   <m:mtd><m:mn>0.04</m:mn></m:mtd>
   <m:mtd><m:mn>0.22</m:mn></m:mtd>
   <m:mtd><m:mn>0.00</m:mn></m:mtd>
</m:mtr><m:mtr>
   <m:mtd><m:mn>0.04</m:mn></m:mtd>
   <m:mtd><m:mn>1.05</m:mn></m:mtd>
   <m:mtd><m:mrow><m:mo>-</m:mo><m:mn>0.79</m:mn></m:mrow></m:mtd>
   <m:mtd><m:mn>1.04</m:mn></m:mtd>
</m:mtr><m:mtr>
   <m:mtd><m:mn>0.22</m:mn></m:mtd>
   <m:mtd><m:mrow><m:mo>-</m:mo><m:mn>0.79</m:mn></m:mrow></m:mtd>
   <m:mtd><m:mrow><m:mo>-</m:mo><m:mn>2.31</m:mn></m:mrow></m:mtd>
   <m:mtd><m:mrow><m:mo>-</m:mo><m:mn>1.30</m:mn></m:mrow></m:mtd>
</m:mtr><m:mtr>
   <m:mtd><m:mn>0.00</m:mn></m:mtd>
   <m:mtd><m:mn>1.04</m:mn></m:mtd>
   <m:mtd><m:mrow><m:mo>-</m:mo><m:mn>1.30</m:mn></m:mrow></m:mtd>
   <m:mtd><m:mrow><m:mo>-</m:mo><m:mn>0.43</m:mn></m:mrow></m:mtd>
</m:mtr>
</m:mtable></m:mfenced>
<m:mtext>.</m:mtext>
</m:math></td><td class="formula2"/></tr></table></div>

Here <m:math><m:mi>A</m:mi></m:math>&#160;is symmetric and is treated as a band matrix.  The program first calls F08HEF (DSBTRD) to reduce <m:math><m:mi>A</m:mi></m:math>&#160;to tridiagonal form <m:math><m:mi>T</m:mi></m:math>, and to form the orthogonal matrix <m:math><m:mi>Q</m:mi></m:math>; the results are then passed to <a class="rout" href="../F08/f08jef.xml">F08JEF (DSTEQR)</a> which computes the eigenvalues and eigenvectors of <m:math><m:mi>A</m:mi></m:math>.</div><h3 class="standard"><a class="sec" name="examtext" id="examtext"/>9.1&#160;&#160;Program Text</h3>
<p><a class="verbatimref" href="../../examples/source/f08hefe.f90">Program Text (f08hefe.f90)</a></p><h3 class="standard"><a class="sec" name="examdata" id="examdata"/>9.2&#160;&#160;Program Data</h3>
<p><a class="verbatimref" href="../../examples/data/f08hefe.d">Program&#160;Data (f08hefe.d)</a></p><h3 class="standard"><a class="sec" name="examresults" id="examresults"/>9.3&#160;&#160;Program Results</h3>
<p><a class="verbatimref" href="../../examples/baseresults/f08hefe.r">Program Results (f08hefe.r)</a></p>
<hr/><div><a class="rout" href="../../pdf/F08/f08hef.pdf">F08HEF (DSBTRD) (PDF version)</a></div><div><a class="chap" href="f08conts.xml">F08 Chapter Contents</a></div><div><a class="chapint" href="f08intro.xml">F08 Chapter Introduction</a></div>
<div><a class="htmltoc" href="../FRONTMATTER/manconts.xml">NAG Library Manual</a></div>
<div><hr/><a class="genint" href="../FRONTMATTER/copyright.xml">&#169; The Numerical Algorithms Group Ltd, Oxford, UK. 2011</a></div></body></html>