/[lmdze]/trunk/Sources/IOIPSL/Histcom/histhori_regular.f
ViewVC logotype

Annotation of /trunk/Sources/IOIPSL/Histcom/histhori_regular.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 178 - (hide annotations)
Fri Mar 11 18:47:26 2016 UTC (8 years, 2 months ago) by guez
File size: 3859 byte(s)
Moved variables date0, deltat, datasz_max, ncvar_ids, point, buff_pos,
buffer, regular from module histcom_var to modules where they are
defined.

Removed procedure ioipslmpp, useless for a sequential program.

Added argument datasz_max to histwrite_real (to avoid circular
dependency with histwrite).

Removed useless variables and computations everywhere.

Changed real litteral constants from default kind to double precision
in lwb, lwu, lwvn, sw1s, swtt, swtt1, swu.

Removed unused arguments: paer of sw, sw1s, sw2s, swclr; pcldsw of
sw1s, sw2s; pdsig, prayl of swr; co2_ppm of clmain, clqh; tsol of
transp_lay; nsrf of screenp; kcrit and kknu of gwstress; pstd of
orosetup.

Added output of relative humidity.

1 guez 61 module histhori_regular_m
2    
3     implicit none
4    
5     contains
6    
7     SUBROUTINE histhori_regular(pfileid, pim, plon, pjm, plat, phname, phtitle, &
8     phid)
9    
10     ! This subroutine is made to declare a new horizontale grid.
11     ! It has to have the same number of points as
12     ! the original Thus in this we routine we will only
13     ! add two variable (longitude and latitude).
14     ! Any variable in the file can thus point to this pair
15     ! through an attribute. This routine is very usefull
16     ! to allow staggered grids.
17    
18     ! INPUT
19    
20     ! pfileid: The id of the file to which the grid should be added
21     ! pim: Size in the longitude direction
22     ! plon: The longitudes
23     ! pjm: Size in the latitude direction
24     ! plat: The latitudes
25     ! phname: The name of grid
26     ! phtitle: The title of the grid
27    
28     ! OUTPUT
29    
30     ! phid: Id of the created grid
31    
32     ! We assume that the grid is rectilinear.
33    
34     USE errioipsl, ONLY: histerr
35 guez 178 USE histcom_var, ONLY: full_size, hax_name, nb_hax, ncdf_ids, slab_sz, &
36     xid, yid
37     USE netcdf, ONLY: nf90_float
38     use netcdf95, only: nf95_redef, nf95_def_var, nf95_enddef, nf95_put_att, &
39     nf95_put_var
40 guez 61
41     INTEGER, INTENT (IN):: pfileid, pim, pjm
42     REAL, INTENT (IN), DIMENSION (pim, pjm):: plon, plat
43     CHARACTER (len=*), INTENT (IN):: phname, phtitle
44     INTEGER, INTENT (OUT):: phid
45    
46     CHARACTER (len=25):: lon_name, lat_name
47     CHARACTER (len=80):: tmp_title, tmp_name
48     INTEGER:: ndim
49 guez 178 INTEGER dims(2)
50 guez 61 INTEGER:: nlonid, nlatid
51 guez 178 INTEGER:: par_szx, par_szy
52     INTEGER ncid
53 guez 61
54     !---------------------------------------------------------------------
55    
56     ! 1.0 Check that all fits in the buffers
57    
58     IF ((pim/=full_size(pfileid, 1)) .OR. (pjm/=full_size(pfileid, 2))) THEN
59     CALL histerr(3, 'histhori', &
60     'The new horizontal grid does not have the same size', &
61     'as the one provided to histbeg. This is not yet ', &
62     'possible in the hist package.')
63     END IF
64    
65     ! 1.1 Create all the variables needed
66    
67     ncid = ncdf_ids(pfileid)
68    
69     ndim = 2
70     dims(1:2) = (/ xid(pfileid), yid(pfileid) /)
71    
72     tmp_name = phname
73     IF (nb_hax(pfileid)==0) THEN
74     lon_name = 'lon'
75     lat_name = 'lat'
76     ELSE
77     lon_name = 'lon_' // trim(tmp_name)
78     lat_name = 'lat_' // trim(tmp_name)
79     END IF
80    
81     ! 1.2 Save the informations
82    
83     phid = nb_hax(pfileid) + 1
84     nb_hax(pfileid) = phid
85    
86     hax_name(pfileid, phid, 1:2) = (/ lon_name, lat_name/)
87     tmp_title = phtitle
88    
89     ! 2.0 Longitude
90    
91     ndim = 1
92     dims(1:1) = (/ xid(pfileid) /)
93    
94 guez 178 call nf95_def_var(ncid, lon_name, nf90_float, dims(1:ndim), nlonid)
95     call nf95_put_att(ncid, nlonid, 'units', 'degrees_east')
96     call nf95_put_att(ncid, nlonid, 'valid_min', real(minval(plon)))
97     call nf95_put_att(ncid, nlonid, 'valid_max', real(maxval(plon)))
98     call nf95_put_att(ncid, nlonid, 'long_name', 'Longitude')
99     call nf95_put_att(ncid, nlonid, 'nav_model', trim(tmp_title))
100 guez 61
101     ! 3.0 Latitude
102    
103     ndim = 1
104     dims(1:1) = (/ yid(pfileid) /)
105    
106 guez 178 call nf95_def_var(ncid, lat_name, nf90_float, dims(1:ndim), nlatid)
107     call nf95_put_att(ncid, nlatid, 'units', 'degrees_north')
108     call nf95_put_att(ncid, nlatid, 'valid_min', real(minval(plat)))
109     call nf95_put_att(ncid, nlatid, 'valid_max', real(maxval(plat)))
110     call nf95_put_att(ncid, nlatid, 'long_name', 'Latitude')
111     call nf95_put_att(ncid, nlatid, 'nav_model', trim(tmp_title))
112 guez 61
113 guez 178 call nf95_enddef(ncid)
114 guez 61
115     ! 4.0 storing the geographical coordinates
116    
117     par_szx = slab_sz(pfileid, 1)
118     par_szy = slab_sz(pfileid, 2)
119    
120     ! Transfer the longitude
121    
122 guez 178 call nf95_put_var(ncid, nlonid, plon(1:par_szx, 1))
123 guez 61
124     ! Transfer the latitude
125    
126 guez 178 call nf95_put_var(ncid, nlatid, plat(1, 1:par_szy))
127 guez 61
128 guez 168 call nf95_redef(ncid)
129 guez 61
130     END SUBROUTINE histhori_regular
131    
132     end module histhori_regular_m

  ViewVC Help
Powered by ViewVC 1.1.21