File: node138.html

package info (click to toggle)
babel 0.10.2-1
  • links: PTS
  • area: contrib
  • in suites: sarge
  • size: 43,932 kB
  • ctags: 29,707
  • sloc: java: 74,695; ansic: 73,142; cpp: 40,649; sh: 18,411; f90: 10,062; fortran: 6,727; python: 6,406; makefile: 3,866; xml: 118; perl: 48
file content (247 lines) | stat: -rw-r--r-- 9,149 bytes parent folder | download | duplicates (2)
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
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML 4.0 Transitional//EN">

<!--Converted with LaTeX2HTML 2002-2-1 (1.70)
original version by:  Nikos Drakos, CBLU, University of Leeds
* revised and updated by:  Marcus Hennecke, Ross Moore, Herb Swan
* with significant contributions from:
  Jens Lippmann, Marek Rouchal, Martin Wilck and others -->
<HTML>
<HEAD>
<TITLE>Accessing SIDL Arrays From FORTRAN 77</TITLE>
<META NAME="description" CONTENT="Accessing SIDL Arrays From FORTRAN 77">
<META NAME="keywords" CONTENT="users_guide">
<META NAME="resource-type" CONTENT="document">
<META NAME="distribution" CONTENT="global">

<META NAME="Generator" CONTENT="LaTeX2HTML v2002-2-1">
<META HTTP-EQUIV="Content-Style-Type" CONTENT="text/css">

<LINK REL="STYLESHEET" HREF="users_guide.css">

<LINK REL="next" HREF="node139.html">
<LINK REL="previous" HREF="node137.html">
<LINK REL="up" HREF="node131.html">
<LINK REL="next" HREF="node139.html">
</HEAD>

<BODY >

<DIV CLASS="navigation"><!--Navigation Panel-->
<A NAME="tex2html2771"
  HREF="node139.html">
<IMG WIDTH="37" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="next" SRC="next.png"></A> 
<A NAME="tex2html2765"
  HREF="node131.html">
<IMG WIDTH="26" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="up" SRC="up.png"></A> 
<A NAME="tex2html2759"
  HREF="node137.html">
<IMG WIDTH="63" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="previous" SRC="prev.png"></A> 
<A NAME="tex2html2767"
  HREF="node14.html">
<IMG WIDTH="65" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="contents" SRC="contents.png"></A> 
<A NAME="tex2html2769"
  HREF="node317.html">
<IMG WIDTH="43" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="index" SRC="index.png"></A> 
<BR>
<B> Next:</B> <A NAME="tex2html2772"
  HREF="node139.html">FORTRAN 77 objects with</A>
<B> Up:</B> <A NAME="tex2html2766"
  HREF="node131.html">FORTRAN 77 Bindings</A>
<B> Previous:</B> <A NAME="tex2html2760"
  HREF="node137.html">Implementing Classes in FORTRAN</A>
 &nbsp; <B>  <A NAME="tex2html2768"
  HREF="node14.html">Contents</A></B> 
 &nbsp; <B>  <A NAME="tex2html2770"
  HREF="node317.html">Index</A></B> 
<BR>
<BR></DIV>
<!--End of Navigation Panel-->

<H1><A NAME="SECTION03370000000000000000"></A>
<A NAME="7324"></A>
<BR>
Accessing SIDL Arrays From FORTRAN 77
</H1>

<P>
For FORTRAN 77, the difference in how you access normal SIDL arrays and
r-arrays is profound. Normal SIDL arrays are passed in as an
<TT>integer*8</TT>, and you either access them using an function API or
by converting the array data to a index into a known array. R-arrays
appear like normal FORTRAN 77 arrays, so there is a big incentive to
use r-arrays unless you cannot.

<P>
The client-side interface for the <TT>solve</TT> example introduced in
Section&nbsp;<A HREF="node60.html#ss:r-arrays">5.4</A> behaves as if it is a FORTRAN 77 function
with the following declarations:

<P>
<BR>
<PRE  CLASS="verbatim">        subroutine num_Linsol_solve_f(self, A, x, b, m, n)
        implicit none
C       in num.Linsol self
        integer*8 self
        integer*4 m, n
C       in rarray&lt;double,2&gt; A(m,n)
        double precision A(0:m-1, 0:n-1)
C       inout rarray&lt;double&gt; x(n)
        double precision x(0:n-1)
C       in rarray&lt;double&gt; b(m)
        double precision b(0:m-1)
        end
</PRE></td></tr></table></blockquote>
<P>
FORTRAN 77 programmers should note that the array indices go from 0 to
<TT>m</TT><SPAN CLASS="MATH"><IMG
 WIDTH="23" HEIGHT="29" ALIGN="MIDDLE" BORDER="0"
 SRC="img24.png"
 ALT="$-1$"></SPAN> instead of the normal 1 to <TT>m</TT>. This is a concession
to the C and C++ programmers who have to deal with the fact that A is
stored in column-major order.

<P>
On the server-side, the interface for <TT>solve</TT> appears as follows:

<P>
<BR>
<PRE  CLASS="verbatim">        subroutine num_Linsol_solve_fi(self, A, x, b, m, n)
        implicit none
C       in num.Linsol self
        integer*8 self
C       in int m
        integer*4 m
C       in int n
        integer*4 n
C       in rarray&lt;double,2&gt; A(m,n)
        double precision A(0:m-1, 0:n-1)
C       inout rarray&lt;double&gt; x(n)
        double precision x(0:n-1)
C       in rarray&lt;double&gt; b(m)
        double precision b(0:m-1)

C       DO-NOT-DELETE splicer.begin(num.Linsol.solve)
C       Insert the implementation here...
C       DO-NOT-DELETE splicer.end(num.Linsol.solve)
        end
</PRE></td></tr></table></blockquote>
<P>
Note again that the array indices go from 0 to <TT>m</TT><SPAN CLASS="MATH"><IMG
 WIDTH="23" HEIGHT="29" ALIGN="MIDDLE" BORDER="0"
 SRC="img24.png"
 ALT="$-1$"></SPAN>. The
implementation should avoid changing the data in <TT>in</TT>
parameters.

<P>
The remainder of this section is dedicated to how you access normal
SIDL arrays. The normal SIDL C function API is available from FORTRAN
77 to create, destroy and access array elements and meta-data.  The
function name for FORTRAN has <TT>_f</TT> appended.

<P>
For SIDL types dcomplex, double, fcomplex , float, int and long, 
SIDL provides a method to get direct access to the array elements. 
For the other types, you must use the functional API to access array elements.

<P>
For type X, there is a FORTRAN 77 function called <TT>sidl_X__array_access_f</TT> to 
provide a method to get direct access. 
An example is given below. Of course, this will
not work if your FORTRAN 77 compiler does array bounds checking.

<P>
<BR>
<PRE  CLASS="verbatim">        integer*4 lower(1), upper(1), stride(1), i, index(1)
        integer*4 value, refindex, refarray(1), modval
        integer*8 nextprime, tmp
        lower(1) = 0
        value = 0
        upper(1) = len - 1
        call sidl_int__array_create_f(1, lower, upper, retval)
        call sidl_int__array_access_f(retval, refarray, lower,
     $       upper, stride, refindex)
        do i = 0, len - 1
           tmp = value
           value = nextprime(tmp)
           modval = mod(i, 3)
           if (modval .eq. 0) then
              call sidl_int__array_set1_f(retval, i, value)
           else
              if (modval .eq. 1) then
                 index(1) = i
                 call sidl_int__array_set_f(retval, index, value)
              else
C this is equivalent to the sidl_int__array_set_f(retval, index, value)
                 refarray(refindex + stride(1)*(i - lower(1))) =
     $                value
              endif
           endif
        enddo
</PRE></td></tr></table></blockquote>
<P>
To access a two dimensional array, the expression referring to element i, j is 

<P>
<BR>
<PRE  CLASS="verbatim">      refarray(refindex + stride(1) * (i - lower(1)) + stride(2) * (j - lower(2))    
</PRE></td></tr></table></blockquote>
<P>
To access a three dimensional array, the expression referring to element i, j k is 

<P>
<BR>
<PRE  CLASS="verbatim">      refarray(refindex + stride(1) * (i - lower(1)) + stride(2) * (j - lower(2))    
</PRE></td></tr></table></blockquote>
<P>
You can call things like LINPACK or BLAS if you want, 
but you should check the stride to make sure the array 
is packed as you need it to be. 
<TT>stride(i)</TT> indicates the distance between elements in dimension <TT>i</TT>. 
A value of 1 means elements are packed densely in dimension <TT>i</TT>. 
Negative stride values are possible, and when an array is
a slice of another array, there may be no dimension with a stride of 1. 

<P>
For a <TT><I CLASS="slanted">dcomplex</I></TT> array, the reference array should a FORTRAN array of 
<TT>REAL*8</TT> instead of a FORTRAN array of double complex to avoid potential 
alignment problems. For a <TT><I CLASS="slanted">fcomplex</I></TT> array, the reference array is a 
<TT>COMPLEX*8</TT> because we don't anticipate an alignment problem in this
case.

<P>

<DIV CLASS="navigation"><HR>
<!--Navigation Panel-->
<A NAME="tex2html2771"
  HREF="node139.html">
<IMG WIDTH="37" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="next" SRC="next.png"></A> 
<A NAME="tex2html2765"
  HREF="node131.html">
<IMG WIDTH="26" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="up" SRC="up.png"></A> 
<A NAME="tex2html2759"
  HREF="node137.html">
<IMG WIDTH="63" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="previous" SRC="prev.png"></A> 
<A NAME="tex2html2767"
  HREF="node14.html">
<IMG WIDTH="65" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="contents" SRC="contents.png"></A> 
<A NAME="tex2html2769"
  HREF="node317.html">
<IMG WIDTH="43" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="index" SRC="index.png"></A> 
<BR>
<B> Next:</B> <A NAME="tex2html2772"
  HREF="node139.html">FORTRAN 77 objects with</A>
<B> Up:</B> <A NAME="tex2html2766"
  HREF="node131.html">FORTRAN 77 Bindings</A>
<B> Previous:</B> <A NAME="tex2html2760"
  HREF="node137.html">Implementing Classes in FORTRAN</A>
 &nbsp; <B>  <A NAME="tex2html2768"
  HREF="node14.html">Contents</A></B> 
 &nbsp; <B>  <A NAME="tex2html2770"
  HREF="node317.html">Index</A></B> </DIV>
<!--End of Navigation Panel-->
<ADDRESS>
<br><br>babel-0.10.2<br>users_guide Last Modified 2005-03-23<br><br><a href="http://www.llnl.gov/CASC/components">http://www.llnl.gov/CASC/components</a><br><a href="mailto:components@llnl.gov">components@llnl.gov</a>
</ADDRESS>
</BODY>
</HTML>