source: CPL/oasis3/trunk/src/mod/oasis3/src/getfld.F @ 1677

Last change on this file since 1677 was 1677, checked in by aclsce, 12 years ago

Imported oasis3 (tag ipslcm5a) from cvs server to svn server (igcmg project).

File size: 25.8 KB
Line 
1      SUBROUTINE getfld (kindex, kfield, kiter)
2C****
3C               *****************************
4C               * OASIS ROUTINE  -  LEVEL 1 *
5C               * -------------     ------- *
6C               *****************************
7C
8C**** *getfld* - reading routine
9C
10C
11C     Purpose:
12C     -------
13C     Read coupling fields for iteration kiter
14C
15C**   Interface:
16C     ---------
17C       *CALL*  *getfld (kindex, kfield, kiter)*
18C
19C     Input:
20C     -----
21C                kindex : current active fields index array
22C                kfield : current active fields total number
23C                kiter  : iteration number
24C
25C     Output:
26C     ------
27C     None
28C
29C     Workspace:
30C     ---------
31C     None
32C
33C     Externals:
34C     ---------
35C     PIPE_Recv, CLIM_Import, SVIPC_read
36C
37C     Reference:
38C     ---------
39C     See OASIS manual (1995)
40C
41C     History:
42C     -------
43C       Version   Programmer     Date      Description
44C       -------   ----------     ----      ----------- 
45C       2.0       L. Terray      95/09/01  created
46C       2.1       L. Terray      96/08/05  modified: correct printing of
47C                                          diagnostic (variable clname)
48C       2.2       S. Valcke      97/08/22  added: introduction of SVIPC
49C       2.2       L.Terray       97/09/20  added: extra test in delay mode
50C                                          + general cleaning + test on
51C                                          info mode
52C       2.3       L. Terray      99/03/01  modified: reading in delay mode
53C                                          (negative values mean one must
54C                                           read files, 0 means channel)
55C       2.3       S. Valcke      99/04/30  added: printing levels
56C       2.3       L. Terray      99/09/15  added: GMEM branch
57C
58C %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
59C
60C* -----------------Include files and USE of modules---------------------------
61C
62      USE mod_kinds_oasis
63#if defined use_comm_MPI1 || defined use_comm_MPI2 || !defined use_comm_MPI1 && !defined use_comm_MPI2 && !defined use_comm_SIPC && !defined use_comm_GMEM && !defined use_comm_PIPE && !defined use_comm_NONE
64      USE mod_clim
65#endif
66      USE mod_parameter
67      USE mod_string
68      USE mod_analysis
69      USE mod_memory
70      USE mod_sipc
71      USE mod_unitncdf
72      USE mod_experiment
73      USE mod_timestep
74      USE mod_unit
75      USE mod_hardware
76      USE mod_label
77      USE mod_printing
78      INCLUDE 'netcdf.inc'
79C
80C* ---------------------------- Argument declarations -------------------
81C
82      INTEGER (kind=ip_intwp_p) kindex(kfield)
83C
84C* ---------------------------- Local declarations ----------------------
85C
86      INTEGER (kind=ip_intwp_p), DIMENSION(:), ALLOCATABLE :: info, 
87     $    iflag
88      CHARACTER*8 clname, clfic
89      CHARACTER*53 clabel
90      INTEGER (kind=ip_intwp_p) itime(3)
91      CHARACTER*3 cljobnam
92      INTEGER (kind=ip_intwp_p) il_varid, il_dimid, il_dimlen
93      INTEGER (kind=ip_intwp_p) il_nvar, il_natt, il_att
94      INTEGER (kind=ip_intwp_p) ist(3), icnt(3)
95      INTEGER (kind=ip_intwp_p) il_indice_since, il_cf_units
96      LOGICAL l_unit_found
97      CHARACTER(len=60) cl_attname
98      CHARACTER(len=17) cl_units, cl_units_cf
99      CHARACTER(len=260) cl_atttext
100      INTEGER (kind=ip_intwp_p) , parameter :: ip_nb_cf_units = 63
101      CHARACTER(len=17), dimension(ip_nb_cf_units) :: 
102     $         cl_cf_time_units = 
103     $         (/'second           ', 'minute           ',
104     $           'hour             ', 'day              ',
105     $           's                ', 'sec              ',
106     $           'shake            ', 'sidereal_day     ',
107     $           'sidereal_hour    ', 'sidereal_minute  ',
108     $           'sidereal_second  ', 'sidereal_year    ',
109     $           'tropical_year    ', 'lunar_month      ',
110     $           'common_year      ', 'leap_year        ',
111     $           'Julian_year      ', 'Gregorian_year   ',
112     $           'sidereal_month   ', 'tropical_month   ',
113     $           'd                ', 'min              ',
114     $           'hr               ', 'h                ',
115     $           'fortnight        ', 'week             ',
116     $           'jiffy            ', 'jiffies          ',
117     $           'year             ', 'yr               ',
118     $           'a                ', 'eon              ',
119     $           'month            ',
120     $           'seconds          ', 'minutes          ',
121     $           'hours            ', 'days             ',
122     $           'secs             ',
123     $           'shakes           ', 'sidereal_days    ',
124     $           'sidereal_hours   ', 'sidereal_minutes ',
125     $           'sidereal_seconds ', 'sidereal_years   ',
126     $           'tropical_years   ', 'lunar_months     ',
127     $           'common_years     ', 'leap_years       ',
128     $           'Julian_years     ', 'Gregorian_years  ',
129     $           'sidereal_months  ', 'tropical_months  ',
130     $           'ds               ', 'mins             ',
131     $           'hrs              ', 'hs               ',
132     $           'fortnights       ', 'weeks            ',
133     $           'years            ', 'yrs              ',
134     $           'as               ', 'eons             ',
135     $           'months           '/)
136C
137C* ---------------------------- Poema verses ----------------------------
138C
139C %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
140C
141C*    1. Allocation and initialization
142C        -----------------------------
143C
144      IF (nlogprt .GE. 2) THEN
145          WRITE (UNIT = nulou,FMT = *) ' '
146          WRITE (UNIT = nulou,FMT = *) ' '
147          WRITE (UNIT = nulou,FMT = *) 
148     $    '           ROUTINE getfld  -  Level 1'
149          WRITE (UNIT = nulou,FMT = *) 
150     $    '           **************     *******'
151          WRITE (UNIT = nulou,FMT = *) ' '
152          WRITE (UNIT = nulou,FMT = *) ' Get coupling fields'
153          WRITE (UNIT = nulou,FMT = *) ' '
154          WRITE (UNIT = nulou,FMT = *) ' '
155      ENDIF
156      imrca = 0
157      imrcb = 0
158      itime(:)=0
159      ist(:)=0
160      icnt(:)=0
161      ALLOCATE(iflag(ig_nfield))
162      ALLOCATE(info(ig_nfield))
163      CALL izero (iflag, ig_nfield)
164      CALL izero (info, ig_nfield)
165C
166C*    2. Loop on active fields for iteration kiter
167C        -----------------------------------------
168!$omp parallel do default (shared)
169!$omp+ private(jf,iloc,iadrold,isizold,ilabel)
170!$omp+ private(clabel,clname,clfic,iunit)
171!$omp+ private(iseqn,il_varid,n_reaty)
172
173C
174      DO 210 jf = 1, kfield
175C
176C* Assign local variables
177C
178         iloc = kindex(jf)
179         iadrold = nadrold(iloc)
180         isizold = nsizold(iloc)
181         ilabel = numlab(iloc)
182         clabel = cfldlab(ilabel)
183         clname = cnaminp(iloc)
184         clfic = cficinp(iloc)
185         iunit = nluinp(iloc)
186         iseqn = nseqn(iloc)
187         ilagn = nlagn(iloc)
188C
189C* Print field name
190C
191         IF (nlogprt .GE. 1) THEN
192            CALL prcout('Reading of field : ', clname, 2)
193         ENDIF
194         IF (nlogprt .GE. 2) THEN
195            CALL prcout('Field long name : ', clabel, 2)
196         ENDIF
197C
198C* - Get coupling fields
199C
200C* Specific treatment for first iteration 
201C
202#if defined use_comm_NONE || defined use_comm_SIPC
203         ilagn = 1
204#endif
205         IF (kiter .EQ. 0 .AND. ilagn .GT. 0) THEN
206C
207C* Reconnect data file to appropriate unit
208C
209            IF (lncdfrst) THEN
210               ist(1)=1 
211               ist(2)=1
212               ist(3)=kiter+1
213               icnt(1)=nlonbf(iloc) 
214               icnt(2)=nlatbf(iloc)
215               icnt(3)=1
216C     
217               CALL hdlerr
218     $          (NF_INQ_VARID(nc_inpid(iloc),clname,il_varid),'getfld')
219#if defined use_comm_NONE
220               nc_invartimeid=0
221               istat = NF_INQ_VARID(nc_inpid(iloc),'time_counter',
222     $              nc_invartimeid ) 
223               IF (istat .ne. NF_NOERR) then
224                   istat = NF_INQ_VARID(nc_inpid(iloc),'time',
225     $                                  nc_invartimeid)
226               ENDIF
227               IF (istat .ne. NF_NOERR) then
228                   istat = NF_INQ_VARID(nc_inpid(iloc),'TIME_COUNTER',
229     $                                  nc_invartimeid)
230               ENDIF
231               IF (istat .ne. NF_NOERR) then
232                   istat = NF_INQ_VARID(nc_inpid(iloc),'TIME',
233     $                                  nc_invartimeid)
234               ENDIF
235               IF (istat .ne. NF_NOERR) then
236                   istat = NF_INQ_UNLIMDIM(nc_inpid(iloc),
237     $                     il_dimid)
238                   IF (istat .eq. NF_NOERR) then
239                     istat = NF_INQ_DIMNAME(nc_inpid(iloc),il_dimid,
240     $                       nc_invartime_name)
241                       IF (istat .eq. NF_NOERR) then
242                           istat = NF_INQ_VARID(nc_inpid(iloc),
243     $                             nc_invartime_name, nc_invartimeid)
244                       ENDIF
245                   ENDIF
246               ENDIF
247               l_unit_found=.FALSE.
248               IF (istat .ne. NF_NOERR) then
249                   CALL hdlerr(NF_INQ_NVARS(nc_inpid(iloc),il_nvar),
250     $                         'getfld')
251                  DO nc_invartimeid=1,il_nvar
252                    CALL hdlerr(NF_INQ_VARNATTS(nc_inpid(iloc),
253     $                             nc_invartimeid, il_natt),'getfld')
254                    DO il_att=1,il_natt
255                       CALL hdlerr(NF_INQ_ATTNAME(nc_inpid(iloc),
256     $                          nc_invartimeid, il_att, cl_attname),
257     $                          'getfld')
258                       IF ( cl_attname .eq. 'units' ) THEN
259                         CALL hdlerr(NF_GET_ATT_TEXT(
260     $                          nc_inpid(iloc),nc_invartimeid ,'units',
261     $                          cl_atttext),'getfld')
262                         il_indice_since=INDEX(cl_atttext,'since')
263                         IF ( il_indice_since .ne. 0 ) THEN
264                            cl_units=cl_atttext(1:il_indice_since-1)
265                            cl_units=ADJUSTL(cl_units)
266                            cl_units=TRIM(cl_units)
267                            DO il_cf_units =1, ip_nb_cf_units
268                               cl_units_cf=cl_cf_time_units(il_cf_units)
269                              IF ( cl_units .EQ. cl_units_cf )
270     $                           l_unit_found=.TRUE.
271                              CALL case_change(cl_units_cf,1)
272                              IF ( cl_units .EQ. cl_units_cf)
273     $                           l_unit_found=.TRUE.
274                            ENDDO
275                            IF ( l_unit_found )  THEN
276                              istat = NF_NOERR
277                              GO TO 230
278                            ENDIF
279                         ENDIF
280                       ENDIF
281                    ENDDO
282                  ENDDO
283 230              CONTINUE
284               ENDIF
285C
286               IF (istat .ne. NF_NOERR) then
287                  WRITE(nulou,*) 'WARNING: Time variable not found'
288                  WRITE(nulou,*) ''
289                  IF ( nitfn .gt. 0) THEN
290                     WRITE(nulou,*) 'Check coherence between'
291                     WRITE(nulou,*) 'namcouple timestep number'
292                     WRITE(nulou,*) 'and input file dimension'
293                     WRITE (nulou, *)
294                     CALL HALTE('STOP in getfld')
295                  END IF
296             WRITE(nulou,*) 'Interpolation only for one time iteration'
297                  nc_invartime_name = 'time_counter'
298               ELSE
299                  CALL hdlerr
300     $              (NF_INQ_VARNAME(nc_inpid(iloc),nc_invartimeid,
301     $               nc_invartime_name),'getfld')
302C
303C                 Check coherence 
304C                 between namcouple timestep number
305C                 and input file dimension
306C   
307                  CALL hdlerr
308     $              (NF_INQ_VARDIMID(nc_inpid(iloc),nc_invartimeid,
309     $               il_dimid),'getfld')
310                  CALL hdlerr
311     $              (NF_INQ_DIMLEN(nc_inpid(iloc),il_dimid,
312     $               il_dimlen),'getfld')
313                  IF ( il_dimlen .le. nitfn ) then
314                      WRITE(nulou,*) 'Incoherence between'
315                      WRITE(nulou,*) 'namcouple timestep number'
316                      WRITE(nulou,*) 'and input file dimension'
317                      WRITE(nulou,*) ''
318                      CALL HALTE('STOP in getfld')
319                  ENDIF
320               ENDIF
321#else
322               nc_invartime_name = 'time_counter'
323#endif
324               CALL hdlerr
325     $              (NF_INQ_VARTYPE(nc_inpid(iloc),il_varid,n_reaty),
326     $              'getfld')
327               IF (n_reaty .eq. NF_FLOAT) THEN
328c                  IF (kind(rl_testvar). ne. 4) THEN
329                  IF ( ip_realwp_p .ne. ip_single_p ) THEN
330                     WRITE (nulou, *)
331     $         'Incoherence: Oasis is compiled in double precision'
332                     WRITE (nulou, *)
333     $      'but Restart auxiliary file contains single precision REAL'
334                     CALL HALTE('STOP in getfld')
335                  ENDIF
336                  CALL hdlerr(NF_GET_VARA_REAL(nc_inpid(iloc), 
337     $                 il_varid, ist, icnt, 
338     $                fldold(iadrold:iadrold+isizold-1)), 'getfld')
339C
340#if defined use_comm_NONE
341                  IF (nc_invartimeid .ne. 0) then
342                      istat = NF_GET_VARA_REAL(nc_inpid(iloc),
343     $                     nc_invartimeid, ist(3), icnt(3), rtime_val)
344                      IF (istat .eq. NF_NOERR) then
345                         istatus=NF_GET_VARA_REAL(nc_inpid(iloc),
346     $                    nc_invartimeid, ist(3), icnt(3), rtime_val)
347                      ENDIF
348                  ENDIF
349C
350#endif
351               ELSE IF (n_reaty .eq. NF_DOUBLE) THEN
352c                  IF (kind(rl_testvar). ne. 8) THEN
353                  IF ( ip_realwp_p .ne. ip_double_p ) THEN
354                     WRITE (nulou, *)
355     $          'Incoherence: Oasis is compiled in single precision'
356                     WRITE (nulou, *)
357     $    'but Restart auxiliary file contains double precision REAL'
358                     CALL HALTE('STOP in getfld')
359                  ENDIF               
360                  CALL hdlerr(NF_GET_VARA_DOUBLE(nc_inpid(iloc), 
361     $                 il_varid, ist, icnt, 
362     $                fldold(iadrold:iadrold+isizold-1) ), 'getfld')
363C
364#if defined use_comm_NONE
365                  IF (nc_invartimeid .ne. 0) then
366                      istat = NF_GET_VARA_DOUBLE(nc_inpid(iloc),
367     $                     nc_invartimeid, ist(3), icnt(3), rtime_val)
368                      IF (istat .eq. NF_NOERR) then
369                         istatus=NF_GET_VARA_DOUBLE(nc_inpid(iloc),
370     $                    nc_invartimeid, ist(3), icnt(3), rtime_val)
371                      ENDIF
372                  ENDIF
373C
374#endif
375               ENDIF
376            ELSE
377cSL
378C
379C* Reconnect data file to appropriate unit; this first close is necessary !
380C
381                    CLOSE(iunit)
382                    OPEN (UNIT=iunit,FILE=clfic,STATUS ='UNKNOWN',
383     $                      FORM ='UNFORMATTED',IOSTAT = iost)
384                    IF (iost .eq. 0) then
385                        IF (nlogprt .GE. 2) THEN
386                            WRITE(UNIT = nulou, FMT = 2100)iunit,clfic
387                        ENDIF
388                    ELSE
389                        CALL prtout
390     $                  ('No binary transfert file for field ',iloc,1)
391                        CALL HALTE('STOP in getfld') 
392                    ENDIF
393cSL
394               IF (lmodinf) THEN
395                  CALL locreadh (clname, cljobnam, itime, 
396     $                 fldold(iadrold), isizold, iunit, iflag(iloc)) 
397               ELSE
398                  CALL locread (clname, fldold(iadrold), isizold,
399     $                 iunit, iflag(iloc))
400               ENDIF
401               IF (iflag(iloc) .NE. 0) THEN
402                  CALL prcout ('WARNING: problem in reading field',
403     $                 clname, 1)
404                  CALL prtout ('Error reading logical unit', iunit, 2)
405                  CALL HALTE('STOP in getfld')
406               ENDIF
407C
408            ENDIF
409C
410C* Next iterations
411C
412         ELSE
413C
414C* PIPE and NONE case
415C
416#ifdef use_comm_PIPE
417C
418C* Wait for message
419C
420            CALL PIPE_Recv (clname, kiter)
421#endif
422#if defined use_comm_PIPE || defined use_comm_NONE
423C
424C* Reconnect data file to appropriate unit IF not already done
425C
426            IF (lncdfrst) THEN
427               CALL hdlerr(NF_OPEN(clfic,NF_NOWRITE,nc_inpid(iloc)),
428     $             'getfld')
429               ist(1)=1 
430               ist(2)=1
431               ist(3)=1+kiter
432               icnt(1)=nlonbf(iloc)
433               icnt(2)=nlatbf(iloc)
434               icnt(3)=1
435               istatus=NF_INQ_VARID(nc_inpid(iloc),clname,il_varid)
436               
437               istatus=NF_INQ_VARID(nc_inpid(iloc),nc_invartime_name,
438     $                        nc_invartimeid)
439               IF (n_reaty .eq. NF_FLOAT) THEN
440                  CALL hdlerr(NF_GET_VARA_REAL(nc_inpid(iloc), 
441     $                 il_varid,ist,icnt,
442     $                 fldold(iadrold:iadrold+isizold-1)),'getfld') 
443                  IF (nc_invartimeid .ne. 0) 
444     $            CALL hdlerr(NF_GET_VARA_REAL(nc_inpid(iloc),
445     $                 nc_invartimeid, ist(3), icnt(3), rtime_val)
446     $                 , 'getfld') 
447               ELSE IF (n_reaty .eq. NF_DOUBLE) THEN
448                  CALL hdlerr(NF_GET_VARA_DOUBLE(nc_inpid(iloc), 
449     $                 il_varid,ist,icnt,
450     $                 fldold(iadrold:iadrold+isizold-1)),'getfld')
451                  IF (nc_invartimeid .ne. 0) 
452     $               CALL hdlerr(NF_GET_VARA_DOUBLE(nc_inpid(iloc),
453     $                 nc_invartimeid, ist(3), icnt(3), rtime_val)
454     $                 , 'getfld')
455               ENDIF
456               CALL hdlerr(NF_CLOSE(nc_inpid(iloc)), 'getfld')
457            ELSE
458               OPEN (UNIT=iunit,FILE=clfic,STATUS ='UNKNOWN',
459     $              FORM ='UNFORMATTED',IOSTAT = iost)
460               IF (iost .eq. 0) then
461                  IF (nlogprt .GE. 2) THEN
462                     WRITE(UNIT = nulou, FMT = 2100)iunit,clfic
463                  ENDIF
464               ELSE
465                  CALL prtout
466     $                 ('No binary transfert file for field ',iloc,1)
467                  CALL HALTE('STOP in getfld') 
468               ENDIF
469               IF (lmodinf) THEN
470                  CALL locreadh (clname, cljobnam, itime, 
471     $                 fldold(iadrold),isizold,iunit,iflag(iloc)) 
472               ELSE
473                  CALL locread (clname, fldold(iadrold), 
474     $                 isizold, iunit, iflag(iloc))
475               ENDIF
476               IF (iflag(iloc) .NE. 0) THEN
477                  clname = cnaminp(iloc)
478                  CALL prcout('Problem reading field',clname,1)
479                  CALL prtout ('Error reading lunit', iunit, 2)
480                  CALL HALTE ('STOP in getfld')
481               ENDIF
482               CLOSE (UNIT = iunit, IOSTAT = ios)
483               IF (ios .EQ. 0) THEN
484                  IF (nlogprt .GE. 2) THEN
485                     WRITE(UNIT=nulou,FMT=2200) iunit, clfic
486                  ENDIF
487               ELSE
488                  CALL prtout('Problem in closing unit',iunit,2)
489                  CALL prtout('Error message nbrr is= ', ios, 2)
490                  CALL HALTE('STOP in getfld')
491               ENDIF
492            ENDIF
493C
494C* SIPC or GMEM case
495C
496#elif defined use_comm_SIPC || defined use_comm_GMEM
497C
498            ipbytecha=kind('A')
499            ipbyteint=kind(itime)
500            ipbyterea=kind(fldold)
501C
502C* Read encapsulated infos in field-specific shared memory pool
503C  (experiment name, initial date, time since start, iteration number)
504C
505            IF (lmodinf) THEN
506               isizeinp = 3*ipbytecha
507               CALL SVIPC_read(mpoolidin(iloc), cljobnam 
508     $              , isizeinp, imrca)
509               isizeinp = 3*ipbyteint
510               CALL SVIPC_read(mpoolidin(iloc), itime 
511     $              , isizeinp, imrcb)
512C
513C* Find error if any
514C
515               IF (imrca .LT. 0 .OR. imrcb .LT. 0) THEN
516                  CALL prcout
517     $   ('Problem in reading encapsulated infos for field', clname, 1)
518                  CALL HALTE ('STOP in getfld')
519               ELSE
520                  CALL prcout
521     $   ('Read encapsulated infos in SHM pool for field', clname, 1)
522               ENDIF
523            ENDIF
524C
525C* Read part of macro array in field-specific shared memory pool
526C
527            isizeinp = isizold * ipbyterea
528            WRITE(UNIT = nulou, FMT = *) 
529     $           'Reading field data from pool = ',mpoolidin(iloc) 
530            CALL SVIPC_read(mpoolidin(iloc),
531     $           fldold(iadrold), isizeinp, imrc)
532C
533C* Find error if any
534C
535            IF (imrc .LT. 0) THEN
536               CALL prcout
537     $           ('Problem in reading the field in SHM pool:',clname,1)
538               CALL HALTE ('STOP in getfld')
539            ELSE
540               CALL prcout
541     $            ('Read field in SHM pool:', clname, 1)
542            ENDIF
543C
544C* CLIM case
545C
546#elif defined use_comm_MPI1 || defined use_comm_MPI2 || !defined use_comm_MPI1 && !defined use_comm_MPI2 && !defined use_comm_SIPC && !defined use_comm_GMEM && !defined use_comm_PIPE && !defined use_comm_NONE
547!$omp critical
548            CALL CLIM_Import 
549     $           (ig_portin_id(iloc), kiter*nstep, fldold(iadrold), 
550     $           info(iloc))
551!$omp end critical
552            IF (info(iloc) .NE. CLIM_Ok) THEN
553               CALL prcout('Problem in reading field on port',
554     $              clname, 1)
555               CALL prtout('error code number', info(iloc), 2)
556               CALL HALTE ('STOP in getfld')
557            ENDIF
558#else
559            CALL prcout 
560     $           ('Wrong channel option for field', clname, 1)
561#endif
562         IF (lmodinf .AND. nlogprt .GE. 2) THEN
563            WRITE (UNIT = nulou,FMT = *) 
564     $           ' Encapsulated data for current field is :' 
565            CALL prcout ('Experiment name', cljobnam, 1)
566            CALL prtout ('Initial date', itime(1), 2)
567            CALL prtout ('Iteration number', itime(2), 2)
568            CALL prtout ('Time since start', itime(3), 2)
569         ENDIF
570      ENDIF
571 210  CONTINUE
572C
573C* Close units or netcdf restart files if first iteration
574C
575      IF (kiter .EQ. 0) THEN
576          DO 220 jf = 1, kfield
577            iloc = kindex(jf)
578            ilagn = ig_lag(iloc)
579            IF (ilagn .GT. 0) THEN
580            isamefic=0
581            DO 225 jj = 1, jf-1
582              IF (nluinp(jf) .eq. nluinp(jj)) isamefic=isamefic+1
583 225        END DO
584            IF (isamefic .lt. 1) THEN
585                IF (lncdfrst) THEN
586                    CALL hdlerr(NF_CLOSE(nc_inpid(jf)), 'getfld')
587                ELSE
588                    CLOSE (UNIT = nluinp(jf), IOSTAT = ios)
589                    IF (ios .NE. 0) THEN
590                        CALL prtout('Problem in closing unit',iunit,2)
591                        CALL prtout('Error message number is = ',ios,2)
592                        CALL HALTE('STOP in getfld')
593                    ELSE
594                        IF (nlogprt .GE. 2) THEN
595                            WRITE(UNIT=nulou, FMT=2200) iunit, clfic
596                        ENDIF
597                    ENDIF
598                ENDIF
599            ENDIF
600            ENDIF
601 220      CONTINUE
602      ENDIF
603C
604C* Formats
605C
606 2100 FORMAT(/,5X,' Unit ',I2,' has been connected to file ',A8)
607 2200 FORMAT(/,5X,' Unit ',I2,' has been disconnected from file ',A8)
608C
609C
610C*    3. End of routine
611C        --------------
612C
613      IF (nlogprt .GE. 2) THEN
614          WRITE (UNIT = nulou,FMT = *) ' '
615          WRITE (UNIT = nulou,FMT = *) 
616     $    '          --------- End of routine getfld ---------'
617          CALL FLUSH (nulou)
618      ENDIF
619      DEALLOCATE(iflag)
620      DEALLOCATE(info)
621      RETURN
622      END
623C
624C====================================================================
625C
626          subroutine case_change( string, type )
627C
628C Purpose
629C _______
630C
631C   This converts each lower case alphabetic letter in STRING to upper
632C   case, or vice versa.
633C
634C
635C Arguments
636C _________
637C
638C   string  (INPUT/OUTPUT) character(len=*)
639C           The string to convert
640C
641C   TYPE    (INPUT/OUTPUT) integer(i1b)
642C           Define the conversion. Specifically,
643C
644C               If TYPE = 1 = TOUPPER,    conversion is lower to upper
645C               If TYPE = 2 = TOLOWER,    conversion is upper to lower
646C               If TYPE = 3 = CAPITALIZE, use upper for first letter; lower for rest
647C
648C           Definitions of TOUPPER, TOLOWER and CAPITALIZE may be obtained from
649C           the host module Strings.
650C
651C
652C Further Details
653C _______________
654C
655C   All non-alphabetic characters are left unchanged.
656C   It uses the underlying machine collating sequence.
657C
658C     History:
659C     -------
660C       Version   Programmer         Date      Description
661C       -------   ----------         ----      -----------
662C       3.0       P. Terray (LODyC)  03/12/19  created
663C
664C _________________________________________________________________________________
665C
666C
667C
668C SPECIFICATIONS FOR ARGUMENTS
669C ____________________________
670C
671          character(len=*), intent(inout) :: string
672C
673          integer,     intent(in)    :: type
674C
675C
676C SPECIFICATIONS FOR LOCAL VARIABLES
677C __________________________________
678C
679         integer :: i
680         integer, parameter :: shift = ichar( 'a' ) - ichar( 'A' )
681C
682C
683C EXECUTABLE STATEMENTS
684C _____________________
685C
686        if ( type==1 ) then
687            do i = 1, len( string )
688                select case( string(i:i) )
689                  case( 'a':'i' , 'j':'r' , 's':'z' )
690                    string(i:i) = char( ichar(string(i:i)) - shift  )
691                end select
692            end do
693        else
694            do i = 1, len( string )
695                select case( string(i:i) )
696                  case( 'A':'I' , 'J':'R' , 'S':'Z' )
697                    string(i:i) = char( ichar(string(i:i)) + shift  )
698                end select
699            end do
700            if ( type==3 ) then
701                select case( string(1:1) )
702                  case( 'a':'i' , 'j':'r' , 's':'z' )
703                    string(1:1) = char( ichar(string(1:1)) - shift  )
704                end select
705            end if
706        end if
707        return
708C
709C
710C END OF SUBROUTINE case_change
711C _____________________________
712C
713        end
714C
715
Note: See TracBrowser for help on using the repository browser.