/[lmdze]/trunk/libf/IOIPSL/Histcom/histbeg_totreg.f90
ViewVC logotype

Contents of /trunk/libf/IOIPSL/Histcom/histbeg_totreg.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 61 - (show annotations)
Fri Apr 20 14:58:43 2012 UTC (12 years ago) by guez
File size: 6073 byte(s)
No more included file in LMDZE, not even "netcdf.inc".

Created a variable containing the list of common source files in
GNUmakefile. So we now also see clearly files that are specific to
each program.

Split module "histcom". Assembled resulting files in directory
"Histcom".

Removed aliasing in calls to "laplacien".

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 which is i, j indices
8 ! and thus it is stored in a 2D matrix.
9 ! This is opposed to an IRREGULAR grid which is only in a vector
10 ! and where we do not know which neighbors we have.
11 ! As a consequence we need the bounds for each grid-cell.
12
13 ! A RECTILINEAR grid is a special case of a regular grid
14 ! in which all longitudes for i constant are equal
15 ! and all latitudes for j constant.
16 ! In other words we do not need the full 2D matrix
17 ! to describe the grid, just two vectors.
18
19 IMPLICIT NONE
20
21 CONTAINS
22
23 SUBROUTINE histbeg_totreg(pfilename, plon_1d, plat_1d, par_orix, par_szx, &
24 par_oriy, par_szy, pitau0, pdate0, pdeltat, phoriid, pfileid)
25
26 ! The user provides "plon_1d" and "plat_1d" as vectors. Obviously
27 ! this can only be used for very regular grids.
28 ! This subroutine initializes a netcdf file and returns the ID.
29 ! It will set up the geographical space on which the data will be
30 ! stored and offers the possibility of seting a zoom.
31 ! It also gets the global parameters into the I/O subsystem.
32
33 ! INPUT
34
35 ! pfilename: Name of the netcdf file to be created
36 ! pim: Size of arrays in longitude direction
37 ! plon_1d: Coordinates of points in longitude
38 ! pjm: Size of arrays in latitude direction
39 ! plat_1d: Coordinates of points in latitude
40
41 ! The next 4 arguments allow to define a horizontal zoom
42 ! for this file. It is assumed that all variables to come
43 ! have the same index space. This can not be assumed for
44 ! the z axis and thus we define the zoom in histdef.
45
46 ! par_orix: Origin of the slab of data within the X axis (pim)
47 ! par_szx: Size of the slab of data in X
48 ! par_oriy: Origin of the slab of data within the Y axis (pjm)
49 ! par_szy: Size of the slab of data in Y
50
51 ! pitau0: time step at which the history tape starts
52 ! pdate0: The Julian date at which the itau was equal to 0
53 ! pdeltat: Time step in seconds. Time step of the counter itau
54 ! used in histwrt for instance
55
56 ! OUTPUT
57
58 ! phoriid: ID of the horizontal grid
59 ! pfileid: ID of the netcdf file
60
61 ! We assume the grid is rectilinear.
62
63 USE ioipslmpp, ONLY: ioipslmpp_file
64 USE errioipsl, ONLY: histerr
65 USE histcom_var, ONLY: assc_file, date0, deltat, full_size, itau0, &
66 lock_modname, model_name, nb_files, nb_files_max, nb_hax, nb_tax, &
67 nb_var, nb_zax, ncdf_ids, regular, slab_ori, slab_sz, xid, yid, zoom
68 use histhori_regular_m, only: histhori_regular
69 USE netcdf, ONLY: nf90_clobber, nf90_create, nf90_def_dim, &
70 nf90_global, nf90_put_att
71
72 CHARACTER (len=*), INTENT (IN):: pfilename
73 REAL, DIMENSION (:), INTENT (IN):: plon_1d
74 REAL, DIMENSION (:), INTENT (IN):: plat_1d
75 INTEGER, INTENT (IN):: par_orix, par_szx, par_oriy, par_szy
76 INTEGER, INTENT (IN):: pitau0
77 REAL, INTENT (IN):: pdate0, pdeltat
78 INTEGER, INTENT (OUT):: pfileid, phoriid
79
80 ! Variables local to the procedure:
81 REAL, DIMENSION (size(plon_1d), size(plat_1d)):: plon, plat
82 INTEGER:: pim, pjm
83 INTEGER:: ncid, iret
84 INTEGER:: lengf, lenga
85 CHARACTER (len=120):: file, tfile
86
87 !---------------------------------------------------------------------
88
89 pim = size(plon_1d)
90 pjm = size(plat_1d)
91
92 plon = spread(plon_1d, 2, pjm)
93 plat = spread(plat_1d, 1, pim)
94
95 nb_files = nb_files + 1
96 pfileid = nb_files
97
98 ! 1.0 Transfering into the common for future use
99
100 itau0(pfileid) = pitau0
101 date0(pfileid) = pdate0
102 deltat(pfileid) = pdeltat
103
104 ! 2.0 Initializes all variables for this file
105
106 IF (nb_files>nb_files_max) THEN
107 CALL histerr(3, 'histbeg', &
108 'Table of files too small. You should increase nb_files_max', &
109 'in M_HISTCOM.f90 in order to accomodate all these files', ' ')
110 END IF
111
112 nb_var(pfileid) = 0
113 nb_tax(pfileid) = 0
114 nb_hax(pfileid) = 0
115 nb_zax(pfileid) = 0
116
117 slab_ori(pfileid, 1:2) = (/ par_orix, par_oriy/)
118 slab_sz(pfileid, 1:2) = (/ par_szx, par_szy/)
119
120 ! 3.0 Opening netcdf file and defining dimensions
121
122 tfile = pfilename
123 lengf = len_trim(tfile)
124 IF (tfile(lengf-2:lengf)/='.nc') THEN
125 file = tfile(1:lengf) // '.nc'
126 ELSE
127 file = tfile(1:lengf)
128 END IF
129
130 ! Add PE number in file name on MPP
131
132 CALL ioipslmpp_file(file)
133
134 ! Keep track of the name of the files opened
135
136 lengf = len_trim(file)
137 lenga = len_trim(assc_file)
138 IF (nb_files==1) THEN
139 assc_file = file(1:lengf)
140 ELSE IF ((lenga+lengf)<500) THEN
141 assc_file = assc_file(1:lenga) // ' ' // file(1:lengf)
142 ELSE IF (((lenga+7)<500) .AND. (index(assc_file(1:lenga), &
143 'et.al.')<1)) THEN
144 assc_file = assc_file(1:lenga) // ' et.al.'
145 ELSE
146 CALL histerr(2, 'histbeg', &
147 'The file names do not fit into the associate_file attribute.', &
148 'Use shorter names if you wish to keep the information.', ' ')
149 END IF
150
151 iret = nf90_create(file, nf90_clobber, ncid)
152 iret = nf90_def_dim(ncid, 'lon', par_szx, xid(nb_files))
153 iret = nf90_def_dim(ncid, 'lat', par_szy, yid(nb_files))
154
155 ! 4.0 Declaring the geographical coordinates and other attributes
156
157 ! 4.3 Global attributes
158
159 iret = nf90_put_att(ncid, nf90_global, 'Conventions', 'GDT 1.3')
160 iret = nf90_put_att(ncid, nf90_global, 'file_name', trim(file))
161 iret = nf90_put_att(ncid, nf90_global, 'production', trim(model_name))
162 lock_modname = .TRUE.
163
164 ! 5.0 Saving some important information on this file in the common
165
166 ncdf_ids(pfileid) = ncid
167 full_size(pfileid, 1:2) = (/ pim, pjm/)
168
169 ! 6.0 storing the geographical coordinates
170
171 IF ((pim/=par_szx) .OR. (pjm/=par_szy)) zoom(pfileid) = .TRUE.
172 regular(pfileid) = .TRUE.
173
174 CALL histhori_regular(pfileid, pim, plon, pjm, plat, ' ', 'Default grid', &
175 phoriid)
176
177 END SUBROUTINE histbeg_totreg
178
179 end MODULE histbeg_totreg_m

  ViewVC Help
Powered by ViewVC 1.1.21