source: trunk/LMDZ.MARS/libf/phymars/writediagsoil.F90 @ 1528

Last change on this file since 1528 was 1528, checked in by emillour, 9 years ago

Mars GCM:

  • Got rid of references to "dimensions.h" from physics packages: use nbp_lon (=iim), nbp_lat (==jjp1) and nbp_lev from module mod_grid_phy_lmdz (in phy_common) instead.
  • Added "ioipsl_getin_p_mod.F90" (getin_p routine) in phy_common to correctly read in parameters from *.def files in a parallel environment.

EM

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