1 | !------------------------------------------------------------------------- |
---|
2 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
3 | !------------------------------------------------------------------------- |
---|
4 | ! CVS m_inpak90.F90,v 1.8 2006-12-19 00:22:35 jacob Exp |
---|
5 | ! CVS MCT_2_8_0 |
---|
6 | !------------------------------------------------------------------------- |
---|
7 | !BOI |
---|
8 | ! |
---|
9 | ! !TITLE: Inpak 90 Documentation \\ Version 1.01 |
---|
10 | ! |
---|
11 | ! !AUTHORS: Arlindo da Silva |
---|
12 | ! |
---|
13 | ! !AFFILIATION: Data Assimilation Office, NASA/GSFC, Greenbelt, MD 20771 |
---|
14 | ! |
---|
15 | ! !DATE: June 20, 1996 |
---|
16 | ! |
---|
17 | ! !INTRODUCTION: Package Overview |
---|
18 | ! |
---|
19 | ! Inpak 90 is a Fortran (77/90) collection of |
---|
20 | ! routines/functions for accessing {\em Resource Files} |
---|
21 | ! in ASCII format. The package is optimized |
---|
22 | ! for minimizing formatted I/O, performing all of its string |
---|
23 | ! operations in memory using Fortran intrinsic functions. |
---|
24 | ! |
---|
25 | ! \subsection{Resource Files} |
---|
26 | ! |
---|
27 | ! A {\em Resource File} is a text file consisting of variable |
---|
28 | ! length lines (records), each possibly starting with a {\em label} |
---|
29 | ! (or {\em key}), followed by some data. A simple resource file |
---|
30 | ! looks like this: |
---|
31 | ! |
---|
32 | ! \begin{verbatim} |
---|
33 | ! # Lines starting with # are comments which are |
---|
34 | ! # ignored during processing. |
---|
35 | ! my_file_names: jan87.dat jan88.dat jan89.dat |
---|
36 | ! radius_of_the_earth: 6.37E6 # these are comments too |
---|
37 | ! constants: 3.1415 25 |
---|
38 | ! my_favourite_colors: green blue 022 # text & number are OK |
---|
39 | ! \end{verbatim} |
---|
40 | ! |
---|
41 | ! In this example, {\tt my\_file\_names:} and {\tt constants:} |
---|
42 | ! are labels, while {\tt jan87.dat, jan88.dat} and {\tt jan89.dat} are |
---|
43 | ! data associated with label {\tt my\_file\_names:}. |
---|
44 | ! Resource files can also contain simple tables of the form, |
---|
45 | ! |
---|
46 | ! \begin{verbatim} |
---|
47 | ! my_table_name:: |
---|
48 | ! 1000 3000 263.0 |
---|
49 | ! 925 3000 263.0 |
---|
50 | ! 850 3000 263.0 |
---|
51 | ! 700 3000 269.0 |
---|
52 | ! 500 3000 287.0 |
---|
53 | ! 400 3000 295.8 |
---|
54 | ! 300 3000 295.8 |
---|
55 | ! :: |
---|
56 | ! \end{verbatim} |
---|
57 | ! |
---|
58 | ! Resource files are random access, the particular order of the |
---|
59 | ! records are not important (except between ::s in a table definition). |
---|
60 | ! |
---|
61 | ! \subsection{A Quick Stroll} |
---|
62 | ! |
---|
63 | ! The first step is to load the ASCII resource (rc) file into |
---|
64 | ! memory\footnote{See next section for a complete description |
---|
65 | ! of parameters for each routine/function}: |
---|
66 | ! |
---|
67 | ! \begin{verbatim} |
---|
68 | ! call i90_LoadF ( 'my_file.rc', iret ) |
---|
69 | ! \end{verbatim} |
---|
70 | ! |
---|
71 | ! The next step is to select the label (record) of interest, say |
---|
72 | ! |
---|
73 | ! \begin{verbatim} |
---|
74 | ! call i90_label ( 'constants:', iret ) |
---|
75 | ! \end{verbatim} |
---|
76 | ! |
---|
77 | ! The 2 constants above can be retrieved with the following code |
---|
78 | ! fragment: |
---|
79 | ! \begin{verbatim} |
---|
80 | ! real r |
---|
81 | ! integer i |
---|
82 | ! call i90_label ( 'constants:', iret ) |
---|
83 | ! r = i90_gfloat(iret) ! results in r = 3.1415 |
---|
84 | ! i = i90_gint(iret) ! results in i = 25 |
---|
85 | ! \end{verbatim} |
---|
86 | ! |
---|
87 | ! The file names above can be retrieved with the following |
---|
88 | ! code fragment: |
---|
89 | ! \begin{verbatim} |
---|
90 | ! character*20 fn1, fn2, fn3 |
---|
91 | ! integer iret |
---|
92 | ! call i90_label ( 'my_file_names:', iret ) |
---|
93 | ! call i90_Gtoken ( fn1, iret ) ! ==> fn1 = 'jan87.dat' |
---|
94 | ! call i90_Gtoken ( fn2, iret ) ! ==> fn1 = 'jan88.dat' |
---|
95 | ! call i90_Gtoken ( fn3, iret ) ! ==> fn1 = 'jan89.dat' |
---|
96 | ! \end{verbatim} |
---|
97 | ! |
---|
98 | ! To access the table above, the user first must use {\tt i90\_label()} to |
---|
99 | ! locate the beginning of the table, e.g., |
---|
100 | ! |
---|
101 | ! \begin{verbatim} |
---|
102 | ! call i90_label ( 'my_table_name::', iret ) |
---|
103 | ! \end{verbatim} |
---|
104 | ! |
---|
105 | ! Subsequently, {\tt i90\_gline()} can be used to gain access to each |
---|
106 | ! row of the table. Here is a code fragment to read the above |
---|
107 | ! table (7 rows, 3 columns): |
---|
108 | ! |
---|
109 | ! \begin{verbatim} |
---|
110 | ! real table(7,3) |
---|
111 | ! character*20 word |
---|
112 | ! integer iret |
---|
113 | ! call i90_label ( 'my_table_name::', iret ) |
---|
114 | ! do i = 1, 7 |
---|
115 | ! call i90_gline ( iret ) |
---|
116 | ! do j = 1, 3 |
---|
117 | ! table(i,j) = i90_gfloat ( iret ) |
---|
118 | ! end do |
---|
119 | ! end do |
---|
120 | ! \end{verbatim} |
---|
121 | ! |
---|
122 | ! Get the idea? |
---|
123 | ! |
---|
124 | ! \newpage |
---|
125 | ! \subsection{Main Routine/Functions} |
---|
126 | ! |
---|
127 | ! \begin{verbatim} |
---|
128 | ! ------------------------------------------------------------------ |
---|
129 | ! Routine/Function Description |
---|
130 | ! ------------------------------------------------------------------ |
---|
131 | ! I90_LoadF ( filen, iret ) loads resource file into memory |
---|
132 | ! I90_Label ( label, iret ) selects a label (key) |
---|
133 | ! I90_GLine ( iret ) selects next line (for tables) |
---|
134 | ! I90_Gtoken ( word, iret ) get next token |
---|
135 | ! I90_Gfloat ( iret ) returns next float number (function) |
---|
136 | ! I90_GInt ( iret ) returns next integer number (function) |
---|
137 | ! i90_AtoF ( string, iret ) ASCII to float (function) |
---|
138 | ! i90_AtoI ( string, iret ) ASCII to integer (function) |
---|
139 | ! I90_Len ( string ) string length without trailing blanks |
---|
140 | ! LabLin ( label ) similar to i90_label (no iret) |
---|
141 | ! FltGet ( default ) returns next float number (function) |
---|
142 | ! IntGet ( default ) returns next integer number (function) |
---|
143 | ! ChrGet ( default ) returns next character (function) |
---|
144 | ! TokGet ( string, default ) get next token |
---|
145 | ! ------------------------------------------------------------------ |
---|
146 | ! \end{verbatim} |
---|
147 | ! |
---|
148 | ! {\em Common Arguments:} |
---|
149 | ! |
---|
150 | ! \begin{verbatim} |
---|
151 | ! character*(*) filen file name |
---|
152 | ! integer iret error return code (0 is OK) |
---|
153 | ! character*(*) label label (key) to locate record |
---|
154 | ! character*(*) word blank delimited string |
---|
155 | ! character*(*) string a sequence of characters |
---|
156 | ! \end{verbatim} |
---|
157 | ! |
---|
158 | ! See the Prologues in the next section for additional details. |
---|
159 | ! |
---|
160 | ! |
---|
161 | ! \subsection{Package History} |
---|
162 | ! Back in the 70s Eli Isaacson wrote IOPACK in Fortran |
---|
163 | ! 66. In June of 1987 I wrote Inpak77 using |
---|
164 | ! Fortran 77 string functions; Inpak 77 is a vastly |
---|
165 | ! simplified IOPACK, but has its own goodies not found in |
---|
166 | ! IOPACK. Inpak 90 removes some obsolete functionality in |
---|
167 | ! Inpak77, and parses the whole resource file in memory for |
---|
168 | ! performance. Despite its name, Inpak 90 compiles fine |
---|
169 | ! under any modern Fortran 77 compiler. |
---|
170 | ! |
---|
171 | ! \subsection{Bugs} |
---|
172 | ! Inpak 90 is not very gracious with error messages. |
---|
173 | ! The interactive functionality of Inpak77 has not been implemented. |
---|
174 | ! The comment character \# cannot be escaped. |
---|
175 | ! |
---|
176 | ! \subsection{Availability} |
---|
177 | ! |
---|
178 | ! This software is available at |
---|
179 | ! \begin{verbatim} |
---|
180 | ! ftp://niteroi.gsfc.nasa.gov/pub/packages/i90/ |
---|
181 | ! \end{verbatim} |
---|
182 | ! There you will find the following files: |
---|
183 | ! \begin{verbatim} |
---|
184 | ! i90.f Fortran 77/90 source code |
---|
185 | ! i90.h Include file needed by i90.f |
---|
186 | ! ti90.f Test code |
---|
187 | ! i90.ps Postscript documentation |
---|
188 | ! \end{verbatim} |
---|
189 | ! An on-line version of this document is available at |
---|
190 | ! \begin{verbatim} |
---|
191 | ! ftp://niteroi.gsfc.nasa.gov/www/packages/i90/i90.html |
---|
192 | ! \end{verbatim} |
---|
193 | ! |
---|
194 | !EOI |
---|
195 | !------------------------------------------------------------------------- |
---|
196 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
197 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
198 | !----------------------------------------------------------------------- |
---|
199 | ! |
---|
200 | ! !REVISION HISTORY: |
---|
201 | ! 03Jul96 - J. Guo - evolved to Fortran 90 module. The |
---|
202 | ! modifications include 1) additional subroutines to |
---|
203 | ! dynamically manage the memory, 2) privatized most |
---|
204 | ! entries, 3) included "i90.h" into the module source |
---|
205 | ! with better initializations, 4) removed blockdata, 5) |
---|
206 | ! used a portable opntext() call to avoid I/O portability |
---|
207 | ! problems. |
---|
208 | ! |
---|
209 | ! See I90_page() I90_Release(), and I90_LoadF() for |
---|
210 | ! details. |
---|
211 | ! |
---|
212 | ! 05Aug98 - Jing Guo - |
---|
213 | ! Removed i90_page() and its references. |
---|
214 | ! Added internal subroutines push_() and pop_(). |
---|
215 | ! Modified i90_release(). |
---|
216 | ! Added i90_fullrelease(). |
---|
217 | ! Removed %loaded. Check i90_depth instead. |
---|
218 | ! 06Aug98 - Todling - made I90_gstr public |
---|
219 | ! 20Dec98 - Jing Guo - replaced the description of I90_Gstr |
---|
220 | ! 28Sep99 - Jing Guo - Merged with the MPI version with |
---|
221 | ! some addtional changes based on |
---|
222 | ! merging decisions. |
---|
223 | ! 12Oct99 - Larson/Guo - Overloaded fltget() to new routines |
---|
224 | ! getfltsp() and fltgetdp(), providing better support |
---|
225 | ! for 32 and 64 bit platforms, respectively. |
---|
226 | !_______________________________________________________________________ |
---|
227 | |
---|
228 | module m_inpak90 |
---|
229 | use m_stdio, only : stderr,stdout |
---|
230 | use m_realkinds, only: FP, SP, DP,kind_r8 |
---|
231 | implicit none |
---|
232 | private |
---|
233 | public :: I90_LoadF ! loads a resource file into memory |
---|
234 | public :: I90_allLoadF! loads/populates a resource file to all PEs |
---|
235 | public :: I90_Release ! Releases one cached resource file |
---|
236 | public :: I90_fullRelease ! Releases the whole stack |
---|
237 | public :: I90_Label ! selects a label (key) |
---|
238 | public :: I90_GLine ! selects the next line (for tables) |
---|
239 | public :: I90_Gtoken ! gets the next token |
---|
240 | public :: I90_Gstr ! get a string upto to a "$" or EOL |
---|
241 | |
---|
242 | public :: I90_AtoF ! ASCII to float (function) |
---|
243 | public :: I90_AtoI ! ASCII to integer (function) |
---|
244 | |
---|
245 | public :: I90_Gfloat ! returns next float number (function) |
---|
246 | public :: I90_GInt ! returns next integer number (function) |
---|
247 | |
---|
248 | public :: lablin,rdnext,fltget,intget,getwrd,str2rn,chrget,getstr |
---|
249 | public :: strget |
---|
250 | |
---|
251 | interface fltget; module procedure & |
---|
252 | fltgetsp, & |
---|
253 | fltgetdp |
---|
254 | end interface |
---|
255 | |
---|
256 | |
---|
257 | !----------------------------------------------------------------------- |
---|
258 | ! |
---|
259 | ! This part was originally in "i90.h", but included for module. |
---|
260 | ! |
---|
261 | |
---|
262 | ! revised parameter table to fit Fortran 90 standard |
---|
263 | |
---|
264 | integer, parameter :: LSZ = 256 |
---|
265 | |
---|
266 | !ams |
---|
267 | ! On Linux with the Fujitsu compiler, I needed to reduce NBUF_MAX |
---|
268 | !ams |
---|
269 | ! integer, parameter :: NBUF_MAX = 400*(LSZ) ! max size of buffer |
---|
270 | ! integer, parameter :: NBUF_MAX = 200*(LSZ) ! max size of buffer |
---|
271 | ! Further reduction of NBUF_MAX was necessary for the Fujitsu VPP: |
---|
272 | integer, parameter :: NBUF_MAX = 128*(LSZ)-1 ! Maximum buffer size |
---|
273 | ! that works with the |
---|
274 | ! Fujitsu-VPP platform. |
---|
275 | |
---|
276 | |
---|
277 | character, parameter :: BLK = achar(32) ! blank (space) |
---|
278 | character, parameter :: TAB = achar(09) ! TAB |
---|
279 | character, parameter :: EOL = achar(10) ! end of line mark (newline) |
---|
280 | character, parameter :: EOB = achar(00) ! end of buffer mark (null) |
---|
281 | character, parameter :: NULL= achar(00) ! what it says |
---|
282 | |
---|
283 | type inpak90 |
---|
284 | ! May be easily paged for extentable file size (J.G.) |
---|
285 | |
---|
286 | integer :: nbuf ! actual size of buffer |
---|
287 | character(len=NBUF_MAX),pointer :: buffer ! hold the whole file? |
---|
288 | character(len=LSZ), pointer :: this_line ! the current line |
---|
289 | |
---|
290 | integer :: next_line ! index for next line on buffer |
---|
291 | |
---|
292 | type(inpak90),pointer :: last |
---|
293 | end type inpak90 |
---|
294 | |
---|
295 | integer,parameter :: MALLSIZE_=10 ! just an estimation |
---|
296 | |
---|
297 | character(len=*),parameter :: myname='MCT(MPEU)::m_inpak90' |
---|
298 | !----------------------------------------------------------------------- |
---|
299 | |
---|
300 | integer,parameter :: i90_MXDEP = 4 |
---|
301 | integer,save :: i90_depth = 0 |
---|
302 | type(inpak90),save,pointer :: i90_now |
---|
303 | |
---|
304 | !----------------------------------------------------------------------- |
---|
305 | contains |
---|
306 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
307 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
308 | !BOP ------------------------------------------------------------------- |
---|
309 | ! |
---|
310 | ! !IROUTINE: I90_allLoadF - populate a rooted database to all PEs |
---|
311 | ! |
---|
312 | ! !DESCRIPTION: |
---|
313 | ! |
---|
314 | ! !INTERFACE: |
---|
315 | |
---|
316 | subroutine I90_allLoadF(fname,root,comm,istat) |
---|
317 | use m_mpif90, only : MP_perr |
---|
318 | use m_mpif90, only : MP_comm_rank |
---|
319 | use m_mpif90, only : MP_CHARACTER |
---|
320 | use m_mpif90, only : MP_INTEGER |
---|
321 | use m_die, only : perr |
---|
322 | implicit none |
---|
323 | character(len=*),intent(in) :: fname |
---|
324 | integer,intent(in) :: root |
---|
325 | integer,intent(in) :: comm |
---|
326 | integer,intent(out) :: istat |
---|
327 | |
---|
328 | ! !REVISION HISTORY: |
---|
329 | ! 28Jul98 - Jing Guo <guo@thunder> - initial prototype/prolog/code |
---|
330 | !EOP ___________________________________________________________________ |
---|
331 | |
---|
332 | character(len=*),parameter :: myname_=myname//'::I90_allLoadF' |
---|
333 | integer :: myID,ier |
---|
334 | |
---|
335 | istat=0 |
---|
336 | |
---|
337 | call MP_comm_rank(comm,myID,ier) |
---|
338 | if(ier/=0) then |
---|
339 | call MP_perr(myname_,'MP_comm_rank()',ier) |
---|
340 | istat=ier |
---|
341 | return |
---|
342 | endif |
---|
343 | |
---|
344 | if(myID == root) then |
---|
345 | call i90_LoadF(fname,ier) |
---|
346 | if(ier /= 0) then |
---|
347 | call perr(myname_,'i90_LoadF("//trim(fname)//")',ier) |
---|
348 | istat=ier |
---|
349 | return |
---|
350 | endif |
---|
351 | else |
---|
352 | call push_(ier) |
---|
353 | if(ier /= 0) then |
---|
354 | call perr(myname_,'push_()',ier) |
---|
355 | istat=ier |
---|
356 | return |
---|
357 | endif |
---|
358 | endif |
---|
359 | |
---|
360 | ! Initialize the buffer on all PEs |
---|
361 | |
---|
362 | call MPI_Bcast(i90_now%buffer,NBUF_MAX,MP_CHARACTER,root,comm,ier) |
---|
363 | if(ier /= 0) then |
---|
364 | call MP_perr(myname_,'MPI_Bcast(%buffer)',ier) |
---|
365 | istat=ier |
---|
366 | return |
---|
367 | endif |
---|
368 | |
---|
369 | call MPI_Bcast(i90_now%nbuf,1,MP_INTEGER,root,comm,ier) |
---|
370 | if(ier /= 0) then |
---|
371 | call MP_perr(myname_,'MPI_Bcast(%nbuf)',ier) |
---|
372 | istat=ier |
---|
373 | return |
---|
374 | endif |
---|
375 | |
---|
376 | i90_now%this_line=' ' |
---|
377 | i90_now%next_line=0 |
---|
378 | |
---|
379 | end subroutine I90_allLoadF |
---|
380 | |
---|
381 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
382 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
383 | !BOP ------------------------------------------------------------------- |
---|
384 | ! |
---|
385 | ! !IROUTINE: push_ - push on a new layer of the internal file _i90_now_ |
---|
386 | ! |
---|
387 | ! !DESCRIPTION: |
---|
388 | ! |
---|
389 | ! !INTERFACE: |
---|
390 | |
---|
391 | subroutine push_(ier) |
---|
392 | use m_die, only : perr |
---|
393 | use m_mall,only : mall_mci,mall_ci,mall_ison |
---|
394 | implicit none |
---|
395 | integer,intent(out) :: ier |
---|
396 | |
---|
397 | ! !REVISION HISTORY: |
---|
398 | ! 05Aug98 - Jing Guo <guo@thunder> - initial prototype/prolog/code |
---|
399 | !EOP ___________________________________________________________________ |
---|
400 | |
---|
401 | character(len=*),parameter :: myname_=myname//'::push_' |
---|
402 | type(inpak90),pointer :: new |
---|
403 | |
---|
404 | if(i90_depth <= 0) nullify(i90_now) ! just an initialization |
---|
405 | |
---|
406 | ! Too many levels |
---|
407 | |
---|
408 | if(i90_depth >= i90_MXDEP) then |
---|
409 | call perr(myname_,'(overflow)',i90_depth) |
---|
410 | ier=1 |
---|
411 | return |
---|
412 | endif |
---|
413 | |
---|
414 | allocate(new,stat=ier) |
---|
415 | if(ier /= 0) then |
---|
416 | call perr(myname_,'allocate(new)',ier) |
---|
417 | return |
---|
418 | endif |
---|
419 | |
---|
420 | if(mall_ison()) call mall_ci(MALLSIZE_,myname) |
---|
421 | |
---|
422 | allocate(new%buffer,new%this_line,stat=ier) |
---|
423 | if(ier /= 0) then |
---|
424 | call perr(myname_,'allocate(new%..)',ier) |
---|
425 | return |
---|
426 | endif |
---|
427 | |
---|
428 | if(mall_ison()) then |
---|
429 | call mall_mci(new%buffer,myname) |
---|
430 | call mall_mci(new%this_line,myname) |
---|
431 | endif |
---|
432 | |
---|
433 | new%last => i90_now |
---|
434 | i90_now => new |
---|
435 | nullify(new) |
---|
436 | |
---|
437 | i90_depth = i90_depth+1 |
---|
438 | end subroutine push_ |
---|
439 | |
---|
440 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
441 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
442 | !BOP ------------------------------------------------------------------- |
---|
443 | ! |
---|
444 | ! !IROUTINE: pop_ - pop off a layer of the internal file _i90_now_ |
---|
445 | ! |
---|
446 | ! !DESCRIPTION: |
---|
447 | ! |
---|
448 | ! !INTERFACE: |
---|
449 | |
---|
450 | subroutine pop_(ier) |
---|
451 | use m_die, only : perr |
---|
452 | use m_mall,only : mall_mco,mall_co,mall_ison |
---|
453 | implicit none |
---|
454 | integer,intent(out) :: ier |
---|
455 | |
---|
456 | ! !REVISION HISTORY: |
---|
457 | ! 05Aug98 - Jing Guo <guo@thunder> - initial prototype/prolog/code |
---|
458 | !EOP ___________________________________________________________________ |
---|
459 | |
---|
460 | character(len=*),parameter :: myname_=myname//'::pop_' |
---|
461 | type(inpak90),pointer :: old |
---|
462 | |
---|
463 | if(i90_depth <= 0) then |
---|
464 | call perr(myname_,'(underflow)',i90_depth) |
---|
465 | ier=1 |
---|
466 | return |
---|
467 | endif |
---|
468 | |
---|
469 | old => i90_now%last |
---|
470 | |
---|
471 | if(mall_ison()) then |
---|
472 | call mall_mco(i90_now%this_line,myname) |
---|
473 | call mall_mco(i90_now%buffer,myname) |
---|
474 | endif |
---|
475 | |
---|
476 | deallocate(i90_now%buffer,i90_now%this_line,stat=ier) |
---|
477 | if(ier /= 0) then |
---|
478 | call perr(myname_,'deallocate(new%..)',ier) |
---|
479 | return |
---|
480 | endif |
---|
481 | |
---|
482 | if(mall_ison()) call mall_co(MALLSIZE_,myname) |
---|
483 | |
---|
484 | deallocate(i90_now,stat=ier) |
---|
485 | if(ier /= 0) then |
---|
486 | call perr(myname_,'deallocate(new)',ier) |
---|
487 | return |
---|
488 | endif |
---|
489 | |
---|
490 | i90_now => old |
---|
491 | nullify(old) |
---|
492 | |
---|
493 | i90_depth = i90_depth - 1 |
---|
494 | end subroutine pop_ |
---|
495 | |
---|
496 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
497 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
498 | !----------------------------------------------------------------------- |
---|
499 | ! |
---|
500 | ! !ROUTINE: I90_Release - deallocate memory used to load a resource file |
---|
501 | ! |
---|
502 | ! !INTERFACE: |
---|
503 | ! |
---|
504 | subroutine I90_Release(stat) |
---|
505 | use m_die,only : perr,die |
---|
506 | implicit none |
---|
507 | integer,optional, intent(out) :: stat |
---|
508 | ! |
---|
509 | ! !DESCRIPTION: |
---|
510 | ! |
---|
511 | ! I90_Release() is used to pair I90_LoadF() to release the memory |
---|
512 | ! used by I90_LoadF() for resourse data input. |
---|
513 | ! |
---|
514 | ! !SEE ALSO: |
---|
515 | ! |
---|
516 | ! !REVISION HISTORY: |
---|
517 | ! 03Jul96 - J. Guo - added to Arlindos inpak90 for its |
---|
518 | ! Fortran 90 revision. |
---|
519 | !_______________________________________________________________________ |
---|
520 | character(len=*),parameter :: myname_=myname//'::i90_Release' |
---|
521 | integer :: ier |
---|
522 | |
---|
523 | if(present(stat)) stat=0 |
---|
524 | |
---|
525 | call pop_(ier) |
---|
526 | if(ier/=0) then |
---|
527 | call perr(myname_,'pop_()',ier) |
---|
528 | if(.not.present(stat)) call die(myname_) |
---|
529 | stat=ier |
---|
530 | return |
---|
531 | endif |
---|
532 | |
---|
533 | end subroutine I90_Release |
---|
534 | |
---|
535 | !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
536 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
537 | !BOP ------------------------------------------------------------------- |
---|
538 | ! |
---|
539 | ! !IROUTINE: i90_fullRelease - releases the whole stack led by _i90_now_ |
---|
540 | ! |
---|
541 | ! !DESCRIPTION: |
---|
542 | ! |
---|
543 | ! !INTERFACE: |
---|
544 | |
---|
545 | subroutine i90_fullRelease(ier) |
---|
546 | use m_die,only : perr |
---|
547 | implicit none |
---|
548 | integer,intent(out) :: ier |
---|
549 | |
---|
550 | ! !REVISION HISTORY: |
---|
551 | ! 05Aug98 - Jing Guo <guo@thunder> - initial prototype/prolog/code |
---|
552 | !EOP ___________________________________________________________________ |
---|
553 | |
---|
554 | character(len=*),parameter :: myname_=myname//'::i90_fullRelease' |
---|
555 | |
---|
556 | do while(i90_depth > 0) |
---|
557 | call pop_(ier) |
---|
558 | if(ier /= 0) then |
---|
559 | call perr(myname_,'pop_()',ier) |
---|
560 | return |
---|
561 | endif |
---|
562 | end do |
---|
563 | ier=0 |
---|
564 | |
---|
565 | end subroutine i90_fullRelease |
---|
566 | !======================================================================= |
---|
567 | subroutine I90_LoadF ( filen, iret ) |
---|
568 | use m_ioutil, only : luavail,opntext,clstext |
---|
569 | use m_die, only : perr |
---|
570 | implicit NONE |
---|
571 | |
---|
572 | !------------------------------------------------------------------------- |
---|
573 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
574 | !------------------------------------------------------------------------- |
---|
575 | !BOP |
---|
576 | ! |
---|
577 | ! !ROUTINE: I90_LoadF() --- Loads resource file into memory. |
---|
578 | ! |
---|
579 | ! !DESCRIPTION: |
---|
580 | ! |
---|
581 | ! Reads resource file, strips out comments, translate TABs into |
---|
582 | ! blanks, and loads the modified file contents into memory. |
---|
583 | ! Must be called only once for each resource file. |
---|
584 | ! |
---|
585 | ! !CALLING SEQUENCE: |
---|
586 | ! |
---|
587 | ! call i90_LoadF ( filen, iret ) |
---|
588 | ! |
---|
589 | ! !INPUT PARAMETERS: |
---|
590 | ! |
---|
591 | character*(*) filen ! file name |
---|
592 | |
---|
593 | ! !OUTPUT PARAMETERS: |
---|
594 | |
---|
595 | integer iret ! Return code: |
---|
596 | ! 0 no error |
---|
597 | ! -98 coult not get unit number |
---|
598 | ! (strange!) |
---|
599 | ! -98 talk to a wizzard |
---|
600 | ! -99 out of memory: increase |
---|
601 | ! NBUF_MAX in 'i90.h' |
---|
602 | ! other iostat from open statement. |
---|
603 | ! |
---|
604 | ! !BUGS: |
---|
605 | ! |
---|
606 | ! It does not perform dynamic allocation, mostly to keep vanilla f77 |
---|
607 | ! compatibility. Overall amount of static memory is small (~100K |
---|
608 | ! for default NBUF_MAX = 400*256). |
---|
609 | ! |
---|
610 | ! !SEE ALSO: |
---|
611 | ! |
---|
612 | ! i90_label() selects a label (key) |
---|
613 | ! |
---|
614 | ! !FILES USED: |
---|
615 | ! |
---|
616 | ! File name supplied on input. The file is opened, read and then closed. |
---|
617 | ! |
---|
618 | ! !REVISION HISTORY: |
---|
619 | ! |
---|
620 | ! 19Jun96 da Silva Original code. |
---|
621 | ! |
---|
622 | !EOP |
---|
623 | !------------------------------------------------------------------------- |
---|
624 | integer lu, ios, loop, ls, ptr |
---|
625 | character*256 line |
---|
626 | character(len=*), parameter :: myname_ = myname//'::i90_loadf' |
---|
627 | |
---|
628 | ! Check to make sure there is not too many levels |
---|
629 | ! of the stacked resource files |
---|
630 | |
---|
631 | if(i90_depth >= i90_MXDEP) then |
---|
632 | call perr(myname_,'(overflow)',i90_depth) |
---|
633 | iret=1 |
---|
634 | return |
---|
635 | endif |
---|
636 | |
---|
637 | ! Open file |
---|
638 | ! --------- |
---|
639 | ! lu = i90_lua() |
---|
640 | |
---|
641 | lu = luavail() ! a more portable version |
---|
642 | if ( lu .lt. 0 ) then |
---|
643 | iret = -97 |
---|
644 | return |
---|
645 | end if |
---|
646 | |
---|
647 | ! A open through an interface to avoid portability problems. |
---|
648 | ! (J.G.) |
---|
649 | |
---|
650 | call opntext(lu,filen,'old',ios) |
---|
651 | if ( ios .ne. 0 ) then |
---|
652 | write(stderr,'(2a,i5)') myname_,': opntext() error, ios =',ios |
---|
653 | iret = ios |
---|
654 | return |
---|
655 | end if |
---|
656 | |
---|
657 | ! Create a dynamic page to store the file. It might be expanded |
---|
658 | ! to allocate memory on requests (a link list) (J.G.) |
---|
659 | |
---|
660 | ! Changed from page_() to push_(), to allow multiple (stacked) |
---|
661 | ! inpak90 buffers. J.G. |
---|
662 | |
---|
663 | call push_(ios) ! to create buffer space |
---|
664 | if ( ios .ne. 0 ) then |
---|
665 | write(stderr,'(2a,i5)') myname_,': push_() error, ios =',ios |
---|
666 | iret = ios |
---|
667 | return |
---|
668 | end if |
---|
669 | |
---|
670 | ! Read to end of file |
---|
671 | ! ------------------- |
---|
672 | i90_now%buffer(1:1) = EOL |
---|
673 | ptr = 2 ! next buffer position |
---|
674 | do loop = 1, NBUF_MAX |
---|
675 | |
---|
676 | ! Read next line |
---|
677 | ! -------------- |
---|
678 | read(lu,'(a)', end=11) line ! read next line |
---|
679 | call i90_trim ( line ) ! remove trailing blanks |
---|
680 | call i90_pad ( line ) ! Pad with # from end of line |
---|
681 | |
---|
682 | ! A non-empty line |
---|
683 | ! ---------------- |
---|
684 | ls = index(line,'#' ) - 1 ! line length |
---|
685 | if ( ls .gt. 0 ) then |
---|
686 | if ( (ptr+ls) .gt. NBUF_MAX ) then |
---|
687 | iret = -99 |
---|
688 | return |
---|
689 | end if |
---|
690 | i90_now%buffer(ptr:ptr+ls) = line(1:ls) // EOL |
---|
691 | ptr = ptr + ls + 1 |
---|
692 | end if |
---|
693 | |
---|
694 | end do |
---|
695 | |
---|
696 | iret = -98 ! good chance i90_now%buffer is not big enough |
---|
697 | return |
---|
698 | |
---|
699 | 11 continue |
---|
700 | |
---|
701 | ! All done |
---|
702 | ! -------- |
---|
703 | ! close(lu) |
---|
704 | call clstext(lu,ios) |
---|
705 | if(ios /= 0) then |
---|
706 | iret=-99 |
---|
707 | return |
---|
708 | endif |
---|
709 | i90_now%buffer(ptr:ptr) = EOB |
---|
710 | i90_now%nbuf = ptr |
---|
711 | i90_now%this_line=' ' |
---|
712 | i90_now%next_line=0 |
---|
713 | iret = 0 |
---|
714 | |
---|
715 | return |
---|
716 | end subroutine I90_LoadF |
---|
717 | |
---|
718 | |
---|
719 | !................................................................... |
---|
720 | |
---|
721 | subroutine i90_label ( label, iret ) |
---|
722 | |
---|
723 | implicit NONE |
---|
724 | |
---|
725 | !------------------------------------------------------------------------- |
---|
726 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
727 | !------------------------------------------------------------------------- |
---|
728 | !BOP |
---|
729 | ! |
---|
730 | ! !ROUTINE: I90_Label() --- Selects a label (record). |
---|
731 | ! |
---|
732 | ! !DESCRIPTION: |
---|
733 | ! |
---|
734 | ! Once the buffer has been loaded with {\tt i90\_loadf()}, this routine |
---|
735 | ! selects a given ``line'' (record/table) associated with ``label''. |
---|
736 | ! Think of ``label'' as a resource name or data base ``key''. |
---|
737 | ! |
---|
738 | ! !CALLING SEQUENCE: |
---|
739 | ! |
---|
740 | ! call i90_Label ( label, iret ) |
---|
741 | ! |
---|
742 | ! !INPUT PARAMETERS: |
---|
743 | ! |
---|
744 | character(len=*),intent(in) :: label ! input label |
---|
745 | |
---|
746 | ! !OUTPUT PARAMETERS: |
---|
747 | |
---|
748 | integer iret ! Return code: |
---|
749 | ! 0 no error |
---|
750 | ! -1 buffer not loaded |
---|
751 | ! -2 could not find label |
---|
752 | ! |
---|
753 | ! !SEE ALSO: |
---|
754 | ! |
---|
755 | ! i90_loadf() load file into buffer |
---|
756 | ! i90_gtoken() get next token |
---|
757 | ! i90_gline() get next line (for tables) |
---|
758 | ! atof() convert word (string) to float |
---|
759 | ! atoi() convert word (string) to integer |
---|
760 | ! |
---|
761 | ! !REVISION HISTORY: |
---|
762 | ! |
---|
763 | ! 19Jun96 da Silva Original code. |
---|
764 | ! 19Jan01 Jay Larson <larson@mcs.anl.gov> - introduced CHARACTER |
---|
765 | ! variable EOL_label, which is used to circumvent pgf90 |
---|
766 | ! problems with passing concatenated characters as an argument |
---|
767 | ! to a function. |
---|
768 | ! |
---|
769 | !EOP |
---|
770 | !------------------------------------------------------------------------- |
---|
771 | |
---|
772 | integer i, j |
---|
773 | |
---|
774 | character(len=(len(label)+len(EOL))) :: EOL_label |
---|
775 | |
---|
776 | ! Make sure that a buffer is defined (JG) |
---|
777 | ! ---------------------------------- |
---|
778 | if(i90_depth <= 0) then |
---|
779 | iret = -1 |
---|
780 | return |
---|
781 | endif |
---|
782 | |
---|
783 | ! Determine whether label exists |
---|
784 | ! ------------------------------ |
---|
785 | EOL_label = EOL // label |
---|
786 | i = index ( i90_now%buffer(1:i90_now%nbuf), EOL_label ) + 1 |
---|
787 | if ( i .le. 1 ) then |
---|
788 | i90_now%this_line = BLK // EOL |
---|
789 | iret = -2 |
---|
790 | return |
---|
791 | end if |
---|
792 | |
---|
793 | ! Extract the line associated with this label |
---|
794 | ! ------------------------------------------- |
---|
795 | i = i + len ( label ) |
---|
796 | j = i + index(i90_now%buffer(i:i90_now%nbuf),EOL) - 2 |
---|
797 | i90_now%this_line = i90_now%buffer(i:j) // BLK // EOL |
---|
798 | |
---|
799 | i90_now%next_line = j + 2 |
---|
800 | |
---|
801 | iret = 0 |
---|
802 | |
---|
803 | return |
---|
804 | end subroutine i90_label |
---|
805 | |
---|
806 | !................................................................... |
---|
807 | |
---|
808 | subroutine i90_gline ( iret ) |
---|
809 | |
---|
810 | implicit NONE |
---|
811 | |
---|
812 | !------------------------------------------------------------------------- |
---|
813 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
814 | !------------------------------------------------------------------------- |
---|
815 | !BOP |
---|
816 | ! |
---|
817 | ! !ROUTINE: I90_GLine() --- Selects next line. |
---|
818 | ! |
---|
819 | ! !DESCRIPTION: |
---|
820 | ! |
---|
821 | ! Selects next line, irrespective of of label. If the next line starts |
---|
822 | ! with :: (end of table mark), then it lets the user know. This sequential |
---|
823 | ! access of the buffer is useful to assess tables, a concept introduced |
---|
824 | ! in Inpak 77 by Jing Guo. A table is a construct like this: |
---|
825 | ! |
---|
826 | ! \begin{verbatim} |
---|
827 | ! my_table_name:: |
---|
828 | ! 1000 3000 263.0 |
---|
829 | ! 925 3000 263.0 |
---|
830 | ! 850 3000 263.0 |
---|
831 | ! 700 3000 269.0 |
---|
832 | ! 500 3000 287.0 |
---|
833 | ! 400 3000 295.8 |
---|
834 | ! 300 3000 295.8 |
---|
835 | ! :: |
---|
836 | ! \end{verbatim} |
---|
837 | ! |
---|
838 | ! To access this table, the user first must use {\tt i90\_label()} to |
---|
839 | ! locate the beginning of the table, e.g., |
---|
840 | ! |
---|
841 | ! \begin{verbatim} |
---|
842 | ! call i90_label ( 'my_table_name::', iret ) |
---|
843 | ! \end{verbatim} |
---|
844 | ! |
---|
845 | ! Subsequently, {\tt i90\_gline()} can be used to gain acess to each |
---|
846 | ! row of the table. Here is a code fragment to read the above |
---|
847 | ! table (7 rows, 3 columns): |
---|
848 | ! |
---|
849 | ! \begin{verbatim} |
---|
850 | ! real table(7,3) |
---|
851 | ! character*20 word |
---|
852 | ! integer iret |
---|
853 | ! call i90_label ( 'my_table_name::', iret ) |
---|
854 | ! do i = 1, 7 |
---|
855 | ! call i90_gline ( iret ) |
---|
856 | ! do j = 1, 3 |
---|
857 | ! table(i,j) = fltget ( 0. ) |
---|
858 | ! end do |
---|
859 | ! end do |
---|
860 | ! \end{verbatim} |
---|
861 | ! |
---|
862 | ! For simplicity we have assumed that the dimensions of table were |
---|
863 | ! known. It is relatively simple to infer the table dimensions |
---|
864 | ! by manipulating ``iret''. |
---|
865 | ! |
---|
866 | ! !CALLING SEQUENCE: |
---|
867 | ! |
---|
868 | ! call i90_gline ( iret ) |
---|
869 | ! |
---|
870 | ! !INPUT PARAMETERS: |
---|
871 | ! |
---|
872 | ! None. |
---|
873 | ! |
---|
874 | ! !OUTPUT PARAMETERS: |
---|
875 | ! |
---|
876 | integer iret ! Return code: |
---|
877 | ! 0 no error |
---|
878 | ! -1 end of buffer reached |
---|
879 | ! +1 end of table reached |
---|
880 | |
---|
881 | ! !SEE ALSO: |
---|
882 | ! |
---|
883 | ! i90_label() selects a line (record/table) |
---|
884 | ! |
---|
885 | ! !REVISION HISTORY: |
---|
886 | ! |
---|
887 | ! 10feb95 Guo Wrote rdnext(), Inpak 77 extension. |
---|
888 | ! 19Jun96 da Silva Original code with functionality of rdnext() |
---|
889 | ! |
---|
890 | !EOP |
---|
891 | !------------------------------------------------------------------------- |
---|
892 | |
---|
893 | integer i, j |
---|
894 | |
---|
895 | ! Make sure that a buffer is defined (JG) |
---|
896 | ! ---------------------------------- |
---|
897 | if(i90_depth <= 0) then |
---|
898 | iret = -1 |
---|
899 | return |
---|
900 | endif |
---|
901 | |
---|
902 | if ( i90_now%next_line .ge. i90_now%nbuf ) then |
---|
903 | iret = -1 |
---|
904 | return |
---|
905 | end if |
---|
906 | |
---|
907 | i = i90_now%next_line |
---|
908 | j = i + index(i90_now%buffer(i:i90_now%nbuf),EOL) - 2 |
---|
909 | i90_now%this_line = i90_now%buffer(i:j) // BLK // EOL |
---|
910 | |
---|
911 | if ( i90_now%this_line(1:2) .eq. '::' ) then |
---|
912 | iret = 1 ! end of table |
---|
913 | i90_now%next_line = i90_now%nbuf + 1 |
---|
914 | return |
---|
915 | end if |
---|
916 | |
---|
917 | i90_now%next_line = j + 2 |
---|
918 | iret = 0 |
---|
919 | |
---|
920 | return |
---|
921 | end subroutine i90_gline |
---|
922 | |
---|
923 | !................................................................... |
---|
924 | |
---|
925 | subroutine i90_GToken ( token, iret ) |
---|
926 | |
---|
927 | implicit NONE |
---|
928 | |
---|
929 | !------------------------------------------------------------------------- |
---|
930 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
931 | !------------------------------------------------------------------------- |
---|
932 | !BOP |
---|
933 | ! |
---|
934 | ! !ROUTINE: I90_GToken() --- Gets next token. |
---|
935 | ! |
---|
936 | ! !DESCRIPTION: |
---|
937 | ! |
---|
938 | ! Get next token from current line. The current line is defined by a |
---|
939 | ! call to {\tt i90\_label()}. Tokens are sequences of characters (including |
---|
940 | ! blanks) which may be enclosed by single or double quotes. |
---|
941 | ! If no quotes are present, the token from the current position to the next |
---|
942 | ! blank of TAB is returned. |
---|
943 | ! |
---|
944 | ! {\em Examples of valid token:} |
---|
945 | ! |
---|
946 | ! \begin{verbatim} |
---|
947 | ! single_token "second token on line" |
---|
948 | ! "this is a token" |
---|
949 | ! 'Another example of a token' |
---|
950 | ! 'this is how you get a " inside a token' |
---|
951 | ! "this is how you get a ' inside a token" |
---|
952 | ! This is valid too # the line ends before the # |
---|
953 | ! \end{verbatim} |
---|
954 | ! The last line has 4 valid tokens: {\tt This, is, valid} and {\tt too}. |
---|
955 | ! |
---|
956 | ! {\em Invalid string constructs:} |
---|
957 | ! |
---|
958 | ! \begin{verbatim} |
---|
959 | ! cannot handle mixed quotes (i.e. single/double) |
---|
960 | ! 'escaping like this \' is not implemented' |
---|
961 | ! 'this # will not work because of the #' |
---|
962 | ! \end{verbatim} |
---|
963 | ! The \# character is reserved for comments and cannot be included |
---|
964 | ! inside quotation marks. |
---|
965 | ! |
---|
966 | ! !CALLING SEQUENCE: |
---|
967 | ! |
---|
968 | ! call i90_GToken ( token, iret ) |
---|
969 | ! |
---|
970 | ! !INPUT PARAMETERS: |
---|
971 | ! |
---|
972 | ! None. |
---|
973 | ! |
---|
974 | ! !OUTPUT PARAMETERS: |
---|
975 | ! |
---|
976 | character*(*) token ! Next token from current line |
---|
977 | integer iret ! Return code: |
---|
978 | ! 0 no error |
---|
979 | ! -1 either nothing left |
---|
980 | ! on line or mismatched |
---|
981 | ! quotation marks. |
---|
982 | |
---|
983 | ! !BUGS: |
---|
984 | ! |
---|
985 | ! Standard Unix escaping is not implemented at the moment. |
---|
986 | ! |
---|
987 | ! |
---|
988 | ! !SEE ALSO: |
---|
989 | ! |
---|
990 | ! i90_label() selects a line (record/table) |
---|
991 | ! i90_gline() get next line (for tables) |
---|
992 | ! atof() convert word (string) to float |
---|
993 | ! atoi() convert word (string) to integer |
---|
994 | ! |
---|
995 | ! |
---|
996 | ! !REVISION HISTORY: |
---|
997 | ! |
---|
998 | ! 19Jun96 da Silva Original code. |
---|
999 | ! |
---|
1000 | !EOP |
---|
1001 | !------------------------------------------------------------------------- |
---|
1002 | |
---|
1003 | character*1 ch |
---|
1004 | integer ib, ie |
---|
1005 | |
---|
1006 | ! Make sure that a buffer is defined (JG) |
---|
1007 | ! ---------------------------------- |
---|
1008 | if(i90_depth <= 0) then |
---|
1009 | iret = -1 |
---|
1010 | return |
---|
1011 | endif |
---|
1012 | |
---|
1013 | call i90_trim ( i90_now%this_line ) |
---|
1014 | |
---|
1015 | ch = i90_now%this_line(1:1) |
---|
1016 | if ( ch .eq. '"' .or. ch .eq. "'" ) then |
---|
1017 | ib = 2 |
---|
1018 | ie = index ( i90_now%this_line(ib:), ch ) |
---|
1019 | else |
---|
1020 | ib = 1 |
---|
1021 | ie = min(index(i90_now%this_line,BLK), & |
---|
1022 | index(i90_now%this_line,EOL)) - 1 |
---|
1023 | |
---|
1024 | end if |
---|
1025 | |
---|
1026 | if ( ie .lt. ib ) then |
---|
1027 | token = BLK |
---|
1028 | iret = -1 |
---|
1029 | return |
---|
1030 | else |
---|
1031 | ! Get the token, and shift the rest of %this_line to |
---|
1032 | ! the left |
---|
1033 | |
---|
1034 | token = i90_now%this_line(ib:ie) |
---|
1035 | i90_now%this_line = i90_now%this_line(ie+2:) |
---|
1036 | iret = 0 |
---|
1037 | end if |
---|
1038 | |
---|
1039 | return |
---|
1040 | end subroutine i90_gtoken |
---|
1041 | !................................................................... |
---|
1042 | subroutine i90_gstr ( string, iret ) |
---|
1043 | |
---|
1044 | implicit NONE |
---|
1045 | |
---|
1046 | !------------------------------------------------------------------------- |
---|
1047 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
1048 | !------------------------------------------------------------------------- |
---|
1049 | ! |
---|
1050 | ! !ROUTINE: I90\_GStr() |
---|
1051 | ! |
---|
1052 | ! !DESCRIPTION: |
---|
1053 | ! |
---|
1054 | ! Get next string from current line. The current line is defined by a |
---|
1055 | ! call to {\tt i90\_label()}. Strings are sequence of characters (including |
---|
1056 | ! blanks) enclosed by single or double quotes. If no quotes |
---|
1057 | ! are present, the string from the current position to the end of |
---|
1058 | ! the line is returned. |
---|
1059 | ! |
---|
1060 | ! NOTE: This routine is defined differently from \verb"i90_GTolen()", |
---|
1061 | ! where a {\sl token} is white-space delimited, but this routine |
---|
1062 | ! will try to fetch a string either terminated by a "$" or by the |
---|
1063 | ! end of the line. |
---|
1064 | ! |
---|
1065 | ! {\em Examples of valid strings:} |
---|
1066 | ! |
---|
1067 | ! \begin{verbatim} |
---|
1068 | ! "this is a string" |
---|
1069 | ! 'Another example of string' |
---|
1070 | ! 'this is how you get a " inside a string' |
---|
1071 | ! "this is how you get a ' inside a string" |
---|
1072 | ! This is valid too # the line ends before the # |
---|
1073 | ! |
---|
1074 | ! \end{verbatim} |
---|
1075 | ! |
---|
1076 | ! {\em Invalid string constructs:} |
---|
1077 | ! |
---|
1078 | ! \begin{verbatim} |
---|
1079 | ! cannot handle mixed quotes |
---|
1080 | ! 'escaping like this \' is not implemented' |
---|
1081 | ! \end{verbatim} |
---|
1082 | ! |
---|
1083 | ! {\em Obsolete feature (for Inpak 77 compatibility):} |
---|
1084 | ! |
---|
1085 | ! \begin{verbatim} |
---|
1086 | ! the string ends after a $ this is another string |
---|
1087 | ! \end{verbatim} |
---|
1088 | ! |
---|
1089 | ! !CALLING SEQUENCE: |
---|
1090 | ! |
---|
1091 | ! \begin{verbatim} |
---|
1092 | ! call i90_Gstr ( string, iret ) |
---|
1093 | ! \end{verbatim} |
---|
1094 | ! |
---|
1095 | ! !INPUT PARAMETERS: |
---|
1096 | ! |
---|
1097 | character*(*) string ! A NULL (char(0)) delimited string. |
---|
1098 | |
---|
1099 | ! !OUTPUT PARAMETERS: |
---|
1100 | ! |
---|
1101 | integer iret ! Return code: |
---|
1102 | ! 0 no error |
---|
1103 | ! -1 either nothing left |
---|
1104 | ! on line or mismatched |
---|
1105 | ! quotation marks. |
---|
1106 | |
---|
1107 | ! !BUGS: |
---|
1108 | ! |
---|
1109 | ! Standard Unix escaping is not implemented at the moment. |
---|
1110 | ! No way to tell sintax error from end of line (same iret). |
---|
1111 | ! |
---|
1112 | ! |
---|
1113 | ! !SEE ALSO: |
---|
1114 | ! |
---|
1115 | ! i90_label() selects a line (record/table) |
---|
1116 | ! i90_gtoken() get next token |
---|
1117 | ! i90_gline() get next line (for tables) |
---|
1118 | ! atof() convert word (string) to float |
---|
1119 | ! atoi() convert word (string) to integer |
---|
1120 | ! |
---|
1121 | ! |
---|
1122 | ! !REVISION HISTORY: |
---|
1123 | ! |
---|
1124 | ! 19Jun96 da Silva Original code. |
---|
1125 | ! 01Oct96 Jing Guo Removed the null terminitor |
---|
1126 | ! |
---|
1127 | !------------------------------------------------------------------------- |
---|
1128 | |
---|
1129 | character*1 ch |
---|
1130 | integer ib, ie |
---|
1131 | |
---|
1132 | ! Make sure that a buffer is defined (JG) |
---|
1133 | ! ---------------------------------- |
---|
1134 | if(i90_depth <= 0) then |
---|
1135 | iret = -1 |
---|
1136 | return |
---|
1137 | endif |
---|
1138 | |
---|
1139 | call i90_trim ( i90_now%this_line ) |
---|
1140 | |
---|
1141 | ch = i90_now%this_line(1:1) |
---|
1142 | if ( ch .eq. '"' .or. ch .eq. "'" ) then |
---|
1143 | ib = 2 |
---|
1144 | ie = index ( i90_now%this_line(ib:), ch ) |
---|
1145 | else |
---|
1146 | ib = 1 |
---|
1147 | ie = index(i90_now%this_line,'$')-1 ! undocumented feature! |
---|
1148 | if ( ie .lt. 1 ) ie = index(i90_now%this_line,EOL)-2 |
---|
1149 | end if |
---|
1150 | |
---|
1151 | if ( ie .lt. ib ) then |
---|
1152 | ! string = NULL |
---|
1153 | iret = -1 |
---|
1154 | return |
---|
1155 | else |
---|
1156 | string = i90_now%this_line(ib:ie) ! // NULL |
---|
1157 | i90_now%this_line = i90_now%this_line(ie+2:) |
---|
1158 | iret = 0 |
---|
1159 | end if |
---|
1160 | |
---|
1161 | return |
---|
1162 | end subroutine i90_gstr |
---|
1163 | |
---|
1164 | !................................................................... |
---|
1165 | |
---|
1166 | real(FP) function i90_GFloat( iret ) |
---|
1167 | |
---|
1168 | implicit NONE |
---|
1169 | |
---|
1170 | !------------------------------------------------------------------------- |
---|
1171 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
1172 | !------------------------------------------------------------------------- |
---|
1173 | !BOP |
---|
1174 | ! |
---|
1175 | ! !ROUTINE: i90_GFloat() --- Returns next float number. |
---|
1176 | ! |
---|
1177 | ! !DESCRIPTION: |
---|
1178 | ! |
---|
1179 | ! Returns next float (real number) from the current line. |
---|
1180 | ! If an error occurs a zero value is returned. |
---|
1181 | ! |
---|
1182 | ! !CALLING SEQUENCE: |
---|
1183 | ! |
---|
1184 | ! real rnumber |
---|
1185 | ! rnumber = i90_gfloat ( default ) |
---|
1186 | ! |
---|
1187 | ! !OUTPUT PARAMETERS: |
---|
1188 | ! |
---|
1189 | integer,intent(out) :: iret ! Return code: |
---|
1190 | ! 0 no error |
---|
1191 | ! -1 either nothing left |
---|
1192 | ! on line or mismatched |
---|
1193 | ! quotation marks. |
---|
1194 | ! -2 parsing error |
---|
1195 | |
---|
1196 | ! |
---|
1197 | ! !REVISION HISTORY: |
---|
1198 | ! |
---|
1199 | ! 19Jun96 da Silva Original code. |
---|
1200 | ! |
---|
1201 | !EOP |
---|
1202 | !------------------------------------------------------------------------- |
---|
1203 | |
---|
1204 | character*256 token |
---|
1205 | integer ios |
---|
1206 | real(FP) x |
---|
1207 | |
---|
1208 | ! Make sure that a buffer is defined (JG) |
---|
1209 | ! ---------------------------------- |
---|
1210 | if(i90_depth <= 0) then |
---|
1211 | iret = -1 |
---|
1212 | return |
---|
1213 | endif |
---|
1214 | |
---|
1215 | call i90_gtoken ( token, iret ) |
---|
1216 | if ( iret .eq. 0 ) then |
---|
1217 | read(token,*,iostat=ios) x ! Does it require an extension? |
---|
1218 | if ( ios .ne. 0 ) iret = -2 |
---|
1219 | end if |
---|
1220 | if ( iret .ne. 0 ) x = 0. |
---|
1221 | i90_GFloat = x |
---|
1222 | |
---|
1223 | return |
---|
1224 | end function i90_GFloat |
---|
1225 | |
---|
1226 | !................................................................... |
---|
1227 | |
---|
1228 | integer function I90_GInt ( iret ) |
---|
1229 | |
---|
1230 | implicit NONE |
---|
1231 | |
---|
1232 | !------------------------------------------------------------------------- |
---|
1233 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
1234 | !------------------------------------------------------------------------- |
---|
1235 | !BOP |
---|
1236 | ! |
---|
1237 | ! !ROUTINE: I90_GInt() --- Returns next integer number. |
---|
1238 | ! |
---|
1239 | ! !DESCRIPTION: |
---|
1240 | ! |
---|
1241 | ! Returns next integer number from the current line. |
---|
1242 | ! If an error occurs a zero value is returned. |
---|
1243 | ! |
---|
1244 | ! !CALLING SEQUENCE: |
---|
1245 | ! |
---|
1246 | ! integer number |
---|
1247 | ! number = i90_gint ( default ) |
---|
1248 | ! |
---|
1249 | ! !OUTPUT PARAMETERS: |
---|
1250 | ! |
---|
1251 | integer iret ! Return code: |
---|
1252 | ! 0 no error |
---|
1253 | ! -1 either nothing left |
---|
1254 | ! on line or mismatched |
---|
1255 | ! quotation marks. |
---|
1256 | ! -2 parsing error |
---|
1257 | |
---|
1258 | ! |
---|
1259 | ! !REVISION HISTORY: |
---|
1260 | ! |
---|
1261 | ! 19Jun96 da Silva Original code. |
---|
1262 | ! 24may00 da Silva delcared x as real*8 in case this module is compiled |
---|
1263 | ! with real*4 |
---|
1264 | ! |
---|
1265 | !EOP |
---|
1266 | !------------------------------------------------------------------------- |
---|
1267 | |
---|
1268 | character*256 token |
---|
1269 | real(kind_r8) x |
---|
1270 | integer ios |
---|
1271 | |
---|
1272 | ! Make sure that a buffer is defined (JG) |
---|
1273 | ! ---------------------------------- |
---|
1274 | if(i90_depth <= 0) then |
---|
1275 | iret = -1 |
---|
1276 | return |
---|
1277 | endif |
---|
1278 | |
---|
1279 | call i90_gtoken ( token, iret ) |
---|
1280 | if ( iret .eq. 0 ) then |
---|
1281 | read(token,*,iostat=ios) x |
---|
1282 | if ( ios .ne. 0 ) iret = -2 |
---|
1283 | end if |
---|
1284 | if ( iret .ne. 0 ) x = 0 |
---|
1285 | i90_gint = nint(x) |
---|
1286 | |
---|
1287 | return |
---|
1288 | end function i90_gint |
---|
1289 | |
---|
1290 | !................................................................... |
---|
1291 | |
---|
1292 | real(FP) function i90_AtoF( string, iret ) |
---|
1293 | |
---|
1294 | implicit NONE |
---|
1295 | |
---|
1296 | !------------------------------------------------------------------------- |
---|
1297 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
1298 | !------------------------------------------------------------------------- |
---|
1299 | !BOP |
---|
1300 | ! |
---|
1301 | ! !ROUTINE: i90_AtoF() --- Translates ASCII (string) to float. |
---|
1302 | ! |
---|
1303 | ! !DESCRIPTION: |
---|
1304 | ! |
---|
1305 | ! Converts string to real number. Same as obsolete {\tt str2rn()}. |
---|
1306 | ! |
---|
1307 | ! !CALLING SEQUENCE: |
---|
1308 | ! |
---|
1309 | ! real rnumber |
---|
1310 | ! rnumber = i90_atof ( string, iret ) |
---|
1311 | ! |
---|
1312 | ! !INPUT PARAMETERS: |
---|
1313 | ! |
---|
1314 | character(len=*),intent(in) :: string ! a string |
---|
1315 | |
---|
1316 | ! !OUTPUT PARAMETERS: |
---|
1317 | ! |
---|
1318 | integer,intent(out) :: iret ! Return code: |
---|
1319 | ! 0 no error |
---|
1320 | ! -1 could not convert, probably |
---|
1321 | ! string is not a number |
---|
1322 | |
---|
1323 | ! |
---|
1324 | ! !REVISION HISTORY: |
---|
1325 | ! |
---|
1326 | ! 19Jun96 da Silva Original code. |
---|
1327 | ! |
---|
1328 | !EOP |
---|
1329 | !------------------------------------------------------------------------- |
---|
1330 | |
---|
1331 | read(string,*,end=11,err=11) i90_AtoF |
---|
1332 | iret = 0 |
---|
1333 | return |
---|
1334 | 11 iret = -1 |
---|
1335 | return |
---|
1336 | end function i90_AtoF |
---|
1337 | |
---|
1338 | !................................................................... |
---|
1339 | |
---|
1340 | integer function i90_atoi ( string, iret ) |
---|
1341 | |
---|
1342 | implicit NONE |
---|
1343 | |
---|
1344 | !------------------------------------------------------------------------- |
---|
1345 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
1346 | !------------------------------------------------------------------------- |
---|
1347 | !BOP |
---|
1348 | ! |
---|
1349 | ! !ROUTINE: I90_AtoI() --- Translates ASCII (strings) to integer. |
---|
1350 | ! |
---|
1351 | ! !DESCRIPTION: |
---|
1352 | ! |
---|
1353 | ! Converts string to integer number. |
---|
1354 | ! |
---|
1355 | ! !CALLING SEQUENCE: |
---|
1356 | ! |
---|
1357 | ! integer number |
---|
1358 | ! number = i90_atoi ( string, iret ) |
---|
1359 | ! |
---|
1360 | ! !INPUT PARAMETERS: |
---|
1361 | ! |
---|
1362 | character*(*) string ! a string |
---|
1363 | |
---|
1364 | ! !OUTPUT PARAMETERS: |
---|
1365 | ! |
---|
1366 | integer iret ! Return code: |
---|
1367 | ! 0 no error |
---|
1368 | ! -1 could not convert, probably |
---|
1369 | ! string is not a number |
---|
1370 | |
---|
1371 | ! |
---|
1372 | ! !REVISION HISTORY: |
---|
1373 | ! |
---|
1374 | ! 19Jun96 da Silva Original code. |
---|
1375 | ! |
---|
1376 | !EOP |
---|
1377 | !------------------------------------------------------------------------- |
---|
1378 | |
---|
1379 | read(string,*,end=11,err=11) i90_atoi |
---|
1380 | iret = 0 |
---|
1381 | return |
---|
1382 | 11 iret = -1 |
---|
1383 | return |
---|
1384 | end function i90_atoi |
---|
1385 | |
---|
1386 | !................................................................... |
---|
1387 | |
---|
1388 | integer function i90_Len ( string ) |
---|
1389 | |
---|
1390 | implicit NONE |
---|
1391 | |
---|
1392 | !------------------------------------------------------------------------- |
---|
1393 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
1394 | !------------------------------------------------------------------------- |
---|
1395 | !BOP |
---|
1396 | ! |
---|
1397 | ! !ROUTINE: I90_Len() --- Returns length of string. |
---|
1398 | ! |
---|
1399 | ! !DESCRIPTION: |
---|
1400 | ! |
---|
1401 | ! Returns the length of a string excluding trailing blanks. |
---|
1402 | ! It follows that |
---|
1403 | ! \begin{verbatim} |
---|
1404 | ! i90_len(string) .le. len(string), |
---|
1405 | ! \end{verbatim} |
---|
1406 | ! where {\tt len} is the intrinsic string length function. |
---|
1407 | ! Example: |
---|
1408 | ! \begin{verbatim} |
---|
1409 | ! ls = len('abc ') ! results in ls = 5 |
---|
1410 | ! ls = i90_len ('abc ') ! results in ls = 3 |
---|
1411 | ! \end{verbatim} |
---|
1412 | ! |
---|
1413 | ! !CALLING SEQUENCE: |
---|
1414 | ! |
---|
1415 | ! integer ls |
---|
1416 | ! ls = i90_len ( string ) |
---|
1417 | ! |
---|
1418 | ! !INPUT PARAMETERS: |
---|
1419 | ! |
---|
1420 | character*(*) string ! a string |
---|
1421 | ! |
---|
1422 | ! !OUTPUT PARAMETERS: |
---|
1423 | ! |
---|
1424 | ! The length of the string, excluding trailing blanks. |
---|
1425 | ! |
---|
1426 | ! !REVISION HISTORY: |
---|
1427 | ! |
---|
1428 | ! 01Apr94 Guo Original code (a.k.a. luavail()) |
---|
1429 | ! 19Jun96 da Silva Minor modification + prologue. |
---|
1430 | ! |
---|
1431 | !EOP |
---|
1432 | !------------------------------------------------------------------------- |
---|
1433 | |
---|
1434 | integer ls, i, l |
---|
1435 | ls = len(string) |
---|
1436 | do i = ls, 1, -1 |
---|
1437 | l = i |
---|
1438 | if ( string(i:i) .ne. BLK ) go to 11 |
---|
1439 | end do |
---|
1440 | l = l - 1 |
---|
1441 | 11 continue |
---|
1442 | i90_len = l |
---|
1443 | return |
---|
1444 | end function i90_len |
---|
1445 | |
---|
1446 | !................................................................... |
---|
1447 | |
---|
1448 | integer function I90_Lua() |
---|
1449 | |
---|
1450 | implicit NONE |
---|
1451 | |
---|
1452 | !------------------------------------------------------------------------- |
---|
1453 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
1454 | !------------------------------------------------------------------------- |
---|
1455 | !BOP |
---|
1456 | ! |
---|
1457 | ! !ROUTINE: I90_Lua() --- Returns available logical unit number. |
---|
1458 | ! |
---|
1459 | ! !DESCRIPTION: |
---|
1460 | ! |
---|
1461 | ! Look for an available (not opened) Fortran logical unit for i/o. |
---|
1462 | ! |
---|
1463 | ! !CALLING SEQUENCE: |
---|
1464 | ! |
---|
1465 | ! integer lu |
---|
1466 | ! lu = i90_lua() |
---|
1467 | ! |
---|
1468 | ! !INPUT PARAMETERS: |
---|
1469 | ! |
---|
1470 | ! None. |
---|
1471 | ! |
---|
1472 | ! !OUTPUT PARAMETERS: |
---|
1473 | ! |
---|
1474 | ! The desired unit number if positive, -1 if unsucessful. |
---|
1475 | ! |
---|
1476 | ! !REVISION HISTORY: |
---|
1477 | ! |
---|
1478 | ! 01Apr94 Guo Original code (a.k.a. luavail()) |
---|
1479 | ! 19Jun96 da Silva Minor modification + prologue. |
---|
1480 | ! |
---|
1481 | !EOP |
---|
1482 | !------------------------------------------------------------------------- |
---|
1483 | |
---|
1484 | |
---|
1485 | integer lu,ios |
---|
1486 | logical opnd |
---|
1487 | lu=7 |
---|
1488 | inquire(unit=lu,opened=opnd,iostat=ios) |
---|
1489 | do while(ios.eq.0.and.opnd) |
---|
1490 | lu=lu+1 |
---|
1491 | inquire(unit=lu,opened=opnd,iostat=ios) |
---|
1492 | end do |
---|
1493 | if(ios.ne.0) lu=-1 |
---|
1494 | i90_lua=lu |
---|
1495 | return |
---|
1496 | end function i90_lua |
---|
1497 | |
---|
1498 | !................................................................... |
---|
1499 | |
---|
1500 | subroutine i90_pad ( string ) |
---|
1501 | |
---|
1502 | implicit NONE |
---|
1503 | |
---|
1504 | !------------------------------------------------------------------------- |
---|
1505 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
1506 | !------------------------------------------------------------------------- |
---|
1507 | !BOP |
---|
1508 | ! |
---|
1509 | ! !ROUTINE: I90_Pad() --- Pad strings. |
---|
1510 | ! |
---|
1511 | ! !DESCRIPTION: |
---|
1512 | ! |
---|
1513 | ! Pads from the right with the comment character (\#). It also |
---|
1514 | ! replaces TABs with blanks for convenience. This is a low level |
---|
1515 | ! i90 routine. |
---|
1516 | ! |
---|
1517 | ! !CALLING SEQUENCE: |
---|
1518 | ! |
---|
1519 | ! call i90_pad ( string ) |
---|
1520 | ! |
---|
1521 | ! !INPUT PARAMETERS: |
---|
1522 | ! |
---|
1523 | character*256 string ! input string |
---|
1524 | |
---|
1525 | ! !OUTPUT PARAMETERS: ! modified string |
---|
1526 | ! |
---|
1527 | ! character*256 string |
---|
1528 | ! |
---|
1529 | ! !BUGS: |
---|
1530 | ! |
---|
1531 | ! It alters TABs even inside strings. |
---|
1532 | ! |
---|
1533 | ! |
---|
1534 | ! !REVISION HISTORY: |
---|
1535 | ! |
---|
1536 | ! 19Jun96 da Silva Original code. |
---|
1537 | ! |
---|
1538 | !EOP |
---|
1539 | !------------------------------------------------------------------------- |
---|
1540 | |
---|
1541 | integer i |
---|
1542 | |
---|
1543 | ! Pad end of string with # |
---|
1544 | ! ------------------------ |
---|
1545 | do i = 256, 1, -1 |
---|
1546 | if ( string(i:i) .ne. ' ' .and. & |
---|
1547 | string(i:i) .ne. '$' ) go to 11 |
---|
1548 | string(i:i) = '#' |
---|
1549 | end do |
---|
1550 | 11 continue |
---|
1551 | |
---|
1552 | ! Replace TABs with blanks |
---|
1553 | ! ------------------------- |
---|
1554 | do i = 1, 256 |
---|
1555 | if ( string(i:i) .eq. TAB ) string(i:i) = BLK |
---|
1556 | if ( string(i:i) .eq. '#' ) go to 21 |
---|
1557 | end do |
---|
1558 | 21 continue |
---|
1559 | |
---|
1560 | return |
---|
1561 | end subroutine i90_pad |
---|
1562 | |
---|
1563 | !................................................................... |
---|
1564 | |
---|
1565 | subroutine I90_Trim ( string ) |
---|
1566 | |
---|
1567 | implicit NONE |
---|
1568 | |
---|
1569 | !------------------------------------------------------------------------- |
---|
1570 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
1571 | !------------------------------------------------------------------------- |
---|
1572 | !BOP |
---|
1573 | ! |
---|
1574 | ! !ROUTINE: I90_Trim() - Removes leading blanks from strings. |
---|
1575 | ! |
---|
1576 | ! !DESCRIPTION: |
---|
1577 | ! |
---|
1578 | ! Removes blanks and TABS from begenning of string. |
---|
1579 | ! This is a low level i90 routine. |
---|
1580 | ! |
---|
1581 | ! !CALLING SEQUENCE: |
---|
1582 | ! |
---|
1583 | ! call i90_Trim ( string ) |
---|
1584 | ! |
---|
1585 | ! !INPUT PARAMETERS: |
---|
1586 | ! |
---|
1587 | character*256 string ! the input string |
---|
1588 | ! |
---|
1589 | ! !OUTPUT PARAMETERS: |
---|
1590 | ! |
---|
1591 | ! character*256 string ! the modified string |
---|
1592 | ! |
---|
1593 | ! |
---|
1594 | ! !REVISION HISTORY: |
---|
1595 | ! |
---|
1596 | ! 19Jun96 da Silva Original code. |
---|
1597 | ! |
---|
1598 | !EOP |
---|
1599 | !------------------------------------------------------------------------- |
---|
1600 | |
---|
1601 | integer ib, i |
---|
1602 | |
---|
1603 | ! Get rid of leading blanks |
---|
1604 | ! ------------------------- |
---|
1605 | ib = 1 |
---|
1606 | do i = 1, 255 |
---|
1607 | if ( string(i:i) .ne. ' ' .and. & |
---|
1608 | string(i:i) .ne. TAB ) go to 21 |
---|
1609 | ib = ib + 1 |
---|
1610 | end do |
---|
1611 | 21 continue |
---|
1612 | |
---|
1613 | ! String without trailling blanks |
---|
1614 | ! ------------------------------- |
---|
1615 | string = string(ib:) |
---|
1616 | |
---|
1617 | return |
---|
1618 | end subroutine i90_trim |
---|
1619 | |
---|
1620 | |
---|
1621 | !========================================================================== |
---|
1622 | |
---|
1623 | |
---|
1624 | ! ----------------------------- |
---|
1625 | ! Inpak 77 Upward Compatibility |
---|
1626 | ! ----------------------------- |
---|
1627 | |
---|
1628 | |
---|
1629 | subroutine lablin ( label ) |
---|
1630 | |
---|
1631 | implicit NONE |
---|
1632 | |
---|
1633 | !------------------------------------------------------------------------- |
---|
1634 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
1635 | !------------------------------------------------------------------------- |
---|
1636 | !BOP |
---|
1637 | ! |
---|
1638 | ! !ROUTINE: Lablin() --- Selects a Label (Inpak 77) |
---|
1639 | ! |
---|
1640 | ! !DESCRIPTION: |
---|
1641 | ! |
---|
1642 | ! Selects a given ``line'' (record/table) associated with ``label''. |
---|
1643 | ! Similar to {\tt i90\_label()}, but prints a message to {\tt stdout} |
---|
1644 | ! if it cannot locate the label. Kept for Inpak 77 upward compatibility. |
---|
1645 | ! |
---|
1646 | ! !CALLING SEQUENCE: |
---|
1647 | ! |
---|
1648 | ! call lablin ( label ) |
---|
1649 | ! |
---|
1650 | ! !INPUT PARAMETERS: |
---|
1651 | |
---|
1652 | character(len=*),intent(in) :: label ! string with label name |
---|
1653 | ! |
---|
1654 | ! !OUTPUT PARAMETERS: |
---|
1655 | ! |
---|
1656 | ! None. |
---|
1657 | ! |
---|
1658 | ! !REVISION HISTORY: |
---|
1659 | ! |
---|
1660 | ! 19Jun96 da Silva Original code. |
---|
1661 | ! |
---|
1662 | !EOP |
---|
1663 | !------------------------------------------------------------------------- |
---|
1664 | |
---|
1665 | integer iret |
---|
1666 | |
---|
1667 | call i90_label ( label, iret ) |
---|
1668 | if ( iret .ne. 0 ) then |
---|
1669 | write(stderr,'(2a)') 'i90/lablin: cannot find label ', label |
---|
1670 | endif |
---|
1671 | |
---|
1672 | end subroutine lablin |
---|
1673 | |
---|
1674 | !................................................................... |
---|
1675 | |
---|
1676 | real(SP) function fltgetsp ( default ) |
---|
1677 | |
---|
1678 | implicit NONE |
---|
1679 | |
---|
1680 | !------------------------------------------------------------------------- |
---|
1681 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
1682 | !------------------------------------------------------------------------- |
---|
1683 | !BOP |
---|
1684 | ! |
---|
1685 | ! !ROUTINE: FltGetsp() --- Returns next float (Inpak 77, single precision) |
---|
1686 | ! |
---|
1687 | ! !DESCRIPTION: |
---|
1688 | ! |
---|
1689 | ! Returns next float (real number, single precision) from the current |
---|
1690 | ! line, or a default value if it fails to obtain the desired number. |
---|
1691 | ! Kept for Inpak 77 upward compatibility. |
---|
1692 | ! |
---|
1693 | ! !CALLING SEQUENCE: |
---|
1694 | ! |
---|
1695 | ! real rnumber, default |
---|
1696 | ! rnumber = fltgetsp ( default ) |
---|
1697 | ! |
---|
1698 | ! !INPUT PARAMETERS: |
---|
1699 | ! |
---|
1700 | real(SP), intent(IN) :: default ! default value. |
---|
1701 | |
---|
1702 | ! |
---|
1703 | ! !REVISION HISTORY: |
---|
1704 | ! |
---|
1705 | ! 19Jun96 da Silva Original code. |
---|
1706 | ! 12Oct99 Guo/Larson - Built from original FltGet() function. |
---|
1707 | ! |
---|
1708 | !EOP |
---|
1709 | !------------------------------------------------------------------------- |
---|
1710 | |
---|
1711 | character*256 token |
---|
1712 | real(FP) x |
---|
1713 | integer iret |
---|
1714 | |
---|
1715 | call i90_gtoken ( token, iret ) |
---|
1716 | if ( iret .eq. 0 ) then |
---|
1717 | read(token,*,iostat=iret) x |
---|
1718 | end if |
---|
1719 | if ( iret .ne. 0 ) x = default |
---|
1720 | !print *, x |
---|
1721 | fltgetsp = x |
---|
1722 | |
---|
1723 | return |
---|
1724 | end function fltgetsp |
---|
1725 | |
---|
1726 | !................................................................... |
---|
1727 | |
---|
1728 | real(DP) function fltgetdp ( default ) |
---|
1729 | |
---|
1730 | implicit NONE |
---|
1731 | |
---|
1732 | !------------------------------------------------------------------------- |
---|
1733 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
1734 | !------------------------------------------------------------------------- |
---|
1735 | !BOP |
---|
1736 | ! |
---|
1737 | ! !ROUTINE: FltGetdp() --- Returns next float (Inpak 77) |
---|
1738 | ! |
---|
1739 | ! !DESCRIPTION: |
---|
1740 | ! |
---|
1741 | ! Returns next float (real number) from the current line, or a |
---|
1742 | ! default value (double precision) if it fails to obtain the desired |
---|
1743 | ! number. Kept for Inpak 77 upward compatibility. |
---|
1744 | ! |
---|
1745 | ! !CALLING SEQUENCE: |
---|
1746 | ! |
---|
1747 | ! real(DP) :: default |
---|
1748 | ! real :: rnumber |
---|
1749 | ! rnumber = FltGetdp(default) |
---|
1750 | ! |
---|
1751 | ! !INPUT PARAMETERS: |
---|
1752 | ! |
---|
1753 | real(DP), intent(IN) :: default ! default value. |
---|
1754 | |
---|
1755 | ! |
---|
1756 | ! !REVISION HISTORY: |
---|
1757 | ! |
---|
1758 | ! 19Jun96 da Silva Original code. |
---|
1759 | ! 12Oct99 Guo/Larson - Built from original FltGet() function. |
---|
1760 | ! |
---|
1761 | !EOP |
---|
1762 | !------------------------------------------------------------------------- |
---|
1763 | |
---|
1764 | character*256 token |
---|
1765 | real(FP) x |
---|
1766 | integer iret |
---|
1767 | |
---|
1768 | call i90_gtoken ( token, iret ) |
---|
1769 | if ( iret .eq. 0 ) then |
---|
1770 | read(token,*,iostat=iret) x |
---|
1771 | end if |
---|
1772 | if ( iret .ne. 0 ) x = default |
---|
1773 | !print *, x |
---|
1774 | fltgetdp = x |
---|
1775 | |
---|
1776 | return |
---|
1777 | end function fltgetdp |
---|
1778 | |
---|
1779 | !................................................................... |
---|
1780 | |
---|
1781 | integer function intget ( default ) |
---|
1782 | |
---|
1783 | implicit NONE |
---|
1784 | |
---|
1785 | !------------------------------------------------------------------------- |
---|
1786 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
1787 | !------------------------------------------------------------------------- |
---|
1788 | !BOP |
---|
1789 | ! |
---|
1790 | ! !ROUTINE: IntGet() --- Returns next integer (Inpak 77). |
---|
1791 | ! |
---|
1792 | ! !DESCRIPTION: |
---|
1793 | ! |
---|
1794 | ! Returns next integer number from the current line, or a default |
---|
1795 | ! value if it fails to obtain the desired number. |
---|
1796 | ! Kept for Inpak 77 upward compatibility. |
---|
1797 | ! |
---|
1798 | ! !CALLING SEQUENCE: |
---|
1799 | ! |
---|
1800 | ! integer number, default |
---|
1801 | ! number = intget ( default ) |
---|
1802 | ! |
---|
1803 | ! !INPUT PARAMETERS: |
---|
1804 | ! |
---|
1805 | integer default ! default value. |
---|
1806 | |
---|
1807 | ! |
---|
1808 | ! !REVISION HISTORY: |
---|
1809 | ! |
---|
1810 | ! 19Jun96 da Silva Original code. |
---|
1811 | ! |
---|
1812 | !EOP |
---|
1813 | !------------------------------------------------------------------------- |
---|
1814 | |
---|
1815 | character*256 token |
---|
1816 | real(FP) x |
---|
1817 | integer iret |
---|
1818 | |
---|
1819 | call i90_gtoken ( token, iret ) |
---|
1820 | if ( iret .eq. 0 ) then |
---|
1821 | read(token,*,iostat=iret) x |
---|
1822 | end if |
---|
1823 | if ( iret .ne. 0 ) x = default |
---|
1824 | intget = nint(x) |
---|
1825 | !print *, intget |
---|
1826 | |
---|
1827 | return |
---|
1828 | end function intget |
---|
1829 | |
---|
1830 | !................................................................... |
---|
1831 | |
---|
1832 | character(len=1) function chrget ( default ) |
---|
1833 | |
---|
1834 | implicit NONE |
---|
1835 | |
---|
1836 | !------------------------------------------------------------------------- |
---|
1837 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
1838 | !------------------------------------------------------------------------- |
---|
1839 | !BOP |
---|
1840 | ! |
---|
1841 | ! !ROUTINE: ChrGet() --- Returns next character (Inpak 77). |
---|
1842 | ! |
---|
1843 | ! !DESCRIPTION: |
---|
1844 | ! |
---|
1845 | ! Returns next non-blank character from the current line, or a default |
---|
1846 | ! character if it fails for whatever reason. |
---|
1847 | ! Kept for Inpak 77 upward compatibility. |
---|
1848 | ! |
---|
1849 | ! !CALLING SEQUENCE: |
---|
1850 | ! |
---|
1851 | ! character*1 ch, default |
---|
1852 | ! ch = chrget ( default ) |
---|
1853 | ! |
---|
1854 | ! !INPUT PARAMETERS: |
---|
1855 | ! |
---|
1856 | character*1 default ! default value. |
---|
1857 | |
---|
1858 | ! |
---|
1859 | ! !REVISION HISTORY: |
---|
1860 | ! |
---|
1861 | ! 19Jun96 da Silva Original code. |
---|
1862 | ! |
---|
1863 | !EOP |
---|
1864 | !------------------------------------------------------------------------- |
---|
1865 | |
---|
1866 | character*256 token |
---|
1867 | integer iret |
---|
1868 | |
---|
1869 | call i90_gtoken ( token, iret ) |
---|
1870 | if ( iret .ne. 0 ) then |
---|
1871 | chrget = default |
---|
1872 | else |
---|
1873 | chrget = token(1:1) |
---|
1874 | end if |
---|
1875 | !print *, chrget |
---|
1876 | |
---|
1877 | return |
---|
1878 | end function chrget |
---|
1879 | |
---|
1880 | !................................................................... |
---|
1881 | |
---|
1882 | subroutine TokGet ( token, default ) |
---|
1883 | |
---|
1884 | implicit NONE |
---|
1885 | |
---|
1886 | |
---|
1887 | !------------------------------------------------------------------------- |
---|
1888 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
1889 | !------------------------------------------------------------------------- |
---|
1890 | !BOP |
---|
1891 | ! |
---|
1892 | ! !ROUTINE: TokGet() --- Gets next token (Inpakk 77 like). |
---|
1893 | ! |
---|
1894 | ! !DESCRIPTION: |
---|
1895 | ! |
---|
1896 | ! Returns next token from the current line, or a default |
---|
1897 | ! word if it fails for whatever reason. |
---|
1898 | ! |
---|
1899 | ! !CALLING SEQUENCE: |
---|
1900 | ! |
---|
1901 | ! call TokGet ( token, default ) |
---|
1902 | ! |
---|
1903 | ! !INPUT PARAMETERS: |
---|
1904 | ! |
---|
1905 | character*(*) default ! default token |
---|
1906 | |
---|
1907 | ! !OUTPUT PARAMETERS: |
---|
1908 | ! |
---|
1909 | character*(*) token ! desired token |
---|
1910 | ! |
---|
1911 | ! !REVISION HISTORY: |
---|
1912 | ! |
---|
1913 | ! 19Jun96 da Silva Original code. |
---|
1914 | ! |
---|
1915 | !EOP |
---|
1916 | !------------------------------------------------------------------------- |
---|
1917 | |
---|
1918 | integer iret |
---|
1919 | |
---|
1920 | call i90_GToken ( token, iret ) |
---|
1921 | if ( iret .ne. 0 ) then |
---|
1922 | token = default |
---|
1923 | end if |
---|
1924 | !print *, token |
---|
1925 | |
---|
1926 | return |
---|
1927 | end subroutine tokget |
---|
1928 | |
---|
1929 | !==================================================================== |
---|
1930 | |
---|
1931 | ! -------------------------- |
---|
1932 | ! Obsolete Inpak 77 Routines |
---|
1933 | ! (Not Documented) |
---|
1934 | ! -------------------------- |
---|
1935 | |
---|
1936 | !................................................................... |
---|
1937 | |
---|
1938 | subroutine iniin() |
---|
1939 | print *, & |
---|
1940 | 'i90: iniin() is obsolete, use i90_loadf() instead!' |
---|
1941 | return |
---|
1942 | end subroutine iniin |
---|
1943 | |
---|
1944 | |
---|
1945 | !................................................................... |
---|
1946 | |
---|
1947 | subroutine iunits ( mifans, moftrm, moferr, miftrm ) |
---|
1948 | integer mifans, moftrm, moferr, miftrm |
---|
1949 | print *, & |
---|
1950 | 'i90: iunits() is obsolete, use i90_loadf() instead!' |
---|
1951 | return |
---|
1952 | end subroutine iunits |
---|
1953 | |
---|
1954 | !................................................................... |
---|
1955 | |
---|
1956 | subroutine getstr ( iret, string ) |
---|
1957 | implicit NONE |
---|
1958 | character*(*) string |
---|
1959 | integer iret !, ls |
---|
1960 | call i90_gstr ( string, iret ) |
---|
1961 | return |
---|
1962 | end subroutine getstr |
---|
1963 | |
---|
1964 | !................................................................... |
---|
1965 | |
---|
1966 | subroutine getwrd ( iret, word ) |
---|
1967 | implicit NONE |
---|
1968 | character*(*) word |
---|
1969 | integer iret |
---|
1970 | call i90_gtoken ( word, iret ) |
---|
1971 | return |
---|
1972 | end subroutine getwrd |
---|
1973 | |
---|
1974 | !................................................................... |
---|
1975 | |
---|
1976 | subroutine rdnext ( iret ) |
---|
1977 | implicit NONE |
---|
1978 | integer iret |
---|
1979 | call i90_gline ( iret ) |
---|
1980 | return |
---|
1981 | end subroutine rdnext |
---|
1982 | |
---|
1983 | !................................................................... |
---|
1984 | |
---|
1985 | real(FP) function str2rn ( string, iret ) |
---|
1986 | implicit NONE |
---|
1987 | character*(*) string |
---|
1988 | integer iret |
---|
1989 | read(string,*,end=11,err=11) str2rn |
---|
1990 | iret = 0 |
---|
1991 | return |
---|
1992 | 11 iret = 1 |
---|
1993 | return |
---|
1994 | end function str2rn |
---|
1995 | |
---|
1996 | !................................................................... |
---|
1997 | |
---|
1998 | subroutine strget ( string, default ) |
---|
1999 | |
---|
2000 | implicit NONE |
---|
2001 | |
---|
2002 | !------------------------------------------------------------------------- |
---|
2003 | ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS ! |
---|
2004 | !------------------------------------------------------------------------- |
---|
2005 | ! |
---|
2006 | ! !ROUTINE: StrGet() |
---|
2007 | ! |
---|
2008 | ! !DESCRIPTION: |
---|
2009 | ! |
---|
2010 | ! Returns next string on the current line, or a default |
---|
2011 | ! string if it fails for whatever reason. Similar to {\tt i90\_gstr()}. |
---|
2012 | ! Kept for Inpak 77 upward compatibility. |
---|
2013 | ! |
---|
2014 | ! NOTE: This is an obsolete routine. The notion of "string" used |
---|
2015 | ! here is not conventional. Please use routine {\tt TokGet()} |
---|
2016 | ! instead. |
---|
2017 | ! |
---|
2018 | ! !CALLING SEQUENCE: |
---|
2019 | ! |
---|
2020 | ! call strget ( string, default ) |
---|
2021 | ! |
---|
2022 | ! !INPUT PARAMETERS: |
---|
2023 | ! |
---|
2024 | character*(*) default ! default string |
---|
2025 | |
---|
2026 | ! !OUTPUT PARAMETERS: |
---|
2027 | |
---|
2028 | character*(*) string ! desired string |
---|
2029 | |
---|
2030 | ! |
---|
2031 | ! !REVISION HISTORY: |
---|
2032 | ! |
---|
2033 | ! 19Jun96 da Silva Original code. |
---|
2034 | ! 01Oct96 Jing Guo Removed the null terminitor |
---|
2035 | ! |
---|
2036 | !------------------------------------------------------------------------- |
---|
2037 | |
---|
2038 | integer iret |
---|
2039 | |
---|
2040 | call i90_gstr ( string, iret ) |
---|
2041 | if ( iret .ne. 0 ) then |
---|
2042 | string = default |
---|
2043 | end if |
---|
2044 | |
---|
2045 | return |
---|
2046 | end subroutine strget |
---|
2047 | |
---|
2048 | |
---|
2049 | end module m_inpak90 |
---|