/[lmdze]/trunk/Sources/IOIPSL/Histcom/histbeg_totreg.f
ViewVC logotype

Contents of /trunk/Sources/IOIPSL/Histcom/histbeg_totreg.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 138 - (show annotations)
Fri May 22 23:13:19 2015 UTC (9 years ago) by guez
File size: 5642 byte(s)
Moved variable nb_files from module histcom_var to module
histbeg_totreg_m.

Removed unused argument q of writehist.

No history file is created in program ce0l so there is no need to call
histclo in etat0.

In phyredem, access variables rlat and rlon directly from module
phyetat0_m instead of having them as arguments. This is clearer for
the program gcm. There are bad side effects for the program ce0l: we
have to modify the module variables rlat and rlon in procedure etat0,
and we need the additional file phyetat0.f to compile ce0l.

1 MODULE histbeg_totreg_m
2
3 ! From histcom.f90, version 2.1 2004/04/21 09:27:10
4
5 ! Some confusing vocabulary in this code!
6
7 ! A regular grid is a grid with i, j indices and thus it is
8 ! stored in a 2D matrix. This is opposed to an irregular grid which
9 ! is in a vector and where we do not know which neighbours we
10 ! have. As a consequence we need the bounds for each grid-cell.
11
12 ! A rectilinear grid is a special case of a regular grid in which
13 ! all longitudes for i constant are equal and all latitudes for j
14 ! constant are equal. In other words we do not need the full 2D
15 ! matrix to describe the grid, just two vectors.
16
17 IMPLICIT NONE
18
19 INTEGER:: nb_files = 0
20
21 CONTAINS
22
23 SUBROUTINE histbeg_totreg(filename, lon_1d, lat_1d, orix, szx, oriy, szy, &
24 pitau0, pdate0, pdeltat, horiid, fileid)
25
26 ! We assume the grid is rectilinear. The user provides "lon_1d"
27 ! and "lat_1d" as vectors. This subroutine initializes a NetCDF
28 ! file and returns the ID. It sets up the geographical space on
29 ! which the data will be stored and offers the possibility of
30 ! setting a zoom. It also gets the global parameters into the
31 ! input-output subsystem.
32
33 USE ioipslmpp, ONLY: ioipslmpp_file
34 USE errioipsl, ONLY: histerr
35 USE histcom_var, ONLY: assc_file, date0, deltat, full_size, itau0, &
36 lock_modname, model_name, nb_files_max, nb_hax, nb_tax, nb_var, &
37 nb_zax, ncdf_ids, regular, slab_ori, slab_sz, xid, yid, zoom
38 use histhori_regular_m, only: histhori_regular
39 USE netcdf, ONLY: nf90_clobber, nf90_global
40 use netcdf95, only: nf95_create, nf95_def_dim, nf95_put_att
41
42 CHARACTER(len=*), INTENT(IN):: filename
43 ! name of the netcdf file to be created
44
45 REAL, INTENT(IN):: lon_1d(:) ! coordinates of points in longitude
46 REAL, INTENT(IN):: lat_1d(:) ! coordinates of points in latitude
47
48 ! The next 4 arguments allow to define a horizontal zoom for this
49 ! file. It is assumed that all variables to come have the same
50 ! index space. This can not be assumed for the z axis and thus we
51 ! define the zoom in histdef.
52 INTEGER, INTENT(IN):: orix ! origin of the slab of data within the X axis
53 INTEGER, INTENT(IN):: szx ! size of the slab of data in X
54 INTEGER, INTENT(IN):: oriy ! origin of the slab of data within the Y axis
55 INTEGER, INTENT(IN):: szy ! size of the slab of data in Y
56
57 INTEGER, INTENT(IN):: pitau0 ! time step at which the history tape starts
58 REAL, INTENT(IN):: pdate0 ! the Julian date at which the itau was equal to 0
59 REAL, INTENT(IN):: pdeltat ! time step of the counter itau, in seconds
60
61 INTEGER, INTENT(OUT):: fileid ! ID of the netcdf file
62 INTEGER, INTENT(OUT):: horiid ! ID of the horizontal grid
63
64 ! Variables local to the procedure:
65 REAL, DIMENSION(size(lon_1d), size(lat_1d)):: lon, lat
66 INTEGER im ! size of arrays in longitude direction
67 integer jm ! size of arrays in latitude direction
68 INTEGER ncid
69 INTEGER lengf, lenga
70 CHARACTER(len=120) file
71
72 !---------------------------------------------------------------------
73
74 im = size(lon_1d)
75 jm = size(lat_1d)
76
77 lon = spread(lon_1d, 2, jm)
78 lat = spread(lat_1d, 1, im)
79
80 nb_files = nb_files + 1
81 fileid = nb_files
82
83 ! 1. Transfer into module variables for future use
84
85 itau0(fileid) = pitau0
86 date0(fileid) = pdate0
87 deltat(fileid) = pdeltat
88
89 ! 2. Initialize all variables for this file
90
91 IF (nb_files > nb_files_max) CALL histerr(3, 'histbeg', &
92 'Table of files too small. You should increase nb_files_max', &
93 'in M_HISTCOM.f90 in order to accomodate all these files', ' ')
94
95 nb_var(fileid) = 0
96 nb_tax(fileid) = 0
97 nb_hax(fileid) = 0
98 nb_zax(fileid) = 0
99
100 slab_ori(fileid, :) = (/orix, oriy/)
101 slab_sz(fileid, :) = (/szx, szy/)
102
103 ! 3. Open NetCDF file and define dimensions
104
105 lengf = len_trim(filename)
106 IF (filename(lengf-2:lengf)/='.nc') THEN
107 file = filename(:lengf) // '.nc'
108 ELSE
109 file = filename(:lengf)
110 END IF
111
112 ! Add PE number in file name on MPP
113
114 CALL ioipslmpp_file(file)
115
116 ! Keep track of the name of the files opened
117
118 lengf = len_trim(file)
119 lenga = len_trim(assc_file)
120 IF (nb_files==1) THEN
121 assc_file = file(:lengf)
122 ELSE IF ((lenga+lengf)<500) THEN
123 assc_file = assc_file(:lenga) // ' ' // file(:lengf)
124 ELSE IF (lenga + 7 < 500 .AND. index(assc_file(:lenga), 'et.al.') < 1) THEN
125 assc_file = assc_file(:lenga) // ' et.al.'
126 ELSE
127 CALL histerr(2, 'histbeg', &
128 'The file names do not fit into the associate_file attribute.', &
129 'Use shorter names if you wish to keep the information.', ' ')
130 END IF
131
132 call nf95_create(file, nf90_clobber, ncid)
133 call nf95_def_dim(ncid, 'lon', szx, xid(nb_files))
134 call nf95_def_dim(ncid, 'lat', szy, yid(nb_files))
135
136 ! 4. Declare the geographical coordinates and other attributes
137
138 ! 4.3 Global attributes
139
140 call nf95_put_att(ncid, nf90_global, 'Conventions', 'GDT 1.3')
141 call nf95_put_att(ncid, nf90_global, 'file_name', trim(file))
142 call nf95_put_att(ncid, nf90_global, 'production', trim(model_name))
143 lock_modname = .TRUE.
144
145 ! 5. Save some important information on this file in the module variables
146 ncdf_ids(fileid) = ncid
147 full_size(fileid, :) = (/im, jm/)
148
149 ! 6. Store the geographical coordinates
150
151 IF ((im /= szx) .OR. (jm /= szy)) zoom(fileid) = .TRUE.
152 regular(fileid) = .TRUE.
153
154 CALL histhori_regular(fileid, im, lon, jm, lat, ' ', 'Default grid', horiid)
155
156 END SUBROUTINE histbeg_totreg
157
158 end MODULE histbeg_totreg_m

  ViewVC Help
Powered by ViewVC 1.1.21