source: trunk/LMDZ.GENERIC/libf/phystd/iniwritesoil.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: 7.6 KB
Line 
1subroutine iniwritesoil(nid,ngrid)
2
3use comsoil_h
4
5! initialization routine for 'writediagoil'. Here we create/define
6! dimensions (longitude, latitude, depth and time) and other fixed
7! (time-independent) parameters.
8
9implicit none
10
11#include"dimensions.h"
12#include"dimphys.h"
13#include"paramet.h"
14#include"comcstfi.h"
15#include"comgeom.h"
16#include"netcdf.inc"
17
18! Arguments:
19integer,intent(in) :: ngrid
20integer,intent(in) :: nid ! NetCDF output file ID
21
22! Local variables:
23
24! NetCDF stuff:
25integer :: ierr ! NetCDF routines return code
26integer :: idim_rlatu ! ID of the 'latitude' dimension
27integer :: idim_rlonv ! ID of the 'longitude' dimension
28integer :: idim_depth ! ID of the 'depth' dimension
29integer :: idim_time  ! ID of the 'time' dimension
30integer :: varid ! to store NetCDF ID of a variable
31integer,dimension(3) :: dimids ! to store IDs of dimensions of a variable
32character(len=60) :: text ! to store some text
33real,dimension(iip1,jjp1,nsoilmx) :: data3 ! to store 3D data
34integer :: i,j,l,ig0
35
36! 1. Define the dimensions
37! Switch to NetCDF define mode
38ierr=NF_REDEF(nid)
39
40! Define the dimensions
41ierr=NF_DEF_DIM(nid,"longitude",iip1,idim_rlonv)
42! iip1 known from paramet.h
43if (ierr.ne.NF_NOERR) then
44  write(*,*)"iniwritesoil: Error, could not define longitude dimension"
45endif
46ierr=NF_DEF_DIM(nid,"latitude",jjp1,idim_rlatu)
47! jjp1 known from paramet.h
48if (ierr.ne.NF_NOERR) then
49  write(*,*)"iniwritesoil: Error, could not define latitude dimension"
50endif
51ierr=NF_DEF_DIM(nid,"depth",nsoilmx,idim_depth)
52! nsoilmx known from dimphys.h
53if (ierr.ne.NF_NOERR) then
54  write(*,*)"iniwritesoil: Error, could not define depth dimension"
55endif
56ierr=NF_DEF_DIM(nid,"time",NF_UNLIMITED,idim_time)
57if (ierr.ne.NF_NOERR) then
58  write(*,*)"iniwritesoil: Error, could not define time dimension"
59endif
60
61! Switch out of NetCDF define mode
62ierr=NF_ENDDEF(nid)
63
64! 2. Define (as variables) and write dimensions, as well as their attributes
65! 2.1. Longitude
66ierr=NF_REDEF(nid) ! switch to NetCDF define mode
67
68! Define the variable
69#ifdef NC_DOUBLE
70ierr=NF_DEF_VAR(nid,"longitude",NF_DOUBLE,1,idim_rlonv,varid)
71#else
72ierr=NF_DEF_VAR(nid,"longitude",NF_FLOAT,1,idim_rlonv,varid)
73#endif
74if (ierr.ne.NF_NOERR) then
75  write(*,*)"iniwritesoil: Error, could not define longitude variable"
76endif
77
78! Longitude attributes
79text="East longitude"
80ierr=NF_PUT_ATT_TEXT(nid,varid,"long_name",len_trim(text),text)
81text="degrees_east"
82ierr=NF_PUT_ATT_TEXT(nid,varid,"units",len_trim(text),text)
83
84! Write longitude to file
85ierr=NF_ENDDEF(nid) ! switch out of NetCDF define mode
86! Write
87#ifdef NC_DOUBLE
88ierr=NF_PUT_VAR_DOUBLE(nid,varid,rlonv*(180./pi))
89#else
90ierr=NF_PUT_VAR_REAL(nid,varid,rlonv*(180./pi))
91#endif
92! Note: rlonv is known from comgeom.h and pi from comcstfi.h
93if (ierr.ne.NF_NOERR) then
94  write(*,*)"iniwritesoil: Error, could not write longitude variable"
95endif
96
97! 2.2. Latitude
98ierr=NF_REDEF(nid) ! switch to NetCDF define mode
99
100! Define the variable
101#ifdef NC_DOUBLE
102ierr=NF_DEF_VAR(nid,"latitude",NF_DOUBLE,1,idim_rlatu,varid)
103#else
104ierr=NF_DEF_VAR(nid,"latitude",NF_FLOAT,1,idim_rlatu,varid)
105#endif
106if (ierr.ne.NF_NOERR) then
107  write(*,*)"iniwritesoil: Error, could not define latitude variable"
108endif
109
110! Latitude attributes
111text="North latitude"
112ierr=NF_PUT_ATT_TEXT(nid,varid,"long_name",len_trim(text),text)
113text="degrees_north"
114ierr=NF_PUT_ATT_TEXT(nid,varid,"units",len_trim(text),text)
115
116! Write latitude to file
117ierr=NF_ENDDEF(nid) ! switch out of NetCDF define mode
118! Write
119#ifdef NC_DOUBLE
120ierr=NF_PUT_VAR_DOUBLE(nid,varid,rlatu*(180./pi))
121#else
122ierr=NF_PUT_VAR_REAL(nid,varid,rlatu*(180./pi))
123#endif
124! Note: rlatu is known from comgeom.h and pi from comcstfi.h
125if (ierr.ne.NF_NOERR) then
126  write(*,*)"iniwritesoil: Error, could not write longitude variable"
127endif
128
129! 2.3. Depth
130ierr=NF_REDEF(nid) ! switch to NetCDF define mode
131
132! Define the variable
133#ifdef NC_DOUBLE
134ierr=NF_DEF_VAR(nid,"depth",NF_DOUBLE,1,idim_depth,varid)
135#else
136ierr=NF_DEF_VAR(nid,"depth",NF_FLOAT,1,idim_depth,varid)
137#endif
138if (ierr.ne.NF_NOERR) then
139  write(*,*)"iniwritesoil: Error, could not define depth variable"
140endif
141
142! Depth attributes
143text="Soil mid-layer depth"
144ierr=NF_PUT_ATT_TEXT(nid,varid,"long_name",len_trim(text),text)
145text="m"
146ierr=NF_PUT_ATT_TEXT(nid,varid,"units",len_trim(text),text)
147text="down"
148ierr=NF_PUT_ATT_TEXT(nid,varid,"positive",len_trim(text),text)
149
150! Write depth to file
151ierr=NF_ENDDEF(nid) ! switch out of NetCDF define mode
152! Write
153#ifdef NC_DOUBLE
154ierr=NF_PUT_VAR_DOUBLE(nid,varid,mlayer)
155#else
156ierr=NF_PUT_VAR_REAL(nid,varid,mlayer)
157#endif
158! Note mlayer(0:nsoilmx-1) known from comsoil.h
159if (ierr.ne.NF_NOERR) then
160  write(*,*)"iniwritesoil: Error, could not write depth variable"
161endif
162
163! 2.4. Time
164ierr=NF_REDEF(nid) ! switch to NetCDF define mode
165
166! Define the variable
167#ifdef NC_DOUBLE
168ierr=NF_DEF_VAR(nid,"time",NF_DOUBLE,1,idim_time,varid)
169#else
170ierr=NF_DEF_VAR(nid,"time",NF_FLOAT,1,idim_time,varid)
171#endif
172if (ierr.ne.NF_NOERR) then
173  write(*,*)"iniwritesoil: Error, could not define depth variable"
174endif
175
176! time attributes
177text="Time"
178ierr=NF_PUT_ATT_TEXT(nid,varid,"long_name",len_trim(text),text)
179text="days since 0000-01-01 00:00:00"
180ierr=NF_PUT_ATT_TEXT(nid,varid,"units",len_trim(text),text)
181
182ierr=NF_ENDDEF(nid) ! switch out of NetCDF define mode
183! Note no need to write time variable here; it is done in writediagsoil.
184
185! 3. Other variables to be included
186
187! 3.1 mesh area surrounding each horizontal point
188ierr=NF_REDEF(nid) ! switch to NetCDF define mode
189
190! Define the variable
191dimids(1)=idim_rlonv ! ID of the 'longitude' dimension
192dimids(2)=idim_rlatu ! ID of the 'latitude' dimension
193#ifdef NC_DOUBLE
194ierr=NF_DEF_VAR(nid,"area",NF_DOUBLE,2,dimids,varid)
195#else
196ierr=NF_DEF_VAR(nid,"area",NF_FLOAT,2,dimids,varid)
197#endif
198if (ierr.ne.NF_NOERR) then
199  write(*,*)"iniwritesoil: Error, could not define area variable"
200endif
201
202! Area attributes
203text="Mesh area"
204ierr=NF_PUT_ATT_TEXT(nid,varid,"long_name",len_trim(text),text)
205text="m2"
206ierr=NF_PUT_ATT_TEXT(nid,varid,"units",len_trim(text),text)
207
208! Write area to file
209ierr=NF_ENDDEF(nid) ! switch out of NetCDF define mode
210! Write
211#ifdef NC_DOUBLE
212ierr=NF_PUT_VAR_DOUBLE(nid,varid,aire)
213#else
214ierr=NF_PUT_VAR_REAL(nid,varid,aire)
215#endif
216! Note: aire is known from comgeom.h
217if (ierr.ne.NF_NOERR) then
218  write(*,*)"iniwritesoil: Error, could not write area variable"
219endif
220
221! 3.2 Thermal inertia
222ierr=NF_REDEF(nid) ! switch to NetCDF define mode
223
224! Define the variable
225dimids(1)=idim_rlonv ! ID of the 'longitude' dimension
226dimids(2)=idim_rlatu ! ID of the 'latitude' dimension
227dimids(3)=idim_depth ! ID of the 'depth' dimension
228#ifdef NC_DOUBLE
229ierr=NF_DEF_VAR(nid,"th_inertia",NF_DOUBLE,3,dimids,varid)
230#else
231ierr=NF_DEF_VAR(nid,"th_inertia",NF_FLOAT,3,dimids,varid)
232#endif
233if (ierr.ne.NF_NOERR) then
234  write(*,*)"iniwritesoil: Error, could not define th_inertia variable"
235endif
236
237! Attributes
238text="Thermal inertia"
239ierr=NF_PUT_ATT_TEXT(nid,varid,"long_name",len_trim(text),text)
240text="J.s-1/2.m-2.K-1"
241ierr=NF_PUT_ATT_TEXT(nid,varid,"units",len_trim(text),text)
242
243! Recast data along 'dynamics' grid
244! Note: inertiedat is known from comsoil.h
245
246do l=1,nsoilmx
247  ! handle the poles
248  do i=1,iip1
249    data3(i,1,l)=inertiedat(1,l)
250    data3(i,jjp1,l)=inertiedat(ngrid,l)
251  enddo
252  !!! THIS WILL NOT WORK IN PARALLEL !!!!
253  ! rest of the grid
254  do j=2,jjm
255    ig0=1+(j-2)*iim
256    do i=1,iim
257      data3(i,j,l)=inertiedat(ig0+i,l)
258    enddo
259    data3(iip1,j,l)=data3(1,j,l) ! extra (modulo) longitude
260  enddo
261enddo ! of do l=1,nsoilmx
262
263! Write data2 to file
264ierr=NF_ENDDEF(nid) ! switch out of NetCDF define mode
265! Write
266#ifdef NC_DOUBLE
267ierr=NF_PUT_VAR_DOUBLE(nid,varid,data3)
268#else
269ierr=NF_PUT_VAR_REAL(nid,varid,data3)
270#endif
271if (ierr.ne.NF_NOERR) then
272  write(*,*)"iniwritesoil: Error, could not write th_inertia variable"
273endif
274
275end subroutine iniwritesoil
Note: See TracBrowser for help on using the repository browser.