/[lmdze]/trunk/libf/IOIPSL/Histcom/histhori_regular.f90
ViewVC logotype

Annotation of /trunk/libf/IOIPSL/Histcom/histhori_regular.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 61 - (hide annotations)
Fri Apr 20 14:58:43 2012 UTC (12 years, 1 month ago) by guez
File size: 3979 byte(s)
No more included file in LMDZE, not even "netcdf.inc".

Created a variable containing the list of common source files in
GNUmakefile. So we now also see clearly files that are specific to
each program.

Split module "histcom". Assembled resulting files in directory
"Histcom".

Removed aliasing in calls to "laplacien".

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     USE histcom_var, ONLY: full_size, hax_name, nb_hax, ncdf_ids, &
36     slab_ori, slab_sz, xid, yid
37     USE netcdf, ONLY: nf90_def_var, nf90_enddef, nf90_float, &
38     nf90_put_att, nf90_put_var, nf90_redef
39    
40     INTEGER, INTENT (IN):: pfileid, pim, pjm
41     REAL, INTENT (IN), DIMENSION (pim, pjm):: plon, plat
42     CHARACTER (len=*), INTENT (IN):: phname, phtitle
43     INTEGER, INTENT (OUT):: phid
44    
45     CHARACTER (len=25):: lon_name, lat_name
46     CHARACTER (len=80):: tmp_title, tmp_name
47     INTEGER:: ndim
48     INTEGER, DIMENSION (2):: dims(2)
49     INTEGER:: nlonid, nlatid
50     INTEGER:: orix, oriy, par_szx, par_szy
51     INTEGER:: iret, ncid
52    
53     !---------------------------------------------------------------------
54    
55     ! 1.0 Check that all fits in the buffers
56    
57     IF ((pim/=full_size(pfileid, 1)) .OR. (pjm/=full_size(pfileid, 2))) THEN
58     CALL histerr(3, 'histhori', &
59     'The new horizontal grid does not have the same size', &
60     'as the one provided to histbeg. This is not yet ', &
61     'possible in the hist package.')
62     END IF
63    
64     ! 1.1 Create all the variables needed
65    
66     ncid = ncdf_ids(pfileid)
67    
68     ndim = 2
69     dims(1:2) = (/ xid(pfileid), yid(pfileid) /)
70    
71     tmp_name = phname
72     IF (nb_hax(pfileid)==0) THEN
73     lon_name = 'lon'
74     lat_name = 'lat'
75     ELSE
76     lon_name = 'lon_' // trim(tmp_name)
77     lat_name = 'lat_' // trim(tmp_name)
78     END IF
79    
80     ! 1.2 Save the informations
81    
82     phid = nb_hax(pfileid) + 1
83     nb_hax(pfileid) = phid
84    
85     hax_name(pfileid, phid, 1:2) = (/ lon_name, lat_name/)
86     tmp_title = phtitle
87    
88     ! 2.0 Longitude
89    
90     ndim = 1
91     dims(1:1) = (/ xid(pfileid) /)
92    
93     iret = nf90_def_var(ncid, lon_name, nf90_float, dims(1:ndim), nlonid)
94     iret = nf90_put_att(ncid, nlonid, 'units', 'degrees_east')
95     iret = nf90_put_att(ncid, nlonid, 'valid_min', real(minval(plon)))
96     iret = nf90_put_att(ncid, nlonid, 'valid_max', real(maxval(plon)))
97     iret = nf90_put_att(ncid, nlonid, 'long_name', 'Longitude')
98     iret = nf90_put_att(ncid, nlonid, 'nav_model', trim(tmp_title))
99    
100     ! 3.0 Latitude
101    
102     ndim = 1
103     dims(1:1) = (/ yid(pfileid) /)
104    
105     iret = nf90_def_var(ncid, lat_name, nf90_float, dims(1:ndim), nlatid)
106     iret = nf90_put_att(ncid, nlatid, 'units', 'degrees_north')
107     iret = nf90_put_att(ncid, nlatid, 'valid_min', real(minval(plat)))
108     iret = nf90_put_att(ncid, nlatid, 'valid_max', real(maxval(plat)))
109     iret = nf90_put_att(ncid, nlatid, 'long_name', 'Latitude')
110     iret = nf90_put_att(ncid, nlatid, 'nav_model', trim(tmp_title))
111    
112     iret = nf90_enddef(ncid)
113    
114     ! 4.0 storing the geographical coordinates
115    
116     orix = slab_ori(pfileid, 1)
117     oriy = slab_ori(pfileid, 2)
118     par_szx = slab_sz(pfileid, 1)
119     par_szy = slab_sz(pfileid, 2)
120    
121     ! Transfer the longitude
122    
123     iret = nf90_put_var(ncid, nlonid, plon(1:par_szx, 1))
124    
125     ! Transfer the latitude
126    
127     iret = nf90_put_var(ncid, nlatid, plat(1, 1:par_szy))
128    
129     iret = nf90_redef(ncid)
130    
131     END SUBROUTINE histhori_regular
132    
133     end module histhori_regular_m

  ViewVC Help
Powered by ViewVC 1.1.21