1 | subroutine writediagsoil(ngrid,name,title,units,dimpx,px) |
---|
2 | |
---|
3 | ! Write variable 'name' to NetCDF file 'diagsoil.nc'. |
---|
4 | ! The variable may be 3D (lon,lat,depth) subterranean field, |
---|
5 | ! a 2D (lon,lat) surface field, or a simple scalar (0D variable). |
---|
6 | ! |
---|
7 | ! Calls to 'writediagsoil' can originate from anywhere in the program; |
---|
8 | ! An initialisation of variable 'name' is done if it is the first time |
---|
9 | ! that this routine is called with given 'name'; otherwise data is appended |
---|
10 | ! (yielding the sought time series of the variable) |
---|
11 | |
---|
12 | ! Modifs: Aug.2010 Ehouarn: enforce outputs to be real*4 |
---|
13 | |
---|
14 | implicit none |
---|
15 | |
---|
16 | #include"dimensions.h" |
---|
17 | #include"dimphys.h" |
---|
18 | #include"paramet.h" |
---|
19 | #include"control.h" |
---|
20 | #include"comsoil.h" |
---|
21 | #include"netcdf.inc" |
---|
22 | |
---|
23 | ! Arguments: |
---|
24 | integer,intent(in) :: ngrid ! number of (horizontal) points of physics grid |
---|
25 | ! i.e. ngrid = 2+(jjm-1)*iim - 1/jjm |
---|
26 | character(len=*),intent(in) :: name ! 'name' of the variable |
---|
27 | character(len=*),intent(in) :: title ! 'long_name' attribute of the variable |
---|
28 | character(len=*),intent(in) :: units ! 'units' attribute of the variable |
---|
29 | integer,intent(in) :: dimpx ! dimension of the variable (3,2 or 0) |
---|
30 | real,dimension(ngrid,nsoilmx),intent(in) :: px ! variable |
---|
31 | ! Note: nsoilmx is a common parameter set in 'dimphys.h' |
---|
32 | |
---|
33 | ! Local variables: |
---|
34 | real*4,dimension(iip1,jjp1,nsoilmx) :: data3 ! to store 3D data |
---|
35 | ! Note iip1,jjp1 known from paramet.h; nsoilmx known from dimphys.h |
---|
36 | real*4,dimension(iip1,jjp1) :: data2 ! to store 2D data |
---|
37 | real*4 :: data0 ! to store 0D data |
---|
38 | integer :: i,j,l ! for loops |
---|
39 | integer :: ig0 |
---|
40 | |
---|
41 | real*4,save :: date ! time counter (in elapsed days) |
---|
42 | integer,save :: isample ! sample rate at which data is to be written to output |
---|
43 | integer,save :: ntime=0 ! counter to internally store time steps |
---|
44 | character(len=20),save :: firstname="1234567890" |
---|
45 | integer,save :: zitau=0 |
---|
46 | |
---|
47 | character(len=30) :: filename="diagsoil.nc" |
---|
48 | |
---|
49 | ! NetCDF stuff: |
---|
50 | integer :: nid ! NetCDF output file ID |
---|
51 | integer :: varid ! NetCDF ID of a variable |
---|
52 | integer :: ierr ! NetCDF routines return code |
---|
53 | integer,dimension(4) :: id ! NetCDF IDs of the dimensions of the variable |
---|
54 | integer,dimension(4) :: edges,corners |
---|
55 | |
---|
56 | ! 1. Initialization step |
---|
57 | if (firstname.eq."1234567890") then |
---|
58 | ! Store 'name' as 'firstname' |
---|
59 | firstname=name |
---|
60 | ! From now on, if 'name'.eq.'firstname', then it is a new time cycle |
---|
61 | |
---|
62 | ! just to be sure, check that firstnom is large enough to hold nom |
---|
63 | if (len_trim(firstname).lt.len_trim(name)) then |
---|
64 | write(*,*) "writediagsoil: Error !!!" |
---|
65 | write(*,*) " firstname string not long enough!!" |
---|
66 | write(*,*) " increase its size to at least ",len_trim(name) |
---|
67 | stop |
---|
68 | endif |
---|
69 | |
---|
70 | ! Set output sample rate |
---|
71 | isample=int(ecritphy) ! same as for diagfi outputs |
---|
72 | ! Note ecritphy is known from control.h |
---|
73 | |
---|
74 | ! Create output NetCDF file |
---|
75 | ierr=NF_CREATE(filename,NF_CLOBBER,nid) |
---|
76 | if (ierr.ne.NF_NOERR) then |
---|
77 | write(*,*)'writediagsoil: Error, failed creating file '//trim(filename) |
---|
78 | stop |
---|
79 | endif |
---|
80 | |
---|
81 | ! Define dimensions and axis attributes |
---|
82 | call iniwritesoil(nid) |
---|
83 | |
---|
84 | ! set zitau to -1 to be compatible with zitau incrementation step below |
---|
85 | zitau=-1 |
---|
86 | |
---|
87 | else |
---|
88 | ! If not an initialization call, simply open the NetCDF file |
---|
89 | ierr=NF_OPEN(filename,NF_WRITE,nid) |
---|
90 | endif ! of if (firstname.eq."1234567890") |
---|
91 | |
---|
92 | ! 2. Increment local time counter, if necessary |
---|
93 | if (name.eq.firstname) then |
---|
94 | ! if we run across 'firstname', then it is a new time step |
---|
95 | zitau=zitau+iphysiq |
---|
96 | ! Note iphysiq is known from control.h |
---|
97 | endif |
---|
98 | |
---|
99 | ! 3. Write data, if the time index matches the sample rate |
---|
100 | if (mod(zitau+1,isample).eq.0) then |
---|
101 | |
---|
102 | ! 3.1 If first call at this date, update 'time' variable |
---|
103 | if (name.eq.firstname) then |
---|
104 | ntime=ntime+1 |
---|
105 | date=real(zitau+1)/real(day_step) |
---|
106 | ! Note: day_step is known from control.h |
---|
107 | |
---|
108 | ! Get NetCDF ID for "time" |
---|
109 | ierr=NF_INQ_VARID(nid,"time",varid) |
---|
110 | ! Add the current value of date to the "time" array |
---|
111 | !#ifdef NC_DOUBLE |
---|
112 | ! ierr=NF_PUT_VARA_DOUBLE(nid,varid,ntime,1,date) |
---|
113 | !#else |
---|
114 | ierr=NF_PUT_VARA_REAL(nid,varid,ntime,1,date) |
---|
115 | !#endif |
---|
116 | if (ierr.ne.NF_NOERR) then |
---|
117 | write(*,*)"writediagsoil: Failed writing date to time variable" |
---|
118 | stop |
---|
119 | endif |
---|
120 | endif ! of if (name.eq.firstname) |
---|
121 | |
---|
122 | ! 3.2 Write the variable to the NetCDF file |
---|
123 | if (dimpx.eq.3) then ! Case of a 3D variable |
---|
124 | ! A. Recast data along 'dynamics' grid |
---|
125 | do l=1,nsoilmx |
---|
126 | ! handle the poles |
---|
127 | do i=1,iip1 |
---|
128 | data3(i,1,l)=px(1,l) |
---|
129 | data3(i,jjp1,l)=px(ngrid,l) |
---|
130 | enddo |
---|
131 | ! rest of the grid |
---|
132 | do j=2,jjm |
---|
133 | ig0=1+(j-2)*iim |
---|
134 | do i=1,iim |
---|
135 | data3(i,j,l)=px(ig0+i,l) |
---|
136 | enddo |
---|
137 | data3(iip1,j,l)=data3(1,j,l) ! extra (modulo) longitude |
---|
138 | enddo |
---|
139 | enddo |
---|
140 | |
---|
141 | ! B. Write (append) the variable to the NetCDF file |
---|
142 | ! B.1. Get the ID of the variable |
---|
143 | ierr=NF_INQ_VARID(nid,name,varid) |
---|
144 | if (ierr.ne.NF_NOERR) then |
---|
145 | ! If we failed geting the variable's ID, we assume it is because |
---|
146 | ! the variable doesn't exist yet and must be created. |
---|
147 | ! Start by obtaining corresponding dimensions IDs |
---|
148 | ierr=NF_INQ_DIMID(nid,"longitude",id(1)) |
---|
149 | ierr=NF_INQ_DIMID(nid,"latitude",id(2)) |
---|
150 | ierr=NF_INQ_DIMID(nid,"depth",id(3)) |
---|
151 | ierr=NF_INQ_DIMID(nid,"time",id(4)) |
---|
152 | ! Tell the world about it |
---|
153 | write(*,*) "=====================" |
---|
154 | write(*,*) "writediagsoil: creating variable "//trim(name) |
---|
155 | call def_var(nid,name,title,units,4,id,varid,ierr) |
---|
156 | endif ! of if (ierr.ne.NF_NOERR) |
---|
157 | |
---|
158 | ! B.2. Prepare things to be able to write/append the variable |
---|
159 | corners(1)=1 |
---|
160 | corners(2)=1 |
---|
161 | corners(3)=1 |
---|
162 | corners(4)=ntime |
---|
163 | |
---|
164 | edges(1)=iip1 |
---|
165 | edges(2)=jjp1 |
---|
166 | edges(3)=nsoilmx |
---|
167 | edges(4)=1 |
---|
168 | |
---|
169 | ! B.3. Write the slab of data |
---|
170 | !#ifdef NC_DOUBLE |
---|
171 | ! ierr=NF_PUT_VARA_DOUBLE(nid,varid,corners,edges,data3) |
---|
172 | !#else |
---|
173 | ierr=NF_PUT_VARA_REAL(nid,varid,corners,edges,data3) |
---|
174 | !#endif |
---|
175 | if (ierr.ne.NF_NOERR) then |
---|
176 | write(*,*) "writediagsoil: Error: Failed writing "//trim(name)//& |
---|
177 | " to file "//trim(filename)//" at time",date |
---|
178 | endif |
---|
179 | |
---|
180 | elseif (dimpx.eq.2) then ! Case of a 2D variable |
---|
181 | ! A. Recast data along 'dynamics' grid |
---|
182 | ! handle the poles |
---|
183 | do i=1,iip1 |
---|
184 | data2(i,1)=px(1,1) |
---|
185 | data2(i,jjp1)=px(ngrid,1) |
---|
186 | enddo |
---|
187 | ! rest of the grid |
---|
188 | do j=2,jjm |
---|
189 | ig0=1+(j-2)*iim |
---|
190 | do i=1,iim |
---|
191 | data2(i,j)=px(ig0+i,1) |
---|
192 | enddo |
---|
193 | data2(iip1,j)=data2(1,j) ! extra (modulo) longitude |
---|
194 | enddo |
---|
195 | |
---|
196 | ! B. Write (append) the variable to the NetCDF file |
---|
197 | ! B.1. Get the ID of the variable |
---|
198 | ierr=NF_INQ_VARID(nid,name,varid) |
---|
199 | if (ierr.ne.NF_NOERR) then |
---|
200 | ! If we failed geting the variable's ID, we assume it is because |
---|
201 | ! the variable doesn't exist yet and must be created. |
---|
202 | ! Start by obtaining corresponding dimensions IDs |
---|
203 | ierr=NF_INQ_DIMID(nid,"longitude",id(1)) |
---|
204 | ierr=NF_INQ_DIMID(nid,"latitude",id(2)) |
---|
205 | ierr=NF_INQ_DIMID(nid,"time",id(3)) |
---|
206 | ! Tell the world about it |
---|
207 | write(*,*) "=====================" |
---|
208 | write(*,*) "writediagsoil: creating variable "//trim(name) |
---|
209 | call def_var(nid,name,title,units,3,id,varid,ierr) |
---|
210 | endif ! of if (ierr.ne.NF_NOERR) |
---|
211 | |
---|
212 | ! B.2. Prepare things to be able to write/append the variable |
---|
213 | corners(1)=1 |
---|
214 | corners(2)=1 |
---|
215 | corners(3)=ntime |
---|
216 | |
---|
217 | edges(1)=iip1 |
---|
218 | edges(2)=jjp1 |
---|
219 | edges(3)=1 |
---|
220 | |
---|
221 | ! B.3. Write the slab of data |
---|
222 | !#ifdef NC_DOUBLE |
---|
223 | ! ierr=NF_PUT_VARA_DOUBLE(nid,varid,corners,edges,data2) |
---|
224 | !#else |
---|
225 | ierr=NF_PUT_VARA_REAL(nid,varid,corners,edges,data2) |
---|
226 | !#endif |
---|
227 | if (ierr.ne.NF_NOERR) then |
---|
228 | write(*,*) "writediagsoil: Error: Failed writing "//trim(name)//& |
---|
229 | " to file "//trim(filename)//" at time",date |
---|
230 | endif |
---|
231 | |
---|
232 | elseif (dimpx.eq.0) then ! Case of a 0D variable |
---|
233 | ! A. Copy data value |
---|
234 | data0=px(1,1) |
---|
235 | |
---|
236 | ! B. Write (append) the variable to the NetCDF file |
---|
237 | ! B.1. Get the ID of the variable |
---|
238 | ierr=NF_INQ_VARID(nid,name,varid) |
---|
239 | if (ierr.ne.NF_NOERR) then |
---|
240 | ! If we failed geting the variable's ID, we assume it is because |
---|
241 | ! the variable doesn't exist yet and must be created. |
---|
242 | ! Start by obtaining corresponding dimensions IDs |
---|
243 | ierr=NF_INQ_DIMID(nid,"time",id(1)) |
---|
244 | ! Tell the world about it |
---|
245 | write(*,*) "=====================" |
---|
246 | write(*,*) "writediagsoil: creating variable "//trim(name) |
---|
247 | call def_var(nid,name,title,units,1,id,varid,ierr) |
---|
248 | endif ! of if (ierr.ne.NF_NOERR) |
---|
249 | |
---|
250 | ! B.2. Prepare things to be able to write/append the variable |
---|
251 | corners(1)=ntime |
---|
252 | |
---|
253 | edges(1)=1 |
---|
254 | |
---|
255 | ! B.3. Write the data |
---|
256 | !#ifdef NC_DOUBLE |
---|
257 | ! ierr=NF_PUT_VARA_DOUBLE(nid,varid,corners,edges,data0) |
---|
258 | !#else |
---|
259 | ierr=NF_PUT_VARA_REAL(nid,varid,corners,edges,data0) |
---|
260 | !#endif |
---|
261 | if (ierr.ne.NF_NOERR) then |
---|
262 | write(*,*) "writediagsoil: Error: Failed writing "//trim(name)//& |
---|
263 | " to file "//trim(filename)//" at time",date |
---|
264 | endif |
---|
265 | |
---|
266 | endif ! of if (dimpx.eq.3) elseif (dimpx.eq.2) ... |
---|
267 | endif ! of if (mod(zitau+1,isample).eq.0) |
---|
268 | |
---|
269 | ! 4. Close the NetCDF file |
---|
270 | ierr=NF_CLOSE(nid) |
---|
271 | |
---|
272 | end subroutine writediagsoil |
---|