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

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

Parent Directory Parent Directory | Revision Log Revision Log


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

  ViewVC Help
Powered by ViewVC 1.1.21