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 | use comsoil_h, only: nsoilmx, inertiedat |
---|
15 | use geometry_mod, only: cell_area |
---|
16 | use time_phylmdz_mod, only: ecritphy, day_step, iphysiq |
---|
17 | use mod_phys_lmdz_para, only : is_mpi_root, is_master, gather |
---|
18 | use mod_grid_phy_lmdz, only : klon_glo, Grid1Dto2D_glo, & |
---|
19 | nbp_lon, nbp_lat |
---|
20 | |
---|
21 | implicit none |
---|
22 | |
---|
23 | include"netcdf.inc" |
---|
24 | |
---|
25 | ! Arguments: |
---|
26 | integer,intent(in) :: ngrid ! number of (horizontal) points of physics grid |
---|
27 | ! i.e. ngrid = 2+(jjm-1)*iim - 1/jjm |
---|
28 | character(len=*),intent(in) :: name ! 'name' of the variable |
---|
29 | character(len=*),intent(in) :: title ! 'long_name' attribute of the variable |
---|
30 | character(len=*),intent(in) :: units ! 'units' attribute of the variable |
---|
31 | integer,intent(in) :: dimpx ! dimension of the variable (3,2 or 0) |
---|
32 | real,dimension(ngrid,nsoilmx),intent(in) :: px ! variable |
---|
33 | |
---|
34 | ! Local variables: |
---|
35 | real*4,dimension(nbp_lon+1,nbp_lat,nsoilmx) :: data3 ! to store 3D data |
---|
36 | real*4,dimension(nbp_lon+1,nbp_lat) :: data2 ! to store 2D data |
---|
37 | real*4 :: data0 ! to store 0D data |
---|
38 | real*4 :: data3_1d(1,nsoilmx) ! to store a profile in 1D model |
---|
39 | real*4 :: data2_1d ! to store surface value with 1D model |
---|
40 | integer :: i,j,l ! for loops |
---|
41 | integer :: ig0 |
---|
42 | |
---|
43 | real*4,save :: date ! time counter (in elapsed days) |
---|
44 | |
---|
45 | real :: inertia((nbp_lon+1),nbp_lat,nsoilmx) |
---|
46 | real :: area((nbp_lon+1),nbp_lat) |
---|
47 | |
---|
48 | real :: inertiafi_glo(klon_glo,nsoilmx) |
---|
49 | real :: areafi_glo(klon_glo) |
---|
50 | |
---|
51 | integer,save :: isample ! sample rate at which data is to be written to output |
---|
52 | integer,save :: ntime=0 ! counter to internally store time steps |
---|
53 | character(len=20),save :: firstname="1234567890" |
---|
54 | integer,save :: zitau=0 |
---|
55 | !$OMP THREADPRIVATE(date,isample,ntime,firstname,zitau) |
---|
56 | |
---|
57 | character(len=30) :: filename="diagsoil.nc" |
---|
58 | |
---|
59 | ! NetCDF stuff: |
---|
60 | integer :: nid ! NetCDF output file ID |
---|
61 | integer :: varid ! NetCDF ID of a variable |
---|
62 | integer :: ierr ! NetCDF routines return code |
---|
63 | integer,dimension(4) :: id ! NetCDF IDs of the dimensions of the variable |
---|
64 | integer,dimension(4) :: edges,corners |
---|
65 | |
---|
66 | #ifdef CPP_PARA |
---|
67 | ! Added to work in parallel mode |
---|
68 | real dx3_glop(klon_glo,nsoilmx) |
---|
69 | real dx3_glo(nbp_lon,nbp_lat,nsoilmx) ! to store a global 3D data set |
---|
70 | real dx2_glop(klon_glo) |
---|
71 | real dx2_glo(nbp_lon,nbp_lat) ! to store a global 2D (surface) data set |
---|
72 | real px2(ngrid) |
---|
73 | #endif |
---|
74 | |
---|
75 | ! 1. Initialization step |
---|
76 | if (firstname.eq."1234567890") then |
---|
77 | ! Store 'name' as 'firstname' |
---|
78 | firstname=name |
---|
79 | ! From now on, if 'name'.eq.'firstname', then it is a new time cycle |
---|
80 | |
---|
81 | ! just to be sure, check that firstnom is large enough to hold nom |
---|
82 | if (len_trim(firstname).lt.len_trim(name)) then |
---|
83 | write(*,*) "writediagsoil: Error !!!" |
---|
84 | write(*,*) " firstname string not long enough!!" |
---|
85 | write(*,*) " increase its size to at least ",len_trim(name) |
---|
86 | call abort_physic("writediagsoil","firstname too short",1) |
---|
87 | endif |
---|
88 | |
---|
89 | ! Set output sample rate |
---|
90 | isample=int(ecritphy) ! same as for diagfi outputs |
---|
91 | ! Note ecritphy is known from control.h |
---|
92 | |
---|
93 | ! Create output NetCDF file |
---|
94 | if (is_master) then |
---|
95 | ierr=NF_CREATE(filename,NF_CLOBBER,nid) |
---|
96 | if (ierr.ne.NF_NOERR) then |
---|
97 | write(*,*)'writediagsoil: Error, failed creating file '//trim(filename) |
---|
98 | call abort_physic("writediagsoil","failed creating"//trim(filename),1) |
---|
99 | endif |
---|
100 | endif |
---|
101 | |
---|
102 | #ifdef CPP_PARA |
---|
103 | ! Gather inertiedat() soil thermal inertia on physics grid |
---|
104 | call Gather(inertiedat,inertiafi_glo) |
---|
105 | ! Gather cell_area() mesh area on physics grid |
---|
106 | call Gather(cell_area,areafi_glo) |
---|
107 | #else |
---|
108 | inertiafi_glo(:,:)=inertiedat(:,:) |
---|
109 | areafi_glo(:)=cell_area(:) |
---|
110 | #endif |
---|
111 | |
---|
112 | if (is_master) then |
---|
113 | ! build inertia() and area() |
---|
114 | if (klon_glo>1) then |
---|
115 | do i=1,nbp_lon+1 ! poles |
---|
116 | inertia(i,1,1:nsoilmx)=inertiafi_glo(1,1:nsoilmx) |
---|
117 | inertia(i,nbp_lat,1:nsoilmx)=inertiafi_glo(klon_glo,1:nsoilmx) |
---|
118 | ! for area, divide at the poles by nbp_lon |
---|
119 | area(i,1)=areafi_glo(1)/nbp_lon |
---|
120 | area(i,nbp_lat)=areafi_glo(klon_glo)/nbp_lon |
---|
121 | enddo |
---|
122 | do j=2,nbp_lat-1 |
---|
123 | ig0= 1+(j-2)*nbp_lon |
---|
124 | do i=1,nbp_lon |
---|
125 | inertia(i,j,1:nsoilmx)=inertiafi_glo(ig0+i,1:nsoilmx) |
---|
126 | area(i,j)=areafi_glo(ig0+i) |
---|
127 | enddo |
---|
128 | ! handle redundant point in longitude |
---|
129 | inertia(nbp_lon+1,j,1:nsoilmx)=inertia(1,j,1:nsoilmx) |
---|
130 | area(nbp_lon+1,j)=area(1,j) |
---|
131 | enddo |
---|
132 | endif |
---|
133 | |
---|
134 | ! write "header" of file (longitudes, latitudes, geopotential, ...) |
---|
135 | if (klon_glo>1) then ! general 3D case |
---|
136 | call iniwritesoil(nid,ngrid,inertia,area,nbp_lon+1,nbp_lat) |
---|
137 | else ! 1D model |
---|
138 | call iniwritesoil(nid,ngrid,inertiafi_glo(1,:),areafi_glo(1),1,1) |
---|
139 | endif |
---|
140 | |
---|
141 | endif ! of if (is_master) |
---|
142 | |
---|
143 | ! set zitau to -1 to be compatible with zitau incrementation step below |
---|
144 | zitau=-1 |
---|
145 | |
---|
146 | else |
---|
147 | ! If not an initialization call, simply open the NetCDF file |
---|
148 | if (is_master) then |
---|
149 | ierr=NF_OPEN(filename,NF_WRITE,nid) |
---|
150 | endif |
---|
151 | endif ! of if (firstname.eq."1234567890") |
---|
152 | |
---|
153 | ! 2. Increment local time counter, if necessary |
---|
154 | if (name.eq.firstname) then |
---|
155 | ! if we run across 'firstname', then it is a new time step |
---|
156 | zitau=zitau+iphysiq |
---|
157 | ! Note iphysiq is known from control.h |
---|
158 | endif |
---|
159 | |
---|
160 | ! 3. Write data, if the time index matches the sample rate |
---|
161 | if (mod(zitau+1,isample).eq.0) then |
---|
162 | |
---|
163 | ! 3.1 If first call at this date, update 'time' variable |
---|
164 | if (name.eq.firstname) then |
---|
165 | ntime=ntime+1 |
---|
166 | date=float(zitau+1)/float(day_step) |
---|
167 | ! Note: day_step is known from control.h |
---|
168 | |
---|
169 | if (is_master) then |
---|
170 | ! Get NetCDF ID for "time" |
---|
171 | ierr=NF_INQ_VARID(nid,"time",varid) |
---|
172 | ! Add the current value of date to the "time" array |
---|
173 | !#ifdef NC_DOUBLE |
---|
174 | ! ierr=NF_PUT_VARA_DOUBLE(nid,varid,ntime,1,date) |
---|
175 | !#else |
---|
176 | ierr=NF_PUT_VARA_REAL(nid,varid,[ntime],[1],[date]) |
---|
177 | !#endif |
---|
178 | if (ierr.ne.NF_NOERR) then |
---|
179 | write(*,*)"writediagsoil: Failed writing date to time variable" |
---|
180 | call abort_physic("writediagsoil","failed writing time",1) |
---|
181 | endif |
---|
182 | endif ! of if (is_master) |
---|
183 | endif ! of if (name.eq.firstname) |
---|
184 | |
---|
185 | ! 3.2 Write the variable to the NetCDF file |
---|
186 | if (dimpx.eq.3) then ! Case of a 3D variable |
---|
187 | ! A. Recast data along 'dynamics' grid |
---|
188 | #ifdef CPP_PARA |
---|
189 | ! gather field on a "global" (without redundant longitude) array |
---|
190 | call Gather(px,dx3_glop) |
---|
191 | !$OMP MASTER |
---|
192 | if (is_mpi_root) then |
---|
193 | call Grid1Dto2D_glo(dx3_glop,dx3_glo) |
---|
194 | ! copy dx3_glo() to dx3(:) and add redundant longitude |
---|
195 | data3(1:nbp_lon,:,:)=dx3_glo(1:nbp_lon,:,:) |
---|
196 | data3(nbp_lon+1,:,:)=data3(1,:,:) |
---|
197 | endif |
---|
198 | !$OMP END MASTER |
---|
199 | !$OMP BARRIER |
---|
200 | #else |
---|
201 | if (klon_glo>1) then ! General case |
---|
202 | do l=1,nsoilmx |
---|
203 | ! handle the poles |
---|
204 | do i=1,nbp_lon+1 |
---|
205 | data3(i,1,l)=px(1,l) |
---|
206 | data3(i,nbp_lat,l)=px(ngrid,l) |
---|
207 | enddo |
---|
208 | ! rest of the grid |
---|
209 | do j=2,nbp_lat-1 |
---|
210 | ig0=1+(j-2)*nbp_lon |
---|
211 | do i=1,nbp_lon |
---|
212 | data3(i,j,l)=px(ig0+i,l) |
---|
213 | enddo |
---|
214 | data3(nbp_lon+1,j,l)=data3(1,j,l) ! extra (modulo) longitude |
---|
215 | enddo |
---|
216 | enddo |
---|
217 | else ! 1D model case |
---|
218 | data3_1d(1,1:nsoilmx)=px(1,1:nsoilmx) |
---|
219 | endif |
---|
220 | #endif |
---|
221 | |
---|
222 | ! B. Write (append) the variable to the NetCDF file |
---|
223 | if (is_master) then |
---|
224 | ! B.1. Get the ID of the variable |
---|
225 | ierr=NF_INQ_VARID(nid,name,varid) |
---|
226 | if (ierr.ne.NF_NOERR) then |
---|
227 | ! If we failed geting the variable's ID, we assume it is because |
---|
228 | ! the variable doesn't exist yet and must be created. |
---|
229 | ! Start by obtaining corresponding dimensions IDs |
---|
230 | ierr=NF_INQ_DIMID(nid,"longitude",id(1)) |
---|
231 | ierr=NF_INQ_DIMID(nid,"latitude",id(2)) |
---|
232 | ierr=NF_INQ_DIMID(nid,"depth",id(3)) |
---|
233 | ierr=NF_INQ_DIMID(nid,"time",id(4)) |
---|
234 | ! Tell the world about it |
---|
235 | write(*,*) "=====================" |
---|
236 | write(*,*) "writediagsoil: creating variable "//trim(name) |
---|
237 | call def_var(nid,name,title,units,4,id,varid,ierr) |
---|
238 | endif ! of if (ierr.ne.NF_NOERR) |
---|
239 | |
---|
240 | ! B.2. Prepare things to be able to write/append the variable |
---|
241 | corners(1)=1 |
---|
242 | corners(2)=1 |
---|
243 | corners(3)=1 |
---|
244 | corners(4)=ntime |
---|
245 | |
---|
246 | if (klon_glo==1) then |
---|
247 | edges(1)=1 |
---|
248 | else |
---|
249 | edges(1)=nbp_lon+1 |
---|
250 | endif |
---|
251 | edges(2)=nbp_lat |
---|
252 | edges(3)=nsoilmx |
---|
253 | edges(4)=1 |
---|
254 | |
---|
255 | ! B.3. Write the slab of data |
---|
256 | !#ifdef NC_DOUBLE |
---|
257 | ! ierr=NF_PUT_VARA_DOUBLE(nid,varid,corners,edges,data3) |
---|
258 | !#else |
---|
259 | if (klon_glo>1) then |
---|
260 | ierr=NF_PUT_VARA_REAL(nid,varid,corners,edges,data3) |
---|
261 | else |
---|
262 | ierr=NF_PUT_VARA_REAL(nid,varid,corners,edges,data3_1d) |
---|
263 | endif |
---|
264 | !#endif |
---|
265 | if (ierr.ne.NF_NOERR) then |
---|
266 | write(*,*) "writediagsoil: Error: Failed writing "//trim(name)//& |
---|
267 | " to file "//trim(filename)//" at time",date |
---|
268 | endif |
---|
269 | endif ! of if (is_master) |
---|
270 | |
---|
271 | elseif (dimpx.eq.2) then ! Case of a 2D variable |
---|
272 | |
---|
273 | ! A. Recast data along 'dynamics' grid |
---|
274 | #ifdef CPP_PARA |
---|
275 | ! gather field on a "global" (without redundant longitude) array |
---|
276 | px2(:)=px(:,1) |
---|
277 | call Gather(px2,dx2_glop) |
---|
278 | !$OMP MASTER |
---|
279 | if (is_mpi_root) then |
---|
280 | call Grid1Dto2D_glo(dx2_glop,dx2_glo) |
---|
281 | ! copy dx3_glo() to dx3(:) and add redundant longitude |
---|
282 | data2(1:nbp_lon,:)=dx2_glo(1:nbp_lon,:) |
---|
283 | data2(nbp_lon+1,:)=data2(1,:) |
---|
284 | endif |
---|
285 | !$OMP END MASTER |
---|
286 | !$OMP BARRIER |
---|
287 | #else |
---|
288 | if (klon_glo>1) then ! general case |
---|
289 | ! handle the poles |
---|
290 | do i=1,nbp_lon+1 |
---|
291 | data2(i,1)=px(1,1) |
---|
292 | data2(i,nbp_lat)=px(ngrid,1) |
---|
293 | enddo |
---|
294 | ! rest of the grid |
---|
295 | do j=2,nbp_lat-1 |
---|
296 | ig0=1+(j-2)*nbp_lon |
---|
297 | do i=1,nbp_lon |
---|
298 | data2(i,j)=px(ig0+i,1) |
---|
299 | enddo |
---|
300 | data2(nbp_lon+1,j)=data2(1,j) ! extra (modulo) longitude |
---|
301 | enddo |
---|
302 | else ! 1D model case |
---|
303 | data2_1d=px(1,1) |
---|
304 | endif |
---|
305 | #endif |
---|
306 | |
---|
307 | ! B. Write (append) the variable to the NetCDF file |
---|
308 | if (is_master) then |
---|
309 | ! B.1. Get the ID of the variable |
---|
310 | ierr=NF_INQ_VARID(nid,name,varid) |
---|
311 | if (ierr.ne.NF_NOERR) then |
---|
312 | ! If we failed geting the variable's ID, we assume it is because |
---|
313 | ! the variable doesn't exist yet and must be created. |
---|
314 | ! Start by obtaining corresponding dimensions IDs |
---|
315 | ierr=NF_INQ_DIMID(nid,"longitude",id(1)) |
---|
316 | ierr=NF_INQ_DIMID(nid,"latitude",id(2)) |
---|
317 | ierr=NF_INQ_DIMID(nid,"time",id(3)) |
---|
318 | ! Tell the world about it |
---|
319 | write(*,*) "=====================" |
---|
320 | write(*,*) "writediagsoil: creating variable "//trim(name) |
---|
321 | call def_var(nid,name,title,units,3,id,varid,ierr) |
---|
322 | endif ! of if (ierr.ne.NF_NOERR) |
---|
323 | |
---|
324 | ! B.2. Prepare things to be able to write/append the variable |
---|
325 | corners(1)=1 |
---|
326 | corners(2)=1 |
---|
327 | corners(3)=ntime |
---|
328 | |
---|
329 | if (klon_glo==1) then |
---|
330 | edges(1)=1 |
---|
331 | else |
---|
332 | edges(1)=nbp_lon+1 |
---|
333 | endif |
---|
334 | edges(2)=nbp_lat |
---|
335 | edges(3)=1 |
---|
336 | |
---|
337 | ! B.3. Write the slab of data |
---|
338 | !#ifdef NC_DOUBLE |
---|
339 | ! ierr=NF_PUT_VARA_DOUBLE(nid,varid,corners,edges,data2) |
---|
340 | !#else |
---|
341 | if (klon_glo>1) then ! General case |
---|
342 | ierr=NF_PUT_VARA_REAL(nid,varid,corners,edges,data2) |
---|
343 | else |
---|
344 | ierr=NF_PUT_VARA_REAL(nid,varid,corners,edges,[data2_1d]) |
---|
345 | endif |
---|
346 | !#endif |
---|
347 | if (ierr.ne.NF_NOERR) then |
---|
348 | write(*,*) "writediagsoil: Error: Failed writing "//trim(name)//& |
---|
349 | " to file "//trim(filename)//" at time",date |
---|
350 | endif |
---|
351 | endif ! of if (is_master) |
---|
352 | |
---|
353 | elseif (dimpx.eq.0) then ! Case of a 0D variable |
---|
354 | #ifdef CPP_PARA |
---|
355 | write(*,*) "writediagsoil: dimps==0 case not implemented in // mode!!" |
---|
356 | call abort_physic("writediagsoil","dimps==0 not implemented",1) |
---|
357 | #endif |
---|
358 | ! A. Copy data value |
---|
359 | data0=px(1,1) |
---|
360 | |
---|
361 | ! B. Write (append) the variable to the NetCDF file |
---|
362 | ! B.1. Get the ID of the variable |
---|
363 | ierr=NF_INQ_VARID(nid,name,varid) |
---|
364 | if (ierr.ne.NF_NOERR) then |
---|
365 | ! If we failed geting the variable's ID, we assume it is because |
---|
366 | ! the variable doesn't exist yet and must be created. |
---|
367 | ! Start by obtaining corresponding dimensions IDs |
---|
368 | ierr=NF_INQ_DIMID(nid,"time",id(1)) |
---|
369 | ! Tell the world about it |
---|
370 | write(*,*) "=====================" |
---|
371 | write(*,*) "writediagsoil: creating variable "//trim(name) |
---|
372 | call def_var(nid,name,title,units,1,id,varid,ierr) |
---|
373 | endif ! of if (ierr.ne.NF_NOERR) |
---|
374 | |
---|
375 | ! B.2. Prepare things to be able to write/append the variable |
---|
376 | corners(1)=ntime |
---|
377 | |
---|
378 | edges(1)=1 |
---|
379 | |
---|
380 | ! B.3. Write the data |
---|
381 | !#ifdef NC_DOUBLE |
---|
382 | ! ierr=NF_PUT_VARA_DOUBLE(nid,varid,corners,edges,data0) |
---|
383 | !#else |
---|
384 | ierr=NF_PUT_VARA_REAL(nid,varid,corners,edges,[data0]) |
---|
385 | !#endif |
---|
386 | if (ierr.ne.NF_NOERR) then |
---|
387 | write(*,*) "writediagsoil: Error: Failed writing "//trim(name)//& |
---|
388 | " to file "//trim(filename)//" at time",date |
---|
389 | endif |
---|
390 | |
---|
391 | endif ! of if (dimpx.eq.3) elseif (dimpx.eq.2) ... |
---|
392 | endif ! of if (mod(zitau+1,isample).eq.0) |
---|
393 | |
---|
394 | ! 4. Close the NetCDF file |
---|
395 | if (is_master) then |
---|
396 | ierr=NF_CLOSE(nid) |
---|
397 | endif |
---|
398 | |
---|
399 | end subroutine writediagsoil |
---|