source: LMDZ6/trunk/libf/dyn3d_common/grilles_gcm_netcdf_sub.f90 @ 5404

Last change on this file since 5404 was 5285, checked in by abarral, 7 weeks ago

As discussed internally, remove generic ONLY: ... for new _mod_h modules

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
File size: 8.2 KB
Line 
1!
2! $Id: $
3!
4! This subroutine creates the grilles_gcm.nc file, containing:
5! -> longitudes and latitudes in degrees for dynamical grids u, v and scalaire,
6! and the following variables added for INCA (informative anyway)
7! -> vertical levels "presnivs"
8! -> mask (land/sea), area (grid), phis=surface geopotential height = phis/g
9!
10! The subroutine is called in dynphy_lonlat/phylmd/ce0l.F90.
11
12SUBROUTINE grilles_gcm_netcdf_sub(masque,phis)
13
14  USE comgeom_mod_h
15  USE comconst_mod, ONLY: cpp, kappa, g, omeg, daysec, rad, pi
16  USE comvert_mod, ONLY: presnivs, preff, pa
17  USE netcdf, ONLY: nf90_def_var, nf90_int, nf90_float, nf90_put_var, nf90_clobber, nf90_64bit_offset, nf90_def_dim, &
18          nf90_put_att, nf90_enddef, nf90_create
19
20  USE dimensions_mod, ONLY: iim, jjm, llm, ndm
21USE paramet_mod_h
22IMPLICIT NONE
23
24
25
26
27!========================
28  REAL,DIMENSION(iip1,jjp1),INTENT(IN)  :: masque ! masque terre/mer
29  REAL,DIMENSION(iip1,jjp1),INTENT(IN)  :: phis   ! geopotentiel au sol
30
31  INTEGER status,i,j
32 
33  ! Attributs netcdf output
34  INTEGER ncid_out,rcode_out
35 
36  INTEGER out_lonuid,out_lonvid,out_latuid,out_latvid
37  INTEGER out_uid,out_vid,out_tempid
38  INTEGER out_lonudim,out_lonvdim
39  INTEGER out_latudim,out_latvdim,out_dim(2)
40  INTEGER out_levdim
41  !
42  INTEGER :: presnivs_id
43  INTEGER :: mask_id,area_id,phis_id
44  !
45  INTEGER start(2),COUNT(2)
46
47  ! Variables
48  REAL rlatudeg(jjp1),rlatvdeg(jjm),rlev(llm)
49  REAL rlonudeg(iip1),rlonvdeg(iip1)
50  REAL uwnd(iip1,jjp1),vwnd(iip1,jjm),temp(iip1,jjp1)
51  !
52  INTEGER masque_int(iip1,jjp1)
53  REAL :: phis_loc(iip1,jjp1)
54 
55  !========================
56  ! CALCULATION of latu, latv, lonu, lonv in deg.
57  ! ---------------------------------------------------
58  rad = 6400000
59  omeg = 7.272205e-05
60  g = 9.8
61  kappa = 0.285716
62  daysec = 86400
63  cpp = 1004.70885
64
65  preff = 101325.
66  pa= 50000.
67
68  CALL conf_gcm( 99, .TRUE. )
69  CALL iniconst
70  CALL inigeom
71
72  DO j=1,jjp1
73     rlatudeg(j)=rlatu(j)*180./pi
74  ENDDO
75 
76  DO j=1,jjm
77     rlatvdeg(j)=rlatv(j)*180./pi
78  ENDDO
79
80  DO i=1,iip1
81     rlonudeg(i)=rlonu(i)*180./pi + 360.
82     rlonvdeg(i)=rlonv(i)*180./pi + 360.
83  ENDDO
84 
85  ! CALCULATION of "false" variables on u, v, s grids
86  ! ---------------------------------------------------
87   DO i=1,iip1
88     DO j=1,jjp1
89        uwnd(i,j)=MOD(i,2)+MOD(j,2)
90        temp(i,j)=MOD(i,2)+MOD(j,2)
91     ENDDO
92     DO j=1,jjm
93        vwnd(i,j)=MOD(i,2)+MOD(j,2)
94     END DO
95  ENDDO 
96
97  ! CALCULATION of local vars for presnivs, mask, sfc. geopot. height
98  ! ---------------------------------------------------
99  rlev(:) = presnivs(:)
100  phis_loc(:,:) = phis(:,:)/g
101  masque_int(:,:) = nINT(masque(:,:))
102
103
104  ! OPEN output netcdf file
105  !-------------------------
106  status=nf90_create('grilles_gcm.nc',IOR(nf90_clobber,nf90_64bit_offset),ncid_out)
107  CALL handle_err(status)
108 
109  ! DEFINE output dimensions
110  !-------------------------
111  status=nf90_def_dim(ncid_out,'lonu',iim+1,out_lonudim)
112  CALL handle_err(status)
113  status=nf90_def_dim(ncid_out,'lonv',iim+1,out_lonvdim)
114  CALL handle_err(status)
115  status=nf90_def_dim(ncid_out,'latu',jjm+1,out_latudim)
116  CALL handle_err(status)
117  status=nf90_def_dim(ncid_out,'latv',jjm,out_latvdim)
118  CALL handle_err(status)
119  !
120  status=nf90_def_dim(ncid_out,'lev',llm,out_levdim)
121  CALL handle_err(status)
122 
123  ! DEFINE output variables
124  !-------------------------
125  !   Longitudes on "u" dynamical grid
126  status=NF90_DEF_VAR(ncid_out,'lonu',NF90_FLOAT,out_lonudim, out_lonuid)
127  CALL handle_err(status)
128  status=nf90_put_att(ncid_out,out_lonuid,'units','degrees_east')
129  status=nf90_put_att(ncid_out,out_lonuid,'long_name','Longitude on u grid')
130  !   Longitudes on "v" dynamical grid
131  status=NF90_DEF_VAR(ncid_out,'lonv',NF90_FLOAT,out_lonvdim, out_lonvid)
132  CALL handle_err(status)
133  status=nf90_put_att(ncid_out,out_lonvid,'units','degrees_east')
134  status=nf90_put_att(ncid_out,out_lonvid,'long_name','Longitude on v grid')
135  !   Latitudes on "u" dynamical grid
136  status=NF90_DEF_VAR(ncid_out,'latu',NF90_FLOAT,out_latudim, out_latuid)
137  CALL handle_err(status)
138  status=nf90_put_att(ncid_out,out_latuid,'units','degrees_north')
139  status=nf90_put_att(ncid_out,out_latuid,'long_name','Latitude on u grid')
140  !  Latitudes on "v" dynamical grid
141  status=NF90_DEF_VAR(ncid_out,'latv',NF90_FLOAT,out_latvdim, out_latvid)
142  CALL handle_err(status)
143  status=nf90_put_att(ncid_out,out_latvid,'units','degrees_north')
144  status=nf90_put_att(ncid_out,out_latvid,'long_name','Latitude on v grid')
145  !  "u" lat/lon dynamical grid
146  out_dim(1)=out_lonudim
147  out_dim(2)=out_latudim
148  status=NF90_DEF_VAR(ncid_out,'grille_u',NF90_FLOAT,out_dim, out_uid)
149  CALL handle_err(status)
150  status=nf90_put_att(ncid_out,out_uid,'units','m/s')
151  status=nf90_put_att(ncid_out,out_uid,'long_name','u-wind dynamical grid')
152  !  "v" lat/lon dynamical grid
153  out_dim(1)=out_lonvdim
154  out_dim(2)=out_latvdim
155  status=NF90_DEF_VAR(ncid_out,'grille_v',NF90_FLOAT,out_dim, out_vid)
156  CALL handle_err(status)
157  status=nf90_put_att(ncid_out,out_vid,'units','m/s')
158  status=nf90_put_att(ncid_out,out_vid,'long_name','v-wind dynamical grid')
159  !  "s" (scalar) lat/lon dynamical grid
160  out_dim(1)=out_lonvdim
161  out_dim(2)=out_latudim
162  status=NF90_DEF_VAR(ncid_out,'grille_s',NF90_FLOAT,out_dim, out_tempid)
163  CALL handle_err(status)
164  status=nf90_put_att(ncid_out,out_tempid,'units','Kelvin')
165  status=nf90_put_att(ncid_out,out_tempid,'long_name','scalar dynamical grid')
166  !
167  ! for INCA :
168  ! vertical levels "presnivs"
169  status=NF90_DEF_VAR(ncid_out,'presnivs',NF90_FLOAT,out_levdim, presnivs_id)
170  CALL handle_err(status)
171  status=nf90_put_att(ncid_out,presnivs_id,'units','Pa')
172  status=nf90_put_att(ncid_out,presnivs_id,'long_name','Vertical levels')
173  ! surface geopotential height: named "phis" as the sfc geopotential, is actually phis/g
174  out_dim(1)=out_lonvdim
175  out_dim(2)=out_latudim
176  status = nf90_def_var(ncid_out,'phis',NF90_FLOAT,out_dim,phis_id)
177  CALL handle_err(status)
178  status=nf90_put_att(ncid_out,phis_id,'units','m')
179  status=nf90_put_att(ncid_out,phis_id,'long_name','surface geopotential height')
180  ! gridcell area
181  status = nf90_def_var(ncid_out,'aire',NF90_FLOAT,out_dim,area_id)
182  CALL handle_err(status)
183  status=nf90_put_att(ncid_out,area_id,'units','m2')
184  status=nf90_put_att(ncid_out,area_id,'long_name','gridcell area')
185  ! land-sea mask (nearest integer approx)
186  status = nf90_def_var(ncid_out,'mask',NF90_INT,out_dim,mask_id)
187  CALL handle_err(status)
188  status=nf90_put_att(ncid_out,mask_id,'long_name','land-sea mask (nINT approx)')
189 
190  ! END the 'define' mode in netCDF file
191  status=nf90_enddef(ncid_out)
192  CALL handle_err(status)
193 
194  ! WRITE the variables
195  !-------------------------
196  ! 1D : lonu, lonv,latu,latv ; INCA : presnivs
197  status=NF90_PUT_VAR(ncid_out,out_lonuid,rlonudeg,[1],[iip1])
198  CALL handle_err(status)
199  status=NF90_PUT_VAR(ncid_out,out_lonvid,rlonvdeg,[1],[iip1])
200  CALL handle_err(status)
201  status=NF90_PUT_VAR(ncid_out,out_latuid,rlatudeg,[1],[jjp1])
202  CALL handle_err(status)
203  status=NF90_PUT_VAR(ncid_out,out_latvid,rlatvdeg,[1],[jjm])
204  CALL handle_err(status)
205  status=NF90_PUT_VAR(ncid_out,presnivs_id,rlev,[1],[llm])
206  CALL handle_err(status)
207
208  ! 2D : grille_u,grille_v,grille_s ; INCA: phis,aire,mask
209  start(:)=1
210  COUNT(1)=iip1
211 
212  COUNT(2)=jjp1  ! for "u" and "s" grids
213  status=NF90_PUT_VAR(ncid_out,out_uid,uwnd,start, count)
214  CALL handle_err(status)
215  COUNT(2)=jjm  ! for "v" grid
216  status=NF90_PUT_VAR(ncid_out,out_vid,vwnd,start, count)
217  CALL handle_err(status)
218  COUNT(2)=jjp1  ! as "s" grid, for all the following vars
219  status=NF90_PUT_VAR(ncid_out,out_tempid,temp,start, count)
220  CALL handle_err(status)
221  status = nf90_put_var(ncid_out, phis_id, phis_loc,start,count)
222  CALL handle_err(status)
223  status = nf90_put_var(ncid_out, area_id, aire,start,count)
224  CALL handle_err(status)
225  status = nf90_put_var(ncid_out, mask_id,masque_int,start,count)
226  CALL handle_err(status)
227 
228  ! CLOSE netcdf file
229  CALL ncclos(ncid_out,rcode_out)
230  write(*,*) "END grilles_gcm_netcdf_sub OK"
231
232END SUBROUTINE grilles_gcm_netcdf_sub
233
234
235SUBROUTINE handle_err(status)
236  USE netcdf, ONLY: nf90_noerr, nf90_strerror
237  IMPLICIT NONE
238
239  INTEGER status
240  IF (status.NE.nf90_noerr) THEN
241     PRINT *,nf90_strerror(status)
242     CALL abort_gcm('grilles_gcm_netcdf','netcdf error',1)
243  ENDIF
244END SUBROUTINE handle_err
245
Note: See TracBrowser for help on using the repository browser.