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

Last change on this file since 5272 was 5272, checked in by abarral, 23 hours ago

Turn paramet.h into a module

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