/usr/share/doc/lorene-doc/refguide/diff__xdsdx2_8C_source.html is in lorene-doc 0.0.0~cvs20161116+dfsg-1.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
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 | <!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "http://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd">
<html xmlns="http://www.w3.org/1999/xhtml">
<head>
<meta http-equiv="Content-Type" content="text/xhtml;charset=UTF-8"/>
<meta http-equiv="X-UA-Compatible" content="IE=9"/>
<meta name="generator" content="Doxygen 1.8.12"/>
<meta name="viewport" content="width=device-width, initial-scale=1"/>
<title>LORENE: diff_xdsdx2.C Source File</title>
<link href="tabs.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="jquery.js"></script>
<script type="text/javascript" src="dynsections.js"></script>
<link href="navtree.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="resize.js"></script>
<script type="text/javascript" src="navtreedata.js"></script>
<script type="text/javascript" src="navtree.js"></script>
<script type="text/javascript">
$(document).ready(initResizable);
</script>
<link href="doxygen.css" rel="stylesheet" type="text/css" />
</head>
<body>
<div id="top"><!-- do not remove this div, it is closed by doxygen! -->
<div id="titlearea">
<table cellspacing="0" cellpadding="0">
<tbody>
<tr style="height: 56px;">
<td id="projectalign" style="padding-left: 0.5em;">
<div id="projectname">LORENE
</div>
</td>
</tr>
</tbody>
</table>
</div>
<!-- end header part -->
<!-- Generated by Doxygen 1.8.12 -->
<script type="text/javascript" src="menudata.js"></script>
<script type="text/javascript" src="menu.js"></script>
<script type="text/javascript">
$(function() {
initMenu('',false,false,'search.php','Search');
});
</script>
<div id="main-nav"></div>
</div><!-- top -->
<div id="side-nav" class="ui-resizable side-nav-resizable">
<div id="nav-tree">
<div id="nav-tree-contents">
<div id="nav-sync" class="sync"></div>
</div>
</div>
<div id="splitbar" style="-moz-user-select:none;"
class="ui-resizable-handle">
</div>
</div>
<script type="text/javascript">
$(document).ready(function(){initNavTree('diff__xdsdx2_8C_source.html','');});
</script>
<div id="doc-content">
<div class="header">
<div class="headertitle">
<div class="title">diff_xdsdx2.C</div> </div>
</div><!--header-->
<div class="contents">
<div class="fragment"><div class="line"><a name="l00001"></a><span class="lineno"> 1</span> <span class="comment">/*</span></div><div class="line"><a name="l00002"></a><span class="lineno"> 2</span> <span class="comment"> * Methods for the class Diff_xdsdx2</span></div><div class="line"><a name="l00003"></a><span class="lineno"> 3</span> <span class="comment"> *</span></div><div class="line"><a name="l00004"></a><span class="lineno"> 4</span> <span class="comment"> * (see file diff.h for documentation).</span></div><div class="line"><a name="l00005"></a><span class="lineno"> 5</span> <span class="comment"> *</span></div><div class="line"><a name="l00006"></a><span class="lineno"> 6</span> <span class="comment"> */</span></div><div class="line"><a name="l00007"></a><span class="lineno"> 7</span> </div><div class="line"><a name="l00008"></a><span class="lineno"> 8</span> <span class="comment">/*</span></div><div class="line"><a name="l00009"></a><span class="lineno"> 9</span> <span class="comment"> * Copyright (c) 2005 Jerome Novak</span></div><div class="line"><a name="l00010"></a><span class="lineno"> 10</span> <span class="comment"> *</span></div><div class="line"><a name="l00011"></a><span class="lineno"> 11</span> <span class="comment"> * This file is part of LORENE.</span></div><div class="line"><a name="l00012"></a><span class="lineno"> 12</span> <span class="comment"> *</span></div><div class="line"><a name="l00013"></a><span class="lineno"> 13</span> <span class="comment"> * LORENE is free software; you can redistribute it and/or modify</span></div><div class="line"><a name="l00014"></a><span class="lineno"> 14</span> <span class="comment"> * it under the terms of the GNU General Public License version 2</span></div><div class="line"><a name="l00015"></a><span class="lineno"> 15</span> <span class="comment"> * as published by the Free Software Foundation.</span></div><div class="line"><a name="l00016"></a><span class="lineno"> 16</span> <span class="comment"> *</span></div><div class="line"><a name="l00017"></a><span class="lineno"> 17</span> <span class="comment"> * LORENE is distributed in the hope that it will be useful,</span></div><div class="line"><a name="l00018"></a><span class="lineno"> 18</span> <span class="comment"> * but WITHOUT ANY WARRANTY; without even the implied warranty of</span></div><div class="line"><a name="l00019"></a><span class="lineno"> 19</span> <span class="comment"> * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the</span></div><div class="line"><a name="l00020"></a><span class="lineno"> 20</span> <span class="comment"> * GNU General Public License for more details.</span></div><div class="line"><a name="l00021"></a><span class="lineno"> 21</span> <span class="comment"> *</span></div><div class="line"><a name="l00022"></a><span class="lineno"> 22</span> <span class="comment"> * You should have received a copy of the GNU General Public License</span></div><div class="line"><a name="l00023"></a><span class="lineno"> 23</span> <span class="comment"> * along with LORENE; if not, write to the Free Software</span></div><div class="line"><a name="l00024"></a><span class="lineno"> 24</span> <span class="comment"> * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA</span></div><div class="line"><a name="l00025"></a><span class="lineno"> 25</span> <span class="comment"> *</span></div><div class="line"><a name="l00026"></a><span class="lineno"> 26</span> <span class="comment"> */</span></div><div class="line"><a name="l00027"></a><span class="lineno"> 27</span> </div><div class="line"><a name="l00028"></a><span class="lineno"> 28</span> <span class="keywordtype">char</span> diff_xdsdx2_C[] = <span class="stringliteral">"$Header: /cvsroot/Lorene/C++/Source/Diff/diff_xdsdx2.C,v 1.4 2014/10/13 08:52:51 j_novak Exp $"</span> ;</div><div class="line"><a name="l00029"></a><span class="lineno"> 29</span> </div><div class="line"><a name="l00030"></a><span class="lineno"> 30</span> <span class="comment">/*</span></div><div class="line"><a name="l00031"></a><span class="lineno"> 31</span> <span class="comment"> * $Id: diff_xdsdx2.C,v 1.4 2014/10/13 08:52:51 j_novak Exp $</span></div><div class="line"><a name="l00032"></a><span class="lineno"> 32</span> <span class="comment"> * $Log: diff_xdsdx2.C,v $</span></div><div class="line"><a name="l00033"></a><span class="lineno"> 33</span> <span class="comment"> * Revision 1.4 2014/10/13 08:52:51 j_novak</span></div><div class="line"><a name="l00034"></a><span class="lineno"> 34</span> <span class="comment"> * Lorene classes and functions now belong to the namespace Lorene.</span></div><div class="line"><a name="l00035"></a><span class="lineno"> 35</span> <span class="comment"> *</span></div><div class="line"><a name="l00036"></a><span class="lineno"> 36</span> <span class="comment"> * Revision 1.3 2014/10/06 15:13:05 j_novak</span></div><div class="line"><a name="l00037"></a><span class="lineno"> 37</span> <span class="comment"> * Modified #include directives to use c++ syntax.</span></div><div class="line"><a name="l00038"></a><span class="lineno"> 38</span> <span class="comment"> *</span></div><div class="line"><a name="l00039"></a><span class="lineno"> 39</span> <span class="comment"> * Revision 1.2 2007/12/11 15:28:11 jl_cornou</span></div><div class="line"><a name="l00040"></a><span class="lineno"> 40</span> <span class="comment"> * Jacobi(0,2) polynomials partially implemented</span></div><div class="line"><a name="l00041"></a><span class="lineno"> 41</span> <span class="comment"> *</span></div><div class="line"><a name="l00042"></a><span class="lineno"> 42</span> <span class="comment"> * Revision 1.1 2005/01/10 16:34:52 j_novak</span></div><div class="line"><a name="l00043"></a><span class="lineno"> 43</span> <span class="comment"> * New class for 1D mono-domain differential operators.</span></div><div class="line"><a name="l00044"></a><span class="lineno"> 44</span> <span class="comment"> *</span></div><div class="line"><a name="l00045"></a><span class="lineno"> 45</span> <span class="comment"> *</span></div><div class="line"><a name="l00046"></a><span class="lineno"> 46</span> <span class="comment"> * $Header: /cvsroot/Lorene/C++/Source/Diff/diff_xdsdx2.C,v 1.4 2014/10/13 08:52:51 j_novak Exp $</span></div><div class="line"><a name="l00047"></a><span class="lineno"> 47</span> <span class="comment"> *</span></div><div class="line"><a name="l00048"></a><span class="lineno"> 48</span> <span class="comment"> */</span></div><div class="line"><a name="l00049"></a><span class="lineno"> 49</span> </div><div class="line"><a name="l00050"></a><span class="lineno"> 50</span> <span class="comment">// C headers</span></div><div class="line"><a name="l00051"></a><span class="lineno"> 51</span> <span class="preprocessor">#include <cassert></span></div><div class="line"><a name="l00052"></a><span class="lineno"> 52</span> <span class="preprocessor">#include <cstdlib></span></div><div class="line"><a name="l00053"></a><span class="lineno"> 53</span> </div><div class="line"><a name="l00054"></a><span class="lineno"> 54</span> <span class="comment">// Lorene headers</span></div><div class="line"><a name="l00055"></a><span class="lineno"> 55</span> <span class="preprocessor">#include "diff.h"</span></div><div class="line"><a name="l00056"></a><span class="lineno"> 56</span> <span class="preprocessor">#include "proto.h"</span></div><div class="line"><a name="l00057"></a><span class="lineno"> 57</span> </div><div class="line"><a name="l00058"></a><span class="lineno"> 58</span> <span class="keyword">namespace </span><a class="code" href="namespaceLorene.html">Lorene</a> {</div><div class="line"><a name="l00059"></a><span class="lineno"> 59</span> <span class="keywordtype">void</span> multxpun_1d(<span class="keywordtype">int</span>, <span class="keywordtype">double</span>**, <span class="keywordtype">int</span>) ;</div><div class="line"><a name="l00060"></a><span class="lineno"> 60</span> </div><div class="line"><a name="l00061"></a><span class="lineno"> 61</span> <span class="keyword">namespace </span>{</div><div class="line"><a name="l00062"></a><span class="lineno"> 62</span>  <span class="keywordtype">int</span> nap = 0 ;</div><div class="line"><a name="l00063"></a><span class="lineno"> 63</span>  Matrice* tab[<a class="code" href="group__def__mac.html#ga32503ae89f3881e280be79eefc4da0d7">MAX_BASE</a>*<a class="code" href="classLorene_1_1Diff.html#ad9529ab3d6813ac08e14689e34e2b9f4">Diff::max_points</a>] ;</div><div class="line"><a name="l00064"></a><span class="lineno"> 64</span>  <span class="keywordtype">int</span> nr_done[<a class="code" href="classLorene_1_1Diff.html#ad9529ab3d6813ac08e14689e34e2b9f4">Diff::max_points</a>] ;</div><div class="line"><a name="l00065"></a><span class="lineno"> 65</span> }</div><div class="line"><a name="l00066"></a><span class="lineno"> 66</span> </div><div class="line"><a name="l00067"></a><span class="lineno"><a class="line" href="classLorene_1_1Diff__xdsdx2.html#a790da66fbf295a3b65bcfe69052e256d"> 67</a></span> <a class="code" href="classLorene_1_1Diff__xdsdx2.html#a790da66fbf295a3b65bcfe69052e256d">Diff_xdsdx2::Diff_xdsdx2</a>(<span class="keywordtype">int</span> base_r, <span class="keywordtype">int</span> nr) : <a class="code" href="classLorene_1_1Diff.html">Diff</a>(base_r, nr) {</div><div class="line"><a name="l00068"></a><span class="lineno"> 68</span>  <a class="code" href="classLorene_1_1Diff__xdsdx2.html#a4adfe4c96630d613d3f1e10af64d65d9">initialize</a>() ;</div><div class="line"><a name="l00069"></a><span class="lineno"> 69</span>  assert ((<a class="code" href="classLorene_1_1Diff.html#afbaa55cedbd7594c4fdc9c085d79ee3c">base</a> != <a class="code" href="group__def__mac.html#gae7dd6de71b3796c971acf0b09982951c">R_CHEBP</a>)&&(<a class="code" href="classLorene_1_1Diff.html#afbaa55cedbd7594c4fdc9c085d79ee3c">base</a> != <a class="code" href="group__def__mac.html#ga610567430ec691ccd8f1120da031c3e0">R_CHEBI</a>)) ;</div><div class="line"><a name="l00070"></a><span class="lineno"> 70</span> }</div><div class="line"><a name="l00071"></a><span class="lineno"> 71</span> </div><div class="line"><a name="l00072"></a><span class="lineno"><a class="line" href="classLorene_1_1Diff__xdsdx2.html#a4f5f891f4d2c2eaaa494d8c870f317a2"> 72</a></span> <a class="code" href="classLorene_1_1Diff__xdsdx2.html#a790da66fbf295a3b65bcfe69052e256d">Diff_xdsdx2::Diff_xdsdx2</a>(<span class="keyword">const</span> <a class="code" href="classLorene_1_1Diff__xdsdx2.html">Diff_xdsdx2</a>& diff_in) : <a class="code" href="classLorene_1_1Diff.html">Diff</a>(diff_in) {</div><div class="line"><a name="l00073"></a><span class="lineno"> 73</span>  assert (nap != 0) ;</div><div class="line"><a name="l00074"></a><span class="lineno"> 74</span>  assert ((<a class="code" href="classLorene_1_1Diff.html#afbaa55cedbd7594c4fdc9c085d79ee3c">base</a> != <a class="code" href="group__def__mac.html#gae7dd6de71b3796c971acf0b09982951c">R_CHEBP</a>)&&(<a class="code" href="classLorene_1_1Diff.html#afbaa55cedbd7594c4fdc9c085d79ee3c">base</a> != <a class="code" href="group__def__mac.html#ga610567430ec691ccd8f1120da031c3e0">R_CHEBI</a>)) ;</div><div class="line"><a name="l00075"></a><span class="lineno"> 75</span> } </div><div class="line"><a name="l00076"></a><span class="lineno"> 76</span> </div><div class="line"><a name="l00077"></a><span class="lineno"><a class="line" href="classLorene_1_1Diff__xdsdx2.html#a97d108d41d21c6131d981683def9c77a"> 77</a></span> <a class="code" href="classLorene_1_1Diff__xdsdx2.html#a97d108d41d21c6131d981683def9c77a">Diff_xdsdx2::~Diff_xdsdx2</a>() {}</div><div class="line"><a name="l00078"></a><span class="lineno"> 78</span> </div><div class="line"><a name="l00079"></a><span class="lineno"><a class="line" href="classLorene_1_1Diff__xdsdx2.html#a4adfe4c96630d613d3f1e10af64d65d9"> 79</a></span> <span class="keywordtype">void</span> <a class="code" href="classLorene_1_1Diff__xdsdx2.html#a4adfe4c96630d613d3f1e10af64d65d9">Diff_xdsdx2::initialize</a>() {</div><div class="line"><a name="l00080"></a><span class="lineno"> 80</span>  <span class="keywordflow">if</span> (nap == 0) {</div><div class="line"><a name="l00081"></a><span class="lineno"> 81</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> i=0; i<<a class="code" href="classLorene_1_1Diff.html#ad9529ab3d6813ac08e14689e34e2b9f4">max_points</a>; i++) {</div><div class="line"><a name="l00082"></a><span class="lineno"> 82</span>  nr_done[i] = -1 ;</div><div class="line"><a name="l00083"></a><span class="lineno"> 83</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> j=0; j<<a class="code" href="group__def__mac.html#ga32503ae89f3881e280be79eefc4da0d7">MAX_BASE</a>; j++) </div><div class="line"><a name="l00084"></a><span class="lineno"> 84</span>  tab[j*max_points+i] = 0x0 ;</div><div class="line"><a name="l00085"></a><span class="lineno"> 85</span>  }</div><div class="line"><a name="l00086"></a><span class="lineno"> 86</span>  nap = 1 ;</div><div class="line"><a name="l00087"></a><span class="lineno"> 87</span>  }</div><div class="line"><a name="l00088"></a><span class="lineno"> 88</span>  return ;</div><div class="line"><a name="l00089"></a><span class="lineno"> 89</span> }</div><div class="line"><a name="l00090"></a><span class="lineno"> 90</span> </div><div class="line"><a name="l00091"></a><span class="lineno"><a class="line" href="classLorene_1_1Diff__xdsdx2.html#a6bd20c3714234d43330e45230f5eab9d"> 91</a></span> <span class="keywordtype">void</span> <a class="code" href="classLorene_1_1Diff__xdsdx2.html#a6bd20c3714234d43330e45230f5eab9d">Diff_xdsdx2::operator=</a>(<span class="keyword">const</span> <a class="code" href="classLorene_1_1Diff__xdsdx2.html">Diff_xdsdx2</a>& diff_in) {</div><div class="line"><a name="l00092"></a><span class="lineno"> 92</span>  assert (nap != 0) ;</div><div class="line"><a name="l00093"></a><span class="lineno"> 93</span>  <a class="code" href="classLorene_1_1Diff.html#a1768a197f28eb6c7735d3eb672e26aab">Diff::operator=</a>(diff_in) ;</div><div class="line"><a name="l00094"></a><span class="lineno"> 94</span> </div><div class="line"><a name="l00095"></a><span class="lineno"> 95</span> }</div><div class="line"><a name="l00096"></a><span class="lineno"> 96</span> </div><div class="line"><a name="l00097"></a><span class="lineno"><a class="line" href="classLorene_1_1Diff__xdsdx2.html#a76832bcbcf932e997767f94fdaddb789"> 97</a></span> <span class="keyword">const</span> <a class="code" href="classLorene_1_1Matrice.html">Matrice</a>& <a class="code" href="classLorene_1_1Diff__xdsdx2.html#a76832bcbcf932e997767f94fdaddb789">Diff_xdsdx2::get_matrice</a>()<span class="keyword"> const </span>{</div><div class="line"><a name="l00098"></a><span class="lineno"> 98</span>  </div><div class="line"><a name="l00099"></a><span class="lineno"> 99</span>  <span class="keywordtype">bool</span> done = false ;</div><div class="line"><a name="l00100"></a><span class="lineno"> 100</span>  <span class="keywordtype">int</span> indice ;</div><div class="line"><a name="l00101"></a><span class="lineno"> 101</span>  <span class="keywordflow">for</span> (indice =0; indice<<a class="code" href="classLorene_1_1Diff.html#ad9529ab3d6813ac08e14689e34e2b9f4">max_points</a>; indice++) {</div><div class="line"><a name="l00102"></a><span class="lineno"> 102</span>  <span class="keywordflow">if</span> (nr_done[indice] == <a class="code" href="classLorene_1_1Diff.html#a14adb98c5740935b2a2771ba605fe4f4">npoints</a>) {</div><div class="line"><a name="l00103"></a><span class="lineno"> 103</span>  <span class="keywordflow">if</span> (tab[<a class="code" href="classLorene_1_1Diff.html#afbaa55cedbd7594c4fdc9c085d79ee3c">base</a>*max_points + indice] != 0x0) done = true ;</div><div class="line"><a name="l00104"></a><span class="lineno"> 104</span>  break ;</div><div class="line"><a name="l00105"></a><span class="lineno"> 105</span>  }</div><div class="line"><a name="l00106"></a><span class="lineno"> 106</span>  <span class="keywordflow">if</span> (nr_done[indice] == -1)</div><div class="line"><a name="l00107"></a><span class="lineno"> 107</span>  break ;</div><div class="line"><a name="l00108"></a><span class="lineno"> 108</span>  }</div><div class="line"><a name="l00109"></a><span class="lineno"> 109</span>  <span class="keywordflow">if</span> (!done) { <span class="comment">//The computation must be done ...</span></div><div class="line"><a name="l00110"></a><span class="lineno"> 110</span>  <span class="keywordflow">if</span> (indice == max_points) {</div><div class="line"><a name="l00111"></a><span class="lineno"> 111</span>  cerr << <span class="stringliteral">"Diff_xdsdx2::get_matrice() : no space left!!"</span> << <span class="charliteral">'\n'</span></div><div class="line"><a name="l00112"></a><span class="lineno"> 112</span>  << <span class="stringliteral">"The value of Diff.max_points must be increased..."</span> << endl ;</div><div class="line"><a name="l00113"></a><span class="lineno"> 113</span>  abort() ;</div><div class="line"><a name="l00114"></a><span class="lineno"> 114</span>  }</div><div class="line"><a name="l00115"></a><span class="lineno"> 115</span>  nr_done[indice] = <a class="code" href="classLorene_1_1Diff.html#a14adb98c5740935b2a2771ba605fe4f4">npoints</a> ;</div><div class="line"><a name="l00116"></a><span class="lineno"> 116</span>  tab[<a class="code" href="classLorene_1_1Diff.html#afbaa55cedbd7594c4fdc9c085d79ee3c">base</a>*max_points + indice] = <span class="keyword">new</span> <a class="code" href="classLorene_1_1Matrice.html">Matrice</a>(<a class="code" href="classLorene_1_1Diff.html#a14adb98c5740935b2a2771ba605fe4f4">npoints</a>, <a class="code" href="classLorene_1_1Diff.html#a14adb98c5740935b2a2771ba605fe4f4">npoints</a>) ;</div><div class="line"><a name="l00117"></a><span class="lineno"> 117</span>  <a class="code" href="classLorene_1_1Matrice.html">Matrice</a>& resu = *tab[<a class="code" href="classLorene_1_1Diff.html#afbaa55cedbd7594c4fdc9c085d79ee3c">base</a>*max_points + indice] ;</div><div class="line"><a name="l00118"></a><span class="lineno"> 118</span>  resu.<a class="code" href="classLorene_1_1Matrice.html#a06ee81317234a9dfecaede8cabbbc166">set_etat_qcq</a>() ;</div><div class="line"><a name="l00119"></a><span class="lineno"> 119</span> </div><div class="line"><a name="l00120"></a><span class="lineno"> 120</span>  <span class="keywordtype">double</span>* vect = <span class="keyword">new</span> <span class="keywordtype">double</span>[<a class="code" href="classLorene_1_1Diff.html#a14adb98c5740935b2a2771ba605fe4f4">npoints</a>] ;</div><div class="line"><a name="l00121"></a><span class="lineno"> 121</span>  <span class="keywordtype">double</span>* cres = <span class="keyword">new</span> <span class="keywordtype">double</span>[<a class="code" href="classLorene_1_1Diff.html#a14adb98c5740935b2a2771ba605fe4f4">npoints</a>] ;</div><div class="line"><a name="l00122"></a><span class="lineno"> 122</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> i=0; i<<a class="code" href="classLorene_1_1Diff.html#a14adb98c5740935b2a2771ba605fe4f4">npoints</a>; i++) {</div><div class="line"><a name="l00123"></a><span class="lineno"> 123</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> j=0; j<<a class="code" href="classLorene_1_1Diff.html#a14adb98c5740935b2a2771ba605fe4f4">npoints</a>; j++)</div><div class="line"><a name="l00124"></a><span class="lineno"> 124</span>  vect[j] = 0. ;</div><div class="line"><a name="l00125"></a><span class="lineno"> 125</span>  vect[i] = 1. ;</div><div class="line"><a name="l00126"></a><span class="lineno"> 126</span>  d2sdx2_1d(npoints, &vect, <a class="code" href="classLorene_1_1Diff.html#afbaa55cedbd7594c4fdc9c085d79ee3c">base</a> << <a class="code" href="group__def__mac.html#gab81de20c8d07e93f26854aa85f376540">TRA_R</a>) ;</div><div class="line"><a name="l00127"></a><span class="lineno"> 127</span>  <span class="keywordflow">if</span> (<a class="code" href="classLorene_1_1Diff.html#afbaa55cedbd7594c4fdc9c085d79ee3c">base</a> == <a class="code" href="group__def__mac.html#ga37e841430236cf206e8ca75a0f67b598">R_CHEBU</a>) {</div><div class="line"><a name="l00128"></a><span class="lineno"> 128</span>  mult_xm1_1d_cheb(npoints, vect, cres) ;</div><div class="line"><a name="l00129"></a><span class="lineno"> 129</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> j=0; j<<a class="code" href="classLorene_1_1Diff.html#a14adb98c5740935b2a2771ba605fe4f4">npoints</a>; j++)</div><div class="line"><a name="l00130"></a><span class="lineno"> 130</span>  resu.<a class="code" href="classLorene_1_1Matrice.html#a6a4f9fee59381b6643d36b0e8d0f0f87">set</a>(j,i) = cres[j] ;</div><div class="line"><a name="l00131"></a><span class="lineno"> 131</span>  }</div><div class="line"><a name="l00132"></a><span class="lineno"> 132</span>  <span class="keywordflow">else</span> <span class="keywordflow">if</span> (<a class="code" href="classLorene_1_1Diff.html#afbaa55cedbd7594c4fdc9c085d79ee3c">base</a> == <a class="code" href="group__def__mac.html#ga3e6b160da8c190c909181f5caaa08c75">R_JACO02</a>) {</div><div class="line"><a name="l00133"></a><span class="lineno"> 133</span>  multxpun_1d(npoints, &vect, <a class="code" href="classLorene_1_1Diff.html#afbaa55cedbd7594c4fdc9c085d79ee3c">base</a> << <a class="code" href="group__def__mac.html#gab81de20c8d07e93f26854aa85f376540">TRA_R</a>) ;</div><div class="line"><a name="l00134"></a><span class="lineno"> 134</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> j=0; j<<a class="code" href="classLorene_1_1Diff.html#a14adb98c5740935b2a2771ba605fe4f4">npoints</a>; j++)</div><div class="line"><a name="l00135"></a><span class="lineno"> 135</span>  resu.<a class="code" href="classLorene_1_1Matrice.html#a6a4f9fee59381b6643d36b0e8d0f0f87">set</a>(j,i) = vect[j] ;</div><div class="line"><a name="l00136"></a><span class="lineno"> 136</span>  }</div><div class="line"><a name="l00137"></a><span class="lineno"> 137</span>  <span class="keywordflow">else</span> {</div><div class="line"><a name="l00138"></a><span class="lineno"> 138</span>  multx_1d(npoints, &vect, <a class="code" href="classLorene_1_1Diff.html#afbaa55cedbd7594c4fdc9c085d79ee3c">base</a> << <a class="code" href="group__def__mac.html#gab81de20c8d07e93f26854aa85f376540">TRA_R</a>) ;</div><div class="line"><a name="l00139"></a><span class="lineno"> 139</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> j=0; j<<a class="code" href="classLorene_1_1Diff.html#a14adb98c5740935b2a2771ba605fe4f4">npoints</a>; j++)</div><div class="line"><a name="l00140"></a><span class="lineno"> 140</span>  resu.<a class="code" href="classLorene_1_1Matrice.html#a6a4f9fee59381b6643d36b0e8d0f0f87">set</a>(j,i) = vect[j] ;</div><div class="line"><a name="l00141"></a><span class="lineno"> 141</span>  }</div><div class="line"><a name="l00142"></a><span class="lineno"> 142</span>  }</div><div class="line"><a name="l00143"></a><span class="lineno"> 143</span>  <span class="keyword">delete</span> [] vect ;</div><div class="line"><a name="l00144"></a><span class="lineno"> 144</span>  <span class="keyword">delete</span> [] cres ;</div><div class="line"><a name="l00145"></a><span class="lineno"> 145</span>  }</div><div class="line"><a name="l00146"></a><span class="lineno"> 146</span>  <span class="keywordflow">return</span> *tab[<a class="code" href="classLorene_1_1Diff.html#afbaa55cedbd7594c4fdc9c085d79ee3c">base</a>*max_points + indice] ;</div><div class="line"><a name="l00147"></a><span class="lineno"> 147</span> }</div><div class="line"><a name="l00148"></a><span class="lineno"> 148</span> </div><div class="line"><a name="l00149"></a><span class="lineno"><a class="line" href="classLorene_1_1Diff__xdsdx2.html#a1bcb8b2db258421275dae6b1d71303da"> 149</a></span> ostream& <a class="code" href="classLorene_1_1Diff__xdsdx2.html#a1bcb8b2db258421275dae6b1d71303da">Diff_xdsdx2::operator>></a>(ostream& ost)<span class="keyword"> const </span>{</div><div class="line"><a name="l00150"></a><span class="lineno"> 150</span> </div><div class="line"><a name="l00151"></a><span class="lineno"> 151</span>  ost << <span class="stringliteral">" x d2 / dx2 "</span> << endl ;</div><div class="line"><a name="l00152"></a><span class="lineno"> 152</span> </div><div class="line"><a name="l00153"></a><span class="lineno"> 153</span>  <span class="keywordflow">return</span> ost ;</div><div class="line"><a name="l00154"></a><span class="lineno"> 154</span> </div><div class="line"><a name="l00155"></a><span class="lineno"> 155</span> }</div><div class="line"><a name="l00156"></a><span class="lineno"> 156</span> }</div><div class="ttc" id="classLorene_1_1Diff_html"><div class="ttname"><a href="classLorene_1_1Diff.html">Lorene::Diff</a></div><div class="ttdoc">Base (abstract) class for 1D spectral differential operators in one domain. </div><div class="ttdef"><b>Definition:</b> <a href="diff_8h_source.html#l00065">diff.h:65</a></div></div>
<div class="ttc" id="namespaceLorene_html"><div class="ttname"><a href="namespaceLorene.html">Lorene</a></div><div class="ttdoc">Lorene prototypes. </div><div class="ttdef"><b>Definition:</b> <a href="app__hor_8h_source.html#l00064">app_hor.h:64</a></div></div>
<div class="ttc" id="classLorene_1_1Diff_html_ad9529ab3d6813ac08e14689e34e2b9f4"><div class="ttname"><a href="classLorene_1_1Diff.html#ad9529ab3d6813ac08e14689e34e2b9f4">Lorene::Diff::max_points</a></div><div class="ttdeci">static const int max_points</div><div class="ttdoc">Maximal number of matrices stored per base. </div><div class="ttdef"><b>Definition:</b> <a href="diff_8h_source.html#l00071">diff.h:71</a></div></div>
<div class="ttc" id="classLorene_1_1Diff__xdsdx2_html_a1bcb8b2db258421275dae6b1d71303da"><div class="ttname"><a href="classLorene_1_1Diff__xdsdx2.html#a1bcb8b2db258421275dae6b1d71303da">Lorene::Diff_xdsdx2::operator>></a></div><div class="ttdeci">virtual ostream & operator>>(ostream &) const</div><div class="ttdoc">Operator >> (virtual function called by the operator <<). </div><div class="ttdef"><b>Definition:</b> <a href="diff__xdsdx2_8C_source.html#l00149">diff_xdsdx2.C:149</a></div></div>
<div class="ttc" id="classLorene_1_1Diff__xdsdx2_html_a4adfe4c96630d613d3f1e10af64d65d9"><div class="ttname"><a href="classLorene_1_1Diff__xdsdx2.html#a4adfe4c96630d613d3f1e10af64d65d9">Lorene::Diff_xdsdx2::initialize</a></div><div class="ttdeci">void initialize()</div><div class="ttdoc">Initializes arrays. </div><div class="ttdef"><b>Definition:</b> <a href="diff__xdsdx2_8C_source.html#l00079">diff_xdsdx2.C:79</a></div></div>
<div class="ttc" id="group__def__mac_html_ga3e6b160da8c190c909181f5caaa08c75"><div class="ttname"><a href="group__def__mac.html#ga3e6b160da8c190c909181f5caaa08c75">R_JACO02</a></div><div class="ttdeci">#define R_JACO02</div><div class="ttdoc">base de Jacobi(0,2) ordinaire (finjac) </div><div class="ttdef"><b>Definition:</b> <a href="type__parite_8h_source.html#l00188">type_parite.h:188</a></div></div>
<div class="ttc" id="classLorene_1_1Diff__xdsdx2_html_a97d108d41d21c6131d981683def9c77a"><div class="ttname"><a href="classLorene_1_1Diff__xdsdx2.html#a97d108d41d21c6131d981683def9c77a">Lorene::Diff_xdsdx2::~Diff_xdsdx2</a></div><div class="ttdeci">virtual ~Diff_xdsdx2()</div><div class="ttdoc">Destructor. </div><div class="ttdef"><b>Definition:</b> <a href="diff__xdsdx2_8C_source.html#l00077">diff_xdsdx2.C:77</a></div></div>
<div class="ttc" id="classLorene_1_1Diff__xdsdx2_html_a790da66fbf295a3b65bcfe69052e256d"><div class="ttname"><a href="classLorene_1_1Diff__xdsdx2.html#a790da66fbf295a3b65bcfe69052e256d">Lorene::Diff_xdsdx2::Diff_xdsdx2</a></div><div class="ttdeci">Diff_xdsdx2(int base_r, int nr)</div><div class="ttdoc">Standard constructor. </div><div class="ttdef"><b>Definition:</b> <a href="diff__xdsdx2_8C_source.html#l00067">diff_xdsdx2.C:67</a></div></div>
<div class="ttc" id="classLorene_1_1Diff_html_a14adb98c5740935b2a2771ba605fe4f4"><div class="ttname"><a href="classLorene_1_1Diff.html#a14adb98c5740935b2a2771ba605fe4f4">Lorene::Diff::npoints</a></div><div class="ttdeci">int npoints</div><div class="ttdoc">Number of coefficients. </div><div class="ttdef"><b>Definition:</b> <a href="diff_8h_source.html#l00075">diff.h:75</a></div></div>
<div class="ttc" id="group__def__mac_html_gab81de20c8d07e93f26854aa85f376540"><div class="ttname"><a href="group__def__mac.html#gab81de20c8d07e93f26854aa85f376540">TRA_R</a></div><div class="ttdeci">#define TRA_R</div><div class="ttdoc">Translation en R, used for a bitwise shift (in hex) </div><div class="ttdef"><b>Definition:</b> <a href="type__parite_8h_source.html#l00158">type_parite.h:158</a></div></div>
<div class="ttc" id="classLorene_1_1Diff__xdsdx2_html_a76832bcbcf932e997767f94fdaddb789"><div class="ttname"><a href="classLorene_1_1Diff__xdsdx2.html#a76832bcbcf932e997767f94fdaddb789">Lorene::Diff_xdsdx2::get_matrice</a></div><div class="ttdeci">virtual const Matrice & get_matrice() const</div><div class="ttdoc">Returns the matrix associated with the operator. </div><div class="ttdef"><b>Definition:</b> <a href="diff__xdsdx2_8C_source.html#l00097">diff_xdsdx2.C:97</a></div></div>
<div class="ttc" id="group__def__mac_html_ga610567430ec691ccd8f1120da031c3e0"><div class="ttname"><a href="group__def__mac.html#ga610567430ec691ccd8f1120da031c3e0">R_CHEBI</a></div><div class="ttdeci">#define R_CHEBI</div><div class="ttdoc">base de Cheb. impaire (rare) seulement </div><div class="ttdef"><b>Definition:</b> <a href="type__parite_8h_source.html#l00170">type_parite.h:170</a></div></div>
<div class="ttc" id="group__def__mac_html_gae7dd6de71b3796c971acf0b09982951c"><div class="ttname"><a href="group__def__mac.html#gae7dd6de71b3796c971acf0b09982951c">R_CHEBP</a></div><div class="ttdeci">#define R_CHEBP</div><div class="ttdoc">base de Cheb. paire (rare) seulement </div><div class="ttdef"><b>Definition:</b> <a href="type__parite_8h_source.html#l00168">type_parite.h:168</a></div></div>
<div class="ttc" id="classLorene_1_1Diff_html_afbaa55cedbd7594c4fdc9c085d79ee3c"><div class="ttname"><a href="classLorene_1_1Diff.html#afbaa55cedbd7594c4fdc9c085d79ee3c">Lorene::Diff::base</a></div><div class="ttdeci">int base</div><div class="ttdoc">Base in radial direction. </div><div class="ttdef"><b>Definition:</b> <a href="diff_8h_source.html#l00074">diff.h:74</a></div></div>
<div class="ttc" id="classLorene_1_1Matrice_html"><div class="ttname"><a href="classLorene_1_1Matrice.html">Lorene::Matrice</a></div><div class="ttdoc">Matrix handling. </div><div class="ttdef"><b>Definition:</b> <a href="matrice_8h_source.html#l00152">matrice.h:152</a></div></div>
<div class="ttc" id="classLorene_1_1Diff__xdsdx2_html_a6bd20c3714234d43330e45230f5eab9d"><div class="ttname"><a href="classLorene_1_1Diff__xdsdx2.html#a6bd20c3714234d43330e45230f5eab9d">Lorene::Diff_xdsdx2::operator=</a></div><div class="ttdeci">void operator=(const Diff_xdsdx2 &)</div><div class="ttdoc">Assignment to another Diff_xdsdx2. </div><div class="ttdef"><b>Definition:</b> <a href="diff__xdsdx2_8C_source.html#l00091">diff_xdsdx2.C:91</a></div></div>
<div class="ttc" id="classLorene_1_1Diff_html_a1768a197f28eb6c7735d3eb672e26aab"><div class="ttname"><a href="classLorene_1_1Diff.html#a1768a197f28eb6c7735d3eb672e26aab">Lorene::Diff::operator=</a></div><div class="ttdeci">void operator=(const Diff &)</div><div class="ttdoc">Assignment to another Diff. </div><div class="ttdef"><b>Definition:</b> <a href="diff_8C_source.html#l00075">diff.C:75</a></div></div>
<div class="ttc" id="classLorene_1_1Matrice_html_a6a4f9fee59381b6643d36b0e8d0f0f87"><div class="ttname"><a href="classLorene_1_1Matrice.html#a6a4f9fee59381b6643d36b0e8d0f0f87">Lorene::Matrice::set</a></div><div class="ttdeci">double & set(int j, int i)</div><div class="ttdoc">Read/write of a particuliar element. </div><div class="ttdef"><b>Definition:</b> <a href="matrice_8h_source.html#l00277">matrice.h:277</a></div></div>
<div class="ttc" id="classLorene_1_1Matrice_html_a06ee81317234a9dfecaede8cabbbc166"><div class="ttname"><a href="classLorene_1_1Matrice.html#a06ee81317234a9dfecaede8cabbbc166">Lorene::Matrice::set_etat_qcq</a></div><div class="ttdeci">void set_etat_qcq()</div><div class="ttdoc">Sets the logical state to ETATQCQ (ordinary state). </div><div class="ttdef"><b>Definition:</b> <a href="matrice_8C_source.html#l00175">matrice.C:175</a></div></div>
<div class="ttc" id="classLorene_1_1Diff__xdsdx2_html"><div class="ttname"><a href="classLorene_1_1Diff__xdsdx2.html">Lorene::Diff_xdsdx2</a></div><div class="ttdoc">Class for the elementary differential operator (see the base class Diff ). </div><div class="ttdef"><b>Definition:</b> <a href="diff_8h_source.html#l00531">diff.h:531</a></div></div>
<div class="ttc" id="group__def__mac_html_ga37e841430236cf206e8ca75a0f67b598"><div class="ttname"><a href="group__def__mac.html#ga37e841430236cf206e8ca75a0f67b598">R_CHEBU</a></div><div class="ttdeci">#define R_CHEBU</div><div class="ttdoc">base de Chebychev ordinaire (fin), dev. en 1/r </div><div class="ttdef"><b>Definition:</b> <a href="type__parite_8h_source.html#l00180">type_parite.h:180</a></div></div>
<div class="ttc" id="group__def__mac_html_ga32503ae89f3881e280be79eefc4da0d7"><div class="ttname"><a href="group__def__mac.html#ga32503ae89f3881e280be79eefc4da0d7">MAX_BASE</a></div><div class="ttdeci">#define MAX_BASE</div><div class="ttdoc">Nombre max. de bases differentes. </div><div class="ttdef"><b>Definition:</b> <a href="type__parite_8h_source.html#l00144">type_parite.h:144</a></div></div>
</div><!-- fragment --></div><!-- contents -->
</div><!-- doc-content -->
<!-- start footer part -->
<div id="nav-path" class="navpath"><!-- id is needed for treeview function! -->
<ul>
<li class="navelem"><a class="el" href="dir_2916f9fbcacdf2e5f0583aa159c707b6.html">C++</a></li><li class="navelem"><a class="el" href="dir_be81ae66a82603ee109e545ee9e7219c.html">Source</a></li><li class="navelem"><a class="el" href="dir_f57e97ea812057147ea5dc093f4832b4.html">Diff</a></li><li class="navelem"><b>diff_xdsdx2.C</b></li>
<li class="footer">Generated by
<a href="http://www.doxygen.org/index.html">
<img class="footer" src="doxygen.png" alt="doxygen"/></a> 1.8.12 </li>
</ul>
</div>
</body>
</html>
|