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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 178 - (hide annotations)
Fri Mar 11 18:47:26 2016 UTC (8 years, 3 months ago) by guez
File size: 5652 byte(s)
Moved variables date0, deltat, datasz_max, ncvar_ids, point, buff_pos,
buffer, regular from module histcom_var to modules where they are
defined.

Removed procedure ioipslmpp, useless for a sequential program.

Added argument datasz_max to histwrite_real (to avoid circular
dependency with histwrite).

Removed useless variables and computations everywhere.

Changed real litteral constants from default kind to double precision
in lwb, lwu, lwvn, sw1s, swtt, swtt1, swu.

Removed unused arguments: paer of sw, sw1s, sw2s, swclr; pcldsw of
sw1s, sw2s; pdsig, prayl of swr; co2_ppm of clmain, clqh; tsol of
transp_lay; nsrf of screenp; kcrit and kknu of gwstress; pstd of
orosetup.

Added output of relative humidity.

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

  ViewVC Help
Powered by ViewVC 1.1.21