source: trunk/LMDZ.GENERIC/libf/phystd/writediagsoil.F90 @ 848

Last change on this file since 848 was 787, checked in by aslmd, 12 years ago

LMDZ.GENERIC. (Sorry for long text but this is a quite major commit)

Paving the path for parallel computations. And moving towards a more flexible code.

Automatic allocation is used within all routines in phystd. No further mention to ngridmx and nqmx.

  1. ngridmx and nqmx are still used in LMDZ.GENERIC in the dyn3d part
  2. if the LMDZ4/LMDZ5 dynamical core is used, there is no more fixed dimensions ngridmx and nqmx --> a fully flexible parallel implementation is now possible (e.g. no need to recompile when changing numbers of processors)

The important stuff :

  • Compilation checked with ifort. OK with and without debug mode. No errors. Checked for: gcm, newstart, rcm1d, kcm1d
  • RUN GCM: Running an Earth test case. Comparison with previous revision --> debug mode : perfect match. bit by bit (diff command). checked with plots --> O1 mode : close match (checked with plots) --> O2 mode : sometimes up to 0.5 K departure.... BUT in this new version O2 and O1 are quite close while in previous version O1 and O2 differed by about, well, typically 0.5 K (pictures available on request)
  • RUN NEWSTART : perfect match (bit-by-bit) in either debug or normal mode.
  • RUN RCM1D : perfect match in normal mode.
  • RUN KCM1D : not tested (I don't know what is the use of kcm1d)

List of main changes :

  • Additional arguments to some subroutines (ngrid and nq)
  • F77 include strategy is obsolete and replaced by F90 module strategy In this new strategy arrays are allocatable and allocated once at first use This has to be done for all common featuring arrays defined with ngridmx or nqmx

surfdat.h >> surfdat_h.F90
tracer.h >> tracer_h.F90
comsaison.h >> comsaison_h.F90
comgeomfi.h >> comgeomfi_h.F90
comsoil.h >> comsoil_h.F90
comdiurn.h >> comdiurn_h.F90
fisice.h >> DELETED. was not used. probably a fossil.
watercap.h >> DELETED. variable put in surfdat_h.F90

  • F77 'save' strategy is obsolete and replaced by F90 'allocatable save' strategy (see previous point and e.g. new version of physiq.F90)
  • Suppressing any mention to advtrac.h which is a common in the dynamics and needs nqmx This was easily solved by adding an argument with tracer names, coming from the dynamics This is probably not a definitive solution, ... but this allows for generic physics to work easily with either LMDZ.GENERIC or LMDZ dynamical cores
  • Removing consistency tests between nq and nqmx ; and ngrid and ngridmx. No use now!
  • Adaptation of rcm1d, kcm1d, newstart given above-mentioned changes

A note on phyetat0 and soil_setting:

  • Now written so that a slice of horizontal size 'ngrid' starting at grid point 'cursor' is read in startfi.nc 'cursor' is defined in dimphys.h and initialized by inifis (or in newstart) this is useful for parallel computations. default behavior is the usual one : sequential runs, cursor is 1, size ngrid is the whole global domain

A note on an additional change :

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