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 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293
|
<!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>Function: first</TITLE>
<META NAME="description" CONTENT="Function: first">
<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="previous" HREF="node96.html">
<LINK REL="up" HREF="node65.html">
<LINK REL="next" HREF="node98.html">
</HEAD>
<BODY >
<DIV CLASS="navigation"><!--Navigation Panel-->
<A NAME="tex2html2100"
HREF="node98.html">
<IMG WIDTH="37" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="next" SRC="next.png"></A>
<A NAME="tex2html2094"
HREF="node65.html">
<IMG WIDTH="26" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="up" SRC="up.png"></A>
<A NAME="tex2html2090"
HREF="node96.html">
<IMG WIDTH="63" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="previous" SRC="prev.png"></A>
<A NAME="tex2html2096"
HREF="node14.html">
<IMG WIDTH="65" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="contents" SRC="contents.png"></A>
<A NAME="tex2html2098"
HREF="node317.html">
<IMG WIDTH="43" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="index" SRC="index.png"></A>
<BR>
<B> Next:</B> <A NAME="tex2html2101"
HREF="node98.html">The C Macro API</A>
<B> Up:</B> <A NAME="tex2html2095"
HREF="node65.html">The Array API</A>
<B> Previous:</B> <A NAME="tex2html2091"
HREF="node96.html">Function: ensure</A>
<B> <A NAME="tex2html2097"
HREF="node14.html">Contents</A></B>
<B> <A NAME="tex2html2099"
HREF="node317.html">Index</A></B>
<BR>
<BR></DIV>
<!--End of Navigation Panel-->
<H3><A NAME="SECTION024463200000000000000"></A><A NAME="3030"></A>
<BR>
Function: first
</H3>
<P>
<BR>
<PRE CLASS="verbatim">/* C */
double *
sidl_double__array_first(const struct sidl_double__array *src);
// C++
double* first() throw();
C FORTRAN 77
subroutine sidl_double__array_access_f(array, ref, lower, upper,
$ stride, index)
integer*8 array
integer*4 lower(), upper(), stride(), index
integer*4 ref()
</PRE></td></tr></table></blockquote>
<P>
This method provides direct access to the element data. Using this
pointer and the stride information, you can perform your own array
accesses without function calls. This method isn't available for
arrays of strings, interface and objects because of memory/reference
management issues. There is no equivalent of this of this
function in Java or Python. To see how to get direct array access in
FORTRAN 90, see Chapter <A HREF="node140.html#c:f90">9</A>.<A NAME="3033"></A><A NAME="3034"></A>
<P>
The FORTRAN versions of the method return the lower, upper and
stride information in three arrays, each with enough elements to hold
an entry for each dimension of <TT>array</TT>. Because FORTRAN 77 does
not have pointers, you must pass in a reference array, <TT>array</TT>.
Upon exit, <TT>ref(index)</TT> is the first element of the array. The
type of <TT>ref</TT> depends on the type of the array.
<P>
<BLOCKQUOTE>
WARNING:<A NAME="3040"></A><A NAME="3041"></A>
While calling the FORTRAN direct access routines, there is a
possibility of an alignment error between your reference pointer,
<TT>ref</TT>, and the pointer to the first element of the array data.
The problem is more likely with arrays of <TT>double</TT> or
<TT>dcomplex</TT>; although, it could occur with any type on some future
platform. If <TT>index</TT> is zero on return, an alignment error
occurred. If an alignment error occurs, you may be able to solve it
by recompiling your FORTRAN files with flags to force doubles to be
aligned on 8 byte boundaries. For example, the <TT>-malign-double</TT>
flag for g77 forces doubles to be aligned on 64-bit boundaries. An
alignment error occurs when <TT>(char *)ref</TT> minus <TT>(char
*)sidl_double__array_first(array)</TT> is not integer divisible by
<TT>sizeof(datatype)</TT> where <TT>ref</TT> refers to the address of the
reference array.
</BLOCKQUOTE>
<BR><P>
Here is an example FORTRAN 77 subroutine to output each element of
a <A NAME="3052"></A>
1-dimensional array of doubles using the direct access
routine. FORTRAN 90 has a pointer in the array derived type when
direct access is possible.
<P>
<BR>
<PRE CLASS="verbatim">C This subroutine will print each element of an array of doubles
subroutine print_array(dblarray)
implicit none
integer*8 dblarray
real*8 refarray(1)
integer*4 lower(1), upper(1), stride(1), index, dimen, i
if (dblarray .ne. 0) then
call sidl_double__array_dimen_f(dblarray, dimen)
if (dimen .eq. 1) then
call sidl_double__array_access_f(dblarray, refarray,
$ lower, upper, stride, index)
if (index .ne. 0) then
do i = lower(1), upper(1)
write(*,*) refarray(index + (i-lower(1))*stride(1))
enddo
else
write(*,*) 'Alignment error occured'
endif
endif
endif
end
</PRE></td></tr></table></blockquote>
<P>
For a 2-dimensional array, the loop and array access is
<P>
<BR>
<PRE CLASS="verbatim"> do i = lower(1), upper(1)
do j = lower(2), upper(2)
write(*,*) refarray(index+(i-lower(1))*stride(1)+
$ (j - lower(2))*stride(2))
enddo
enddo
</PRE></td></tr></table></blockquote>
<P>
Suppose you are wrapping a legacy FORTRAN application and you need to
pass a SIDL array to a FORTRAN subroutine. Further suppose there is
a FORTRAN 77 and FORTRAN 90 version of the subroutine. For example,
the FORTRAN 77 subroutine has a signature such as:
<P>
<BR>
<PRE CLASS="verbatim"> subroutine TriedAndTrue(x, n)
integer n
real*8 x(n)
C insert wonderful, efficient, debugged code here
end
</PRE></td></tr></table></blockquote>
<P>
The FORTRAN 90 subroutine has basically the same signature as follows:
<A NAME="3056"></A>
<BR>
<PRE CLASS="verbatim">subroutine TriedAndTrue(x, n)
integer (selected_int_kind(9)) :: n
real (selected_real_kind(17, 308)) :: x(n)
! insert wonderful, efficient, debugged code here
end subroutine TriedAndTrue
</PRE></td></tr></table></blockquote>
<P>
Here is one way to wrap this method using SIDL. First of all, the
SIDL method definition specifies that the array must be a
1-dimensional, column-major ordered array. This forces the incoming
array to be a dense column.
<P>
<BR>
<PRE CLASS="verbatim"> static void TriedAndTrue(inout array<double,1,column-major> arg);
</PRE></td></tr></table></blockquote>
<P>
Given that method definition in a class named Class and a package
named Pkg, the implementation of the wrapper should look something
like the following for FORTRAN 77:
<P>
<BR>
<PRE CLASS="verbatim"> subroutine Pkg_Class_TriedAndTrue_fi(arg)
implicit none
integer*8 arg
C DO-NOT-DELETE splicer.begin(Pkg.Class.TriedAndTrue)
real*8 refarray(1)
integer*4 lower(1), upper(1), stride(1), index
integer n
call sidl_double__array_access_f(arg, refarray,
$ lower, upper, stride, index)
if (index .ne. 0) then
c we can assume stride(1) = 1 because of column-major specification
n = 1 + upper(1) - lower(1)
call TriedAndTrue(refarray(index), n)
else
write(*,*) 'ERROR: array alignment'
endif
C DO-NOT-DELETE splicer.end(Pkg.Class.TriedAndTrue)
end
</PRE></td></tr></table></blockquote>
<P>
Similarly, it should look something like the following for FORTRAN 90, where
the include statements are required at the top of the Impl file to ensure
proper handling of subroutine names that have automatically been mangled
by the Babel compiler:
<P>
<BR>
<PRE CLASS="verbatim">#include "Pkg_Class_fAbbrev.h"
#include "sidl_BaseClass_fAbbrev.h"
#include "sidl_BaseInterface_fAbbrev.h"
! DO-NOT-DELETE splicer.begin(_miscellaneous_code_start)
#include "sidl_double_fAbbrev.h"
! DO-NOT-DELETE splicer.end(_miscellaneous_code_start)
.
.
.
subroutine Pkg_Class_TriedAndTrue_mi(arg)
! DO-NOT-DELETE splicer.begin(Pkg.Class.TriedAndTrue.use)
use SIDL_double_array
! DO-NOT-DELETE splicer.end(Pkg.Class.TriedAndTrue.use)
implicit none
type(sidl_double_a) :: arg
! DO-NOT-DELETE splicer.begin(Pkg.Class.TriedAndTrue)
real (selected_real_kind(17,308)), dimension(1) :: refarray
integer (selected_int_kind(8)), dimension(1) :: low, up, str
integer (selected_int_kind(8)) :: index, n
call access(arg, refarray, low, up, str, index)
if (index .ne. 0) then
! We can assume stride(1) = 1 because of column-major specification
n = 1 + upper(1) - lower(1)
call TriedAndTrue(refarray(index), n)
else
write(*,*) 'ERROR: array alignment'
endif
! DO-NOT-DELETE splicer.end(Pkg.Class.TriedAndTrue)
end subroutine Pkg_Class_TriedAndTrue_mi
</PRE></td></tr></table></blockquote>
<P>
<DIV CLASS="navigation"><HR>
<!--Navigation Panel-->
<A NAME="tex2html2100"
HREF="node98.html">
<IMG WIDTH="37" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="next" SRC="next.png"></A>
<A NAME="tex2html2094"
HREF="node65.html">
<IMG WIDTH="26" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="up" SRC="up.png"></A>
<A NAME="tex2html2090"
HREF="node96.html">
<IMG WIDTH="63" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="previous" SRC="prev.png"></A>
<A NAME="tex2html2096"
HREF="node14.html">
<IMG WIDTH="65" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="contents" SRC="contents.png"></A>
<A NAME="tex2html2098"
HREF="node317.html">
<IMG WIDTH="43" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="index" SRC="index.png"></A>
<BR>
<B> Next:</B> <A NAME="tex2html2101"
HREF="node98.html">The C Macro API</A>
<B> Up:</B> <A NAME="tex2html2095"
HREF="node65.html">The Array API</A>
<B> Previous:</B> <A NAME="tex2html2091"
HREF="node96.html">Function: ensure</A>
<B> <A NAME="tex2html2097"
HREF="node14.html">Contents</A></B>
<B> <A NAME="tex2html2099"
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>
|