/[lmdze]/trunk/dyn3d/grilles_gcm_netcdf_sub.f
ViewVC logotype

Annotation of /trunk/dyn3d/grilles_gcm_netcdf_sub.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 294 - (hide annotations)
Thu Jul 26 13:11:22 2018 UTC (5 years, 10 months ago) by guez
File size: 3576 byte(s)
Remove variables grille_u and grille_v of file grilles_gcm.nc, which
were not given any value.

In procedure clqh, we only need to define gamt if iflag_pbl == 1.

1 guez 107 module grilles_gcm_netcdf_sub_m
2    
3     IMPLICIT NONE
4    
5     contains
6    
7     SUBROUTINE grilles_gcm_netcdf_sub(phis)
8    
9 guez 118 ! This subroutine creates the file "grilles_gcm.nc" containg
10 guez 107 ! longitudes and latitudes in degrees for grids u and v.
11    
12     USE comconst, ONLY: g
13 guez 139 USE comgeom, ONLY: aire_2d
14 guez 265 USE dimensions, ONLY: iim, jjm
15 guez 139 use dynetat0_m, only: rlatu, rlatv, rlonu, rlonv
16 guez 107 USE netcdf, ONLY: nf90_clobber, nf90_float, nf90_int
17     USE netcdf95, ONLY: nf95_close, nf95_create, nf95_def_dim, nf95_def_var, &
18 guez 178 nf95_enddef, nf95_put_att, nf95_put_var
19 guez 107 USE nr_util, ONLY: pi
20     USE start_init_orog_m, ONLY: mask
21    
22 guez 118 REAL, INTENT(IN):: phis(:, :) ! (iim + 1, jjm + 1)
23     ! surface geopotential, in m2 s-2
24 guez 107
25     ! Local:
26    
27 guez 140 REAL grille_s(iim + 1, jjm + 1)
28 guez 107 INTEGER i, j
29    
30     ! For Netcdf:
31     INTEGER ncid
32 guez 118 INTEGER dimid_lonu, dimid_lonv, dimid_latu, dimid_latv
33 guez 107 INTEGER varid_lonu, varid_lonv, varid_latu, varid_latv
34     INTEGER varid
35 guez 118 INTEGER varid_orog, varid_area, varid_mask
36 guez 107
37     !-----------------------------------------------------------------
38    
39     print *, "Call sequence information: grilles_gcm_netcdf_sub"
40    
41     call NF95_CREATE('grilles_gcm.nc', NF90_CLOBBER, ncid)
42 guez 118
43     call NF95_DEF_DIM(ncid, 'lonu', iim + 1, dimid_lonu)
44     call NF95_DEF_DIM(ncid, 'lonv', iim + 1, dimid_lonv)
45     call NF95_DEF_DIM(ncid, 'latu', jjm + 1, dimid_latu)
46 guez 107 call NF95_DEF_DIM(ncid, 'latv', jjm, dimid_latv)
47    
48     call NF95_DEF_VAR(ncid, 'lonu', NF90_FLOAT, dimid_lonu, varid_lonu)
49     call NF95_PUT_ATT(ncid, varid_lonu, 'units', 'degrees_east')
50     call NF95_PUT_ATT(ncid, varid_lonu, 'long_name', 'Longitude en u')
51    
52     call NF95_DEF_VAR(ncid, 'lonv', NF90_FLOAT, dimid_lonv, varid_lonv)
53     call NF95_PUT_ATT(ncid, varid_lonv, 'units', 'degrees_east')
54     call NF95_PUT_ATT(ncid, varid_lonv, 'long_name', 'Longitude en v')
55    
56     call NF95_DEF_VAR(ncid, 'latu', NF90_FLOAT, dimid_latu, varid_latu)
57     call NF95_PUT_ATT(ncid, varid_latu, 'units', 'degrees_north')
58     call NF95_PUT_ATT(ncid, varid_latu, 'long_name', 'Latitude en u')
59    
60     call NF95_DEF_VAR(ncid, 'latv', NF90_FLOAT, dimid_latv, varid_latv)
61     call NF95_PUT_ATT(ncid, varid_latv, 'units', 'degrees_north')
62     call NF95_PUT_ATT(ncid, varid_latv, 'long_name', 'Latitude en v')
63    
64     call NF95_DEF_VAR(ncid, 'grille_s', NF90_FLOAT, (/dimid_lonv, &
65     dimid_latu/), varid)
66 guez 118 call NF95_PUT_ATT(ncid, varid, 'long_name', 'Grille aux points scalaires')
67 guez 107
68 guez 118 call NF95_DEF_VAR(ncid, 'orog', NF90_FLOAT, (/dimid_lonv, dimid_latu/), &
69     varid_orog)
70     call NF95_PUT_ATT(ncid, varid_orog, 'units', 'm')
71     call NF95_PUT_ATT(ncid, varid_orog, 'standard_name', 'surface_altitude')
72 guez 107
73 guez 118 call NF95_DEF_VAR(ncid, 'aire', NF90_FLOAT, (/dimid_lonv, dimid_latu/), &
74     varid_area)
75     call NF95_PUT_ATT(ncid, varid_area, 'units', 'm2')
76 guez 107
77     call NF95_DEF_VAR(ncid, 'mask', NF90_INT, (/dimid_lonv, dimid_latu/), &
78     varid_mask)
79    
80     call NF95_ENDDEF(ncid)
81    
82 guez 118 call NF95_PUT_VAR(ncid, varid_lonu, rlonu * 180. / pi)
83     call NF95_PUT_VAR(ncid, varid_lonv, rlonv * 180. / pi)
84     call NF95_PUT_VAR(ncid, varid_latu, rlatu * 180. / pi)
85     call NF95_PUT_VAR(ncid, varid_latv, rlatv * 180. / pi)
86 guez 107
87 guez 118 forall (i = 1:iim + 1, j = 1:jjm + 1) grille_s(i, j) = MOD(i, 2) + MOD(j, 2)
88     call NF95_PUT_VAR(ncid, varid, grille_s)
89 guez 107
90 guez 118 call NF95_PUT_VAR(ncid, varid_orog, phis / g)
91 guez 107 call NF95_PUT_VAR(ncid, varid_area, aire_2d)
92     call NF95_PUT_VAR(ncid, varid_mask, nINT(mask))
93    
94     CALL nf95_close(ncid)
95    
96     END SUBROUTINE grilles_gcm_netcdf_sub
97    
98     end module grilles_gcm_netcdf_sub_m

  ViewVC Help
Powered by ViewVC 1.1.21