/[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 134 - (hide annotations)
Wed Apr 29 15:47:56 2015 UTC (9 years, 1 month ago) by guez
Original Path: trunk/Sources/dyn3d/grilles_gcm_netcdf_sub.f
File size: 3897 byte(s)
Sources inside, compilation outside.
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     USE comgeom, ONLY: aire_2d, rlatu, rlatv, rlonu, rlonv
14 guez 118 USE dimens_m, ONLY: iim, jjm
15 guez 107 USE netcdf, ONLY: nf90_clobber, nf90_float, nf90_int
16     USE netcdf95, ONLY: nf95_close, nf95_create, nf95_def_dim, nf95_def_var, &
17     nf95_enddef, nf95_put_att, nf95_put_var, nf95_redef
18     USE nr_util, ONLY: pi
19     USE start_init_orog_m, ONLY: mask
20    
21 guez 118 REAL, INTENT(IN):: phis(:, :) ! (iim + 1, jjm + 1)
22     ! surface geopotential, in m2 s-2
23 guez 107
24     ! Local:
25    
26 guez 118 REAL grille_s(iim+1, jjm+1)
27 guez 107 INTEGER i, j
28    
29     ! For Netcdf:
30     INTEGER ncid
31 guez 118 INTEGER dimid_lonu, dimid_lonv, dimid_latu, dimid_latv
32 guez 107 INTEGER varid_lonu, varid_lonv, varid_latu, varid_latv
33     INTEGER varid
34 guez 118 INTEGER varid_orog, varid_area, varid_mask
35 guez 107
36     !-----------------------------------------------------------------
37    
38     print *, "Call sequence information: grilles_gcm_netcdf_sub"
39    
40     call NF95_CREATE('grilles_gcm.nc', NF90_CLOBBER, ncid)
41 guez 118
42     call NF95_DEF_DIM(ncid, 'lonu', iim + 1, dimid_lonu)
43     call NF95_DEF_DIM(ncid, 'lonv', iim + 1, dimid_lonv)
44     call NF95_DEF_DIM(ncid, 'latu', jjm + 1, dimid_latu)
45 guez 107 call NF95_DEF_DIM(ncid, 'latv', jjm, dimid_latv)
46    
47     call NF95_DEF_VAR(ncid, 'lonu', NF90_FLOAT, dimid_lonu, varid_lonu)
48     call NF95_PUT_ATT(ncid, varid_lonu, 'units', 'degrees_east')
49     call NF95_PUT_ATT(ncid, varid_lonu, 'long_name', 'Longitude en u')
50    
51     call NF95_DEF_VAR(ncid, 'lonv', NF90_FLOAT, dimid_lonv, varid_lonv)
52     call NF95_PUT_ATT(ncid, varid_lonv, 'units', 'degrees_east')
53     call NF95_PUT_ATT(ncid, varid_lonv, 'long_name', 'Longitude en v')
54    
55     call NF95_DEF_VAR(ncid, 'latu', NF90_FLOAT, dimid_latu, varid_latu)
56     call NF95_PUT_ATT(ncid, varid_latu, 'units', 'degrees_north')
57     call NF95_PUT_ATT(ncid, varid_latu, 'long_name', 'Latitude en u')
58    
59     call NF95_DEF_VAR(ncid, 'latv', NF90_FLOAT, dimid_latv, varid_latv)
60     call NF95_PUT_ATT(ncid, varid_latv, 'units', 'degrees_north')
61     call NF95_PUT_ATT(ncid, varid_latv, 'long_name', 'Latitude en v')
62    
63     call NF95_DEF_VAR(ncid, 'grille_u', NF90_FLOAT, (/dimid_lonu, &
64     dimid_latu/), varid)
65     call NF95_PUT_ATT(ncid, varid, 'long_name', 'Grille aux points u')
66    
67     call NF95_DEF_VAR(ncid, 'grille_v', NF90_FLOAT, (/dimid_lonv, &
68     dimid_latv/), varid)
69     call NF95_PUT_ATT(ncid, varid, 'long_name', 'Grille aux points v')
70    
71     call NF95_DEF_VAR(ncid, 'grille_s', NF90_FLOAT, (/dimid_lonv, &
72     dimid_latu/), varid)
73 guez 118 call NF95_PUT_ATT(ncid, varid, 'long_name', 'Grille aux points scalaires')
74 guez 107
75 guez 118 call NF95_DEF_VAR(ncid, 'orog', NF90_FLOAT, (/dimid_lonv, dimid_latu/), &
76     varid_orog)
77     call NF95_PUT_ATT(ncid, varid_orog, 'units', 'm')
78     call NF95_PUT_ATT(ncid, varid_orog, 'standard_name', 'surface_altitude')
79 guez 107
80 guez 118 call NF95_DEF_VAR(ncid, 'aire', NF90_FLOAT, (/dimid_lonv, dimid_latu/), &
81     varid_area)
82     call NF95_PUT_ATT(ncid, varid_area, 'units', 'm2')
83 guez 107
84     call NF95_DEF_VAR(ncid, 'mask', NF90_INT, (/dimid_lonv, dimid_latu/), &
85     varid_mask)
86    
87     call NF95_ENDDEF(ncid)
88    
89 guez 118 call NF95_PUT_VAR(ncid, varid_lonu, rlonu * 180. / pi)
90     call NF95_PUT_VAR(ncid, varid_lonv, rlonv * 180. / pi)
91     call NF95_PUT_VAR(ncid, varid_latu, rlatu * 180. / pi)
92     call NF95_PUT_VAR(ncid, varid_latv, rlatv * 180. / pi)
93 guez 107
94 guez 118 forall (i = 1:iim + 1, j = 1:jjm + 1) grille_s(i, j) = MOD(i, 2) + MOD(j, 2)
95     call NF95_PUT_VAR(ncid, varid, grille_s)
96 guez 107
97 guez 118 call NF95_PUT_VAR(ncid, varid_orog, phis / g)
98 guez 107 call NF95_PUT_VAR(ncid, varid_area, aire_2d)
99     call NF95_PUT_VAR(ncid, varid_mask, nINT(mask))
100    
101     CALL nf95_close(ncid)
102    
103     END SUBROUTINE grilles_gcm_netcdf_sub
104    
105     end module grilles_gcm_netcdf_sub_m

  ViewVC Help
Powered by ViewVC 1.1.21