/[lmdze]/trunk/IOIPSL/Histcom/histclo.f90
ViewVC logotype

Contents of /trunk/IOIPSL/Histcom/histclo.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 328 - (show annotations)
Thu Jun 13 14:40:06 2019 UTC (4 years, 11 months ago) by guez
File size: 1164 byte(s)
Change all `.f` suffixes to `.f90`. (The opposite was done in revision
82.)  Because of change of philosopy in GNUmakefile: we already had a
rewritten rule for `.f`, so it does not make the makefile longer to
replace it by a rule for `.f90`. And it spares us options of
makedepf90 and of the compiler. Also we prepare the way for a simpler
`CMakeLists.txt`.

1 module histclo_m
2
3 implicit none
4
5 contains
6
7 SUBROUTINE histclo(fid)
8
9 ! This subroutine will close the file corresponding
10 ! to the argument, if the argument is present. Else it will close
11 ! all opened files.
12
13 USE errioipsl, ONLY: histerr
14 use histbeg_totreg_m, ONLY: nb_files
15 USE histcom_var, ONLY: ncdf_ids
16 USE netcdf, ONLY: nf90_close, nf90_noerr
17
18 INTEGER, INTENT (IN), OPTIONAL:: fid ! file id
19
20 ! Variables local to the procedure:
21 INTEGER ifile, ncid, iret
22 INTEGER start_loop, end_loop
23 CHARACTER(len=70) str70
24
25 !---------------------------------------------------------------------
26
27 print *, "Call sequence information: histclo"
28
29 IF (present(fid)) THEN
30 start_loop = fid
31 end_loop = fid
32 ELSE
33 start_loop = 1
34 end_loop = nb_files
35 END IF
36
37 DO ifile = start_loop, end_loop
38 ncid = ncdf_ids(ifile)
39 iret = nf90_close(ncid)
40 IF (iret/=nf90_noerr) THEN
41 WRITE(str70, '("This file has already been closed:", I3)') ifile
42 CALL histerr(2, 'histclo', str70, '', '')
43 END IF
44 END DO
45
46 END SUBROUTINE histclo
47
48 end module histclo_m

  ViewVC Help
Powered by ViewVC 1.1.21