1 | PROGRAM teststomate |
---|
2 | !--------------------------------------------------------------------- |
---|
3 | ! This program reads a forcing file for STOMATE |
---|
4 | ! which was created with STOMATE coupled to SECHIBA. |
---|
5 | ! It then integrates STOMATE for a |
---|
6 | ! given number of years using this forcing. |
---|
7 | ! The GPP is scaled to the new calculated LAI... |
---|
8 | ! Nothing is done for soil moisture which should |
---|
9 | ! actually evolve with the vegetation. |
---|
10 | !--------------------------------------------------------------------- |
---|
11 | !- IPSL (2006) |
---|
12 | !- This software is governed by the CeCILL licence see ORCHIDEE/ORCHIDEE_CeCILL.LIC |
---|
13 | USE netcdf |
---|
14 | !- |
---|
15 | USE defprec |
---|
16 | USE constantes |
---|
17 | USE constantes_soil |
---|
18 | USE constantes_veg |
---|
19 | USE constantes_co2 |
---|
20 | USE stomate_constants |
---|
21 | USE ioipsl |
---|
22 | USE grid |
---|
23 | USE slowproc |
---|
24 | USE stomate |
---|
25 | USE intersurf, ONLY: stom_define_history , intsurf_time |
---|
26 | USE parallel |
---|
27 | !- |
---|
28 | IMPLICIT NONE |
---|
29 | !- |
---|
30 | ! Declarations |
---|
31 | !- |
---|
32 | INTEGER(i_std) :: vegax_id |
---|
33 | INTEGER(i_std) :: kjpij,kjpindex |
---|
34 | REAL(r_std) :: dtradia,dt_force |
---|
35 | INTEGER(i_std),DIMENSION(:),ALLOCATABLE :: indices |
---|
36 | INTEGER(i_std),DIMENSION(:),ALLOCATABLE :: indexveg |
---|
37 | REAL(r_std),DIMENSION(:,:),ALLOCATABLE :: soiltype, veget_x |
---|
38 | REAL(r_std),DIMENSION(:),ALLOCATABLE :: totfrac_nobio |
---|
39 | REAL(r_std),DIMENSION(:,:),ALLOCATABLE :: frac_nobio |
---|
40 | REAL(r_std),DIMENSION(:,:),ALLOCATABLE :: veget_max_x |
---|
41 | REAL(r_std),DIMENSION(:,:),ALLOCATABLE :: lai_x |
---|
42 | REAL(r_std),DIMENSION(:,:),ALLOCATABLE :: veget_force_x |
---|
43 | REAL(r_std),DIMENSION(:,:),ALLOCATABLE :: veget_max_force_x |
---|
44 | REAL(r_std),DIMENSION(:,:),ALLOCATABLE :: lai_force_x |
---|
45 | REAL(r_std),DIMENSION(:,:),ALLOCATABLE :: lon,lat |
---|
46 | REAL(r_std),DIMENSION(:),ALLOCATABLE :: t2m,t2m_min,temp_sol |
---|
47 | REAL(r_std),DIMENSION(:,:),ALLOCATABLE :: soiltemp,soilhum |
---|
48 | REAL(r_std),DIMENSION(:,:),ALLOCATABLE :: humrel_x |
---|
49 | REAL(r_std),DIMENSION(:),ALLOCATABLE :: litterhum |
---|
50 | REAL(r_std),DIMENSION(:),ALLOCATABLE :: precip_rain,precip_snow |
---|
51 | REAL(r_std),DIMENSION(:,:),ALLOCATABLE :: gpp_x |
---|
52 | REAL(r_std),DIMENSION(:),ALLOCATABLE :: deadleaf_cover |
---|
53 | REAL(r_std),DIMENSION(:,:,:),ALLOCATABLE :: assim_param_x |
---|
54 | REAL(r_std),DIMENSION(:,:),ALLOCATABLE :: height_x |
---|
55 | REAL(r_std),DIMENSION(:,:),ALLOCATABLE :: qsintmax_x |
---|
56 | REAL(r_std),DIMENSION(:,:),ALLOCATABLE :: co2_flux |
---|
57 | INTEGER :: ier,iret |
---|
58 | INTEGER :: ncfid |
---|
59 | LOGICAL :: a_er |
---|
60 | CHARACTER(LEN=80) :: & |
---|
61 | & dri_restname_in,dri_restname_out, & |
---|
62 | & sec_restname_in,sec_restname_out, & |
---|
63 | & sto_restname_in,sto_restname_out, stom_histname |
---|
64 | INTEGER(i_std) :: iim,jjm |
---|
65 | INTEGER(i_std),PARAMETER :: llm = 1 |
---|
66 | REAL(r_std),DIMENSION(llm) :: lev |
---|
67 | REAL(r_std) :: dt_files |
---|
68 | INTEGER(i_std) :: itau_dep,itau,itau_len,itau_step |
---|
69 | REAL(r_std) :: date0 |
---|
70 | INTEGER(i_std) :: rest_id_sec,rest_id_sto |
---|
71 | INTEGER(i_std) :: hist_id_sec,hist_id_sec2,hist_id_sto,hist_id_stom_IPCC |
---|
72 | CHARACTER(LEN=30) :: time_str |
---|
73 | REAL :: hist_days_stom,hist_dt_stom |
---|
74 | REAL,DIMENSION(nvm) :: hist_PFTaxis |
---|
75 | REAL(r_std),DIMENSION(10) :: hist_pool_10axis |
---|
76 | REAL(r_std),DIMENSION(100) :: hist_pool_100axis |
---|
77 | REAL(r_std),DIMENSION(11) :: hist_pool_11axis |
---|
78 | REAL(r_std),DIMENSION(101) :: hist_pool_101axis |
---|
79 | INTEGER :: hist_PFTaxis_id,hori_id |
---|
80 | INTEGER :: hist_pool_10axis_id |
---|
81 | INTEGER :: hist_pool_100axis_id |
---|
82 | INTEGER :: hist_pool_11axis_id |
---|
83 | INTEGER :: hist_pool_101axis_id |
---|
84 | INTEGER :: hist_level |
---|
85 | INTEGER,PARAMETER :: max_hist_level = 10 |
---|
86 | INTEGER(i_std) :: i,j,iv,id |
---|
87 | CHARACTER*80 :: var_name |
---|
88 | CHARACTER(LEN=40),DIMENSION(10) :: fluxop |
---|
89 | LOGICAL :: ldrestart_read,ldrestart_write |
---|
90 | LOGICAL :: l1d |
---|
91 | INTEGER(i_std),PARAMETER :: nbvarmax=200 |
---|
92 | INTEGER(i_std) :: nbvar |
---|
93 | CHARACTER*50,DIMENSION(nbvarmax) :: varnames |
---|
94 | INTEGER(i_std) :: varnbdim |
---|
95 | INTEGER(i_std),PARAMETER :: varnbdim_max=20 |
---|
96 | INTEGER,DIMENSION(varnbdim_max) :: vardims |
---|
97 | CHARACTER*200 :: taboo_vars |
---|
98 | REAL(r_std),DIMENSION(1) :: xtmp |
---|
99 | INTEGER :: nsfm,nsft |
---|
100 | INTEGER :: iisf |
---|
101 | INTEGER(i_std) :: max_totsize,totsize_1step |
---|
102 | INTEGER(i_std),DIMENSION(0:2) :: ifirst,ilast |
---|
103 | INTEGER(i_std) :: iblocks,nblocks |
---|
104 | INTEGER,PARAMETER :: ndm = 10 |
---|
105 | INTEGER,DIMENSION(ndm) :: start,count |
---|
106 | INTEGER :: ndim,v_id |
---|
107 | INTEGER :: force_id |
---|
108 | CHARACTER(LEN=100) :: forcing_name |
---|
109 | REAL :: x |
---|
110 | REAL(r_std),DIMENSION(:),ALLOCATABLE :: x_indices |
---|
111 | REAL(r_std),DIMENSION(:,:),ALLOCATABLE :: x_neighbours |
---|
112 | REAL(r_std),DIMENSION(:,:),ALLOCATABLE :: var_3d |
---|
113 | !- |
---|
114 | REAL(r_std) :: time_sec,time_step_sec |
---|
115 | REAL(r_std) :: time_dri,time_step_dri |
---|
116 | REAL(r_std),DIMENSION(1) :: r1d |
---|
117 | REAL(r_std) :: julian,djulian |
---|
118 | ! REAL(r_std),DIMENSION(:,:),ALLOCATABLE :: soiltype |
---|
119 | !- |
---|
120 | ! the following variables contain the forcing data |
---|
121 | !- |
---|
122 | REAL(r_std),ALLOCATABLE,SAVE,DIMENSION(:,:) :: clay_fm |
---|
123 | REAL(r_std),ALLOCATABLE,SAVE,DIMENSION(:,:,:) :: humrel_x_fm |
---|
124 | REAL(r_std),ALLOCATABLE,SAVE,DIMENSION(:,:) :: litterhum_fm |
---|
125 | REAL(r_std),ALLOCATABLE,SAVE,DIMENSION(:,:) :: t2m_fm |
---|
126 | REAL(r_std),ALLOCATABLE,SAVE,DIMENSION(:,:) :: t2m_min_fm |
---|
127 | REAL(r_std),ALLOCATABLE,SAVE,DIMENSION(:,:) :: temp_sol_fm |
---|
128 | REAL(r_std),ALLOCATABLE,SAVE,DIMENSION(:,:,:) :: soiltemp_fm |
---|
129 | REAL(r_std),ALLOCATABLE,SAVE,DIMENSION(:,:,:) :: soilhum_fm |
---|
130 | REAL(r_std),ALLOCATABLE,SAVE,DIMENSION(:,:) :: precip_fm |
---|
131 | REAL(r_std),ALLOCATABLE,SAVE,DIMENSION(:,:,:) :: gpp_x_fm |
---|
132 | REAL(r_std),ALLOCATABLE,SAVE,DIMENSION(:,:,:) :: veget_force_x_fm |
---|
133 | REAL(r_std),ALLOCATABLE,SAVE,DIMENSION(:,:,:) :: veget_max_force_x_fm |
---|
134 | REAL(r_std),ALLOCATABLE,SAVE,DIMENSION(:,:,:) :: lai_force_x_fm |
---|
135 | INTEGER(i_std),ALLOCATABLE,SAVE,DIMENSION(:) :: isf |
---|
136 | LOGICAL,ALLOCATABLE,SAVE,DIMENSION(:) :: nf_written |
---|
137 | INTEGER(i_std),ALLOCATABLE,SAVE,DIMENSION(:) :: nf_cumul |
---|
138 | INTEGER(i_std) :: ji,jv |
---|
139 | !--------------------------------------------------------------------- |
---|
140 | !- No parallelisation yet in teststomate ! |
---|
141 | #ifdef CPP_PARA |
---|
142 | CALL ipslerr (3,'teststomate', & |
---|
143 | & 'You try to run testsomate compiled with parallelisation. (CPP_PARA key)', & |
---|
144 | & 'But it wasn''t programmed yet and teststomate will stop.','You must compiled it without CPP_PARA key.') |
---|
145 | #endif |
---|
146 | |
---|
147 | !- |
---|
148 | ! calendar |
---|
149 | !- |
---|
150 | CALL ioconf_calendar ('noleap') |
---|
151 | CALL ioget_calendar (one_year,one_day) |
---|
152 | !- |
---|
153 | ! open STOMATE's forcing file to read some basic info |
---|
154 | !- |
---|
155 | forcing_name = 'stomate_forcing.nc' |
---|
156 | CALL getin ('STOMATE_FORCING_NAME',forcing_name) |
---|
157 | iret = NF90_OPEN (TRIM(forcing_name),NF90_NOWRITE,force_id) |
---|
158 | IF (iret /= NF90_NOERR) THEN |
---|
159 | CALL ipslerr (3,'teststomate', & |
---|
160 | & 'Could not open file : ', & |
---|
161 | & forcing_name,'(Do you have forget it ?)') |
---|
162 | ENDIF |
---|
163 | ier = NF90_GET_ATT (force_id,NF90_GLOBAL,'dtradia',dtradia) |
---|
164 | ier = NF90_GET_ATT (force_id,NF90_GLOBAL,'dt_slow',dt_force) |
---|
165 | ier = NF90_GET_ATT (force_id,NF90_GLOBAL,'nsft',x) |
---|
166 | nsft = NINT(x) |
---|
167 | ier = NF90_GET_ATT (force_id,NF90_GLOBAL,'kjpij',x) |
---|
168 | kjpij = NINT(x) |
---|
169 | ier = NF90_GET_ATT (force_id,NF90_GLOBAL,'kjpindex',x) |
---|
170 | kjpindex = NINT(x) |
---|
171 | CALL init_grid( kjpindex ) |
---|
172 | !- |
---|
173 | write(*,*) 'ATTENTION',dtradia,dt_force |
---|
174 | !- |
---|
175 | ! allocate arrays |
---|
176 | !- |
---|
177 | a_er = .FALSE. |
---|
178 | ALLOCATE (indices(kjpindex),stat=ier) |
---|
179 | a_er = a_er .OR. (ier.NE.0) |
---|
180 | ALLOCATE (indexveg(kjpindex*nvm), stat=ier) |
---|
181 | a_er = a_er .OR. (ier.NE.0) |
---|
182 | ALLOCATE (soiltype(kjpindex,nstm),stat=ier) |
---|
183 | a_er = a_er .OR. (ier.NE.0) |
---|
184 | ALLOCATE (veget_x(kjpindex,nvm),stat=ier) |
---|
185 | a_er = a_er .OR. (ier.NE.0) |
---|
186 | ALLOCATE (totfrac_nobio(kjpindex),stat=ier) |
---|
187 | a_er = a_er .OR. (ier.NE.0) |
---|
188 | ALLOCATE (frac_nobio(kjpindex,nnobio),stat=ier) |
---|
189 | a_er = a_er .OR. (ier.NE.0) |
---|
190 | ALLOCATE (veget_max_x(kjpindex,nvm),stat=ier) |
---|
191 | a_er = a_er .OR. (ier.NE.0) |
---|
192 | ALLOCATE (lai_x(kjpindex,nvm),stat=ier) |
---|
193 | a_er = a_er .OR. (ier.NE.0) |
---|
194 | ALLOCATE (veget_force_x(kjpindex,nvm),stat=ier) |
---|
195 | a_er = a_er .OR. (ier.NE.0) |
---|
196 | ALLOCATE (veget_max_force_x(kjpindex,nvm),stat=ier) |
---|
197 | a_er = a_er .OR. (ier.NE.0) |
---|
198 | ALLOCATE (lai_force_x(kjpindex,nvm),stat=ier) |
---|
199 | a_er = a_er .OR. (ier.NE.0) |
---|
200 | ALLOCATE (t2m(kjpindex),stat=ier) |
---|
201 | a_er = a_er .OR. (ier.NE.0) |
---|
202 | ALLOCATE (t2m_min(kjpindex),stat=ier) |
---|
203 | a_er = a_er .OR. (ier.NE.0) |
---|
204 | ALLOCATE (temp_sol(kjpindex),stat=ier) |
---|
205 | a_er = a_er .OR. (ier.NE.0) |
---|
206 | ALLOCATE (soiltemp(kjpindex,nbdl),stat=ier) |
---|
207 | a_er = a_er .OR. (ier.NE.0) |
---|
208 | ALLOCATE (soilhum(kjpindex,nbdl),stat=ier) |
---|
209 | a_er = a_er .OR. (ier.NE.0) |
---|
210 | ALLOCATE (humrel_x(kjpindex,nvm),stat=ier) |
---|
211 | a_er = a_er .OR. (ier.NE.0) |
---|
212 | ALLOCATE (litterhum(kjpindex),stat=ier) |
---|
213 | a_er = a_er .OR. (ier.NE.0) |
---|
214 | ALLOCATE (precip_rain(kjpindex),stat=ier) |
---|
215 | a_er = a_er .OR. (ier.NE.0) |
---|
216 | ALLOCATE (precip_snow(kjpindex),stat=ier) |
---|
217 | a_er = a_er .OR. (ier.NE.0) |
---|
218 | ALLOCATE (gpp_x(kjpindex,nvm),stat=ier) |
---|
219 | a_er = a_er .OR. (ier.NE.0) |
---|
220 | ALLOCATE (deadleaf_cover(kjpindex),stat=ier) |
---|
221 | a_er = a_er .OR. (ier.NE.0) |
---|
222 | ALLOCATE (assim_param_x(kjpindex,nvm,npco2),stat=ier) |
---|
223 | a_er = a_er .OR. (ier.NE.0) |
---|
224 | ALLOCATE (height_x(kjpindex,nvm),stat=ier) |
---|
225 | a_er = a_er .OR. (ier.NE.0) |
---|
226 | ALLOCATE (qsintmax_x(kjpindex,nvm),stat=ier) |
---|
227 | a_er = a_er .OR. (ier.NE.0) |
---|
228 | ALLOCATE (co2_flux(kjpindex,nvm),stat=ier) |
---|
229 | a_er = a_er .OR. (ier.NE.0) |
---|
230 | IF ( a_er ) STOP 'PROBLEM WITH ALLOCATION' |
---|
231 | ! |
---|
232 | ! prepare forcing |
---|
233 | ! |
---|
234 | max_totsize = 50 |
---|
235 | CALL getin ('STOMATE_FORCING_MEMSIZE',max_totsize) |
---|
236 | max_totsize = max_totsize * 1000000 |
---|
237 | totsize_1step = SIZE(soiltype(:,3))*KIND(soiltype(:,3)) + & |
---|
238 | SIZE(humrel_x)*KIND(humrel_x) + & |
---|
239 | SIZE(litterhum)*KIND(litterhum) + & |
---|
240 | SIZE(t2m)*KIND(t2m) + & |
---|
241 | SIZE(t2m_min)*KIND(t2m_min) + & |
---|
242 | SIZE(temp_sol)*KIND(temp_sol) + & |
---|
243 | SIZE(soiltemp)*KIND(soiltemp) + & |
---|
244 | SIZE(soilhum)*KIND(soilhum) + & |
---|
245 | SIZE(precip_rain)*KIND(precip_rain) + & |
---|
246 | SIZE(precip_snow)*KIND(precip_snow) + & |
---|
247 | SIZE(gpp_x)*KIND(gpp_x) + & |
---|
248 | SIZE(veget_force_x)*KIND(veget_force_x) + & |
---|
249 | SIZE(veget_max_force_x)*KIND(veget_max_force_x) + & |
---|
250 | SIZE(lai_force_x)*KIND(lai_force_x) |
---|
251 | ! total number of forcing steps |
---|
252 | nsft = INT(one_year/(dt_force/one_day)) |
---|
253 | ! number of forcing steps in memory |
---|
254 | nsfm = MIN(nsft,MAX(1,NINT(FLOAT(max_totsize)/FLOAT(totsize_1step)))) |
---|
255 | !- |
---|
256 | WRITE(numout,*) 'Offline forcing of Stomate:' |
---|
257 | WRITE(numout,*) ' Total number of forcing steps:',nsft |
---|
258 | WRITE(numout,*) ' Number of forcing steps in memory:',nsfm |
---|
259 | !- |
---|
260 | a_er = .FALSE. |
---|
261 | ALLOCATE (clay_fm(kjpindex,nsfm),stat=ier) |
---|
262 | a_er = a_er.OR.(ier.NE.0) |
---|
263 | ALLOCATE (humrel_x_fm(kjpindex,nvm,nsfm),stat=ier) |
---|
264 | a_er = a_er.OR.(ier.NE.0) |
---|
265 | ALLOCATE (litterhum_fm(kjpindex,nsfm),stat=ier) |
---|
266 | a_er = a_er.OR.(ier.NE.0) |
---|
267 | ALLOCATE (t2m_fm(kjpindex,nsfm),stat=ier) |
---|
268 | a_er = a_er.OR.(ier.NE.0) |
---|
269 | ALLOCATE (t2m_min_fm(kjpindex,nsfm),stat=ier) |
---|
270 | a_er = a_er.OR.(ier.NE.0) |
---|
271 | ALLOCATE (temp_sol_fm(kjpindex,nsfm),stat=ier) |
---|
272 | a_er = a_er.OR.(ier.NE.0) |
---|
273 | ALLOCATE (soiltemp_fm(kjpindex,nbdl,nsfm),stat=ier) |
---|
274 | a_er = a_er.OR.(ier.NE.0) |
---|
275 | ALLOCATE (soilhum_fm(kjpindex,nbdl,nsfm),stat=ier) |
---|
276 | a_er = a_er.OR.(ier.NE.0) |
---|
277 | ALLOCATE (precip_fm(kjpindex,nsfm),stat=ier) |
---|
278 | a_er = a_er.OR.(ier.NE.0) |
---|
279 | ALLOCATE (gpp_x_fm(kjpindex,nvm,nsfm),stat=ier) |
---|
280 | a_er = a_er.OR.(ier.NE.0) |
---|
281 | ALLOCATE (veget_force_x_fm(kjpindex,nvm,nsfm),stat=ier) |
---|
282 | a_er = a_er.OR.(ier.NE.0) |
---|
283 | ALLOCATE (veget_max_force_x_fm(kjpindex,nvm,nsfm),stat=ier) |
---|
284 | a_er = a_er.OR.(ier.NE.0) |
---|
285 | ALLOCATE (lai_force_x_fm(kjpindex,nvm,nsfm),stat=ier) |
---|
286 | a_er = a_er.OR.(ier.NE.0) |
---|
287 | ALLOCATE (isf(nsfm),stat=ier) |
---|
288 | a_er = a_er.OR.(ier.NE.0) |
---|
289 | ALLOCATE (nf_written(nsft),stat=ier) |
---|
290 | a_er = a_er.OR.(ier.NE.0) |
---|
291 | ALLOCATE (nf_cumul(nsft),stat=ier) |
---|
292 | a_er = a_er.OR.(ier.NE.0) |
---|
293 | IF (a_er) THEN |
---|
294 | STOP 'stomate: error in memory allocation for forcing data' |
---|
295 | ENDIF |
---|
296 | ! ensure that we read all new forcing states |
---|
297 | iisf = nsfm |
---|
298 | ! initialize that table that contains the indices |
---|
299 | ! of the forcing states that will be in memory |
---|
300 | isf(:) = (/ (i,i=1,nsfm) /) |
---|
301 | !- |
---|
302 | ! read info about grids |
---|
303 | !- |
---|
304 | contfrac(:) = 1.0 |
---|
305 | !- |
---|
306 | ALLOCATE (x_indices(kjpindex),stat=ier) |
---|
307 | ier = NF90_INQ_VARID (force_id,'index',v_id) |
---|
308 | ier = NF90_GET_VAR (force_id,v_id,x_indices) |
---|
309 | indices(:) = NINT(x_indices(:)) |
---|
310 | DEALLOCATE (x_indices) |
---|
311 | !- |
---|
312 | DO ji=1,kjpindex |
---|
313 | DO jv=1,nvm |
---|
314 | indexveg((jv-1)*kjpindex+ji) = indices(ji)+(jv-1)*kjpij |
---|
315 | ENDDO |
---|
316 | ENDDO |
---|
317 | !- |
---|
318 | ier = NF90_INQ_VARID (force_id,'lalo',v_id) |
---|
319 | ier = NF90_GET_VAR (force_id,v_id,lalo) |
---|
320 | !- |
---|
321 | ALLOCATE (x_neighbours(kjpindex,8),stat=ier) |
---|
322 | ier = NF90_INQ_VARID (force_id,'neighbours',v_id) |
---|
323 | ier = NF90_GET_VAR (force_id,v_id,x_neighbours) |
---|
324 | neighbours(:,:) = NINT(x_neighbours(:,:)) |
---|
325 | DEALLOCATE (x_neighbours) |
---|
326 | !- |
---|
327 | ier = NF90_INQ_VARID (force_id,'resolution',v_id) |
---|
328 | ier = NF90_GET_VAR (force_id,v_id,resolution) |
---|
329 | !- |
---|
330 | ier = NF90_INQ_VARID (force_id,'contfrac',v_id) |
---|
331 | ier = NF90_GET_VAR (force_id,v_id,contfrac) |
---|
332 | !- |
---|
333 | ! activate CO2, STOMATE, but not sechiba |
---|
334 | !- |
---|
335 | control%river_routing = .FALSE. |
---|
336 | control%hydrol_cwrr = .FALSE. |
---|
337 | control%ok_sechiba = .FALSE. |
---|
338 | ! |
---|
339 | control%stomate_watchout = .TRUE. |
---|
340 | control%ok_co2 = .TRUE. |
---|
341 | control%ok_stomate = .TRUE. |
---|
342 | !- |
---|
343 | ! is DGVM activated? |
---|
344 | !- |
---|
345 | control%ok_dgvm = .FALSE. |
---|
346 | CALL getin('STOMATE_OK_DGVM',control%ok_dgvm) |
---|
347 | WRITE(*,*) 'LPJ is activated: ',control%ok_dgvm |
---|
348 | !- |
---|
349 | ! restart files |
---|
350 | !- |
---|
351 | ! Sechiba's restart files |
---|
352 | sec_restname_in = 'sechiba_start.nc' |
---|
353 | CALL getin('SECHIBA_restart_in',sec_restname_in) |
---|
354 | WRITE(*,*) 'SECHIBA INPUT RESTART_FILE: ',TRIM(sec_restname_in) |
---|
355 | IF ( TRIM(sec_restname_in) .EQ. 'NONE' ) THEN |
---|
356 | STOP 'Need a restart file for Sechiba' |
---|
357 | ENDIF |
---|
358 | sec_restname_out = 'sechiba_restart.nc' |
---|
359 | CALL getin('SECHIBA_rest_out',sec_restname_out) |
---|
360 | WRITE(*,*) 'SECHIBA OUTPUT RESTART_FILE: ',TRIM(sec_restname_out) |
---|
361 | ! Stomate's restart files |
---|
362 | sto_restname_in = 'stomate_start.nc' |
---|
363 | CALL getin ('STOMATE_RESTART_FILEIN',sto_restname_in) |
---|
364 | WRITE(*,*) 'STOMATE INPUT RESTART_FILE: ',TRIM(sto_restname_in) |
---|
365 | sto_restname_out = 'stomate_restart.nc' |
---|
366 | CALL getin ('STOMATE_RESTART_FILEOUT',sto_restname_out) |
---|
367 | WRITE(*,*) 'STOMATE OUTPUT RESTART_FILE: ',TRIM(sto_restname_out) |
---|
368 | !- |
---|
369 | ! We need to know iim, jjm. |
---|
370 | ! Get them from the restart files themselves. |
---|
371 | !- |
---|
372 | iret = NF90_OPEN (sec_restname_in,NF90_NOWRITE,ncfid) |
---|
373 | IF (iret /= NF90_NOERR) THEN |
---|
374 | CALL ipslerr (3,'teststomate', & |
---|
375 | & 'Could not open file : ', & |
---|
376 | & sec_restname_in,'(Do you have forget it ?)') |
---|
377 | ENDIF |
---|
378 | iret = NF90_INQUIRE_DIMENSION (ncfid,1,len=iim) |
---|
379 | iret = NF90_INQUIRE_DIMENSION (ncfid,2,len=jjm) |
---|
380 | iret = NF90_CLOSE (ncfid) |
---|
381 | ! Allocate longitudes and latitudes |
---|
382 | ALLOCATE (lon(iim,jjm),stat=ier) |
---|
383 | a_er = a_er.OR.(ier.NE.0) |
---|
384 | ALLOCATE (lat(iim,jjm),stat=ier) |
---|
385 | a_er = a_er.OR.(ier.NE.0) |
---|
386 | lon(:,:) = 0.0 |
---|
387 | lat(:,:) = 0.0 |
---|
388 | lev(1) = 0.0 |
---|
389 | !- |
---|
390 | CALL restini & |
---|
391 | & (sec_restname_in, iim, jjm, lon, lat, llm, lev, & |
---|
392 | & sec_restname_out, itau_dep, date0, dt_files, rest_id_sec) |
---|
393 | !- |
---|
394 | CALL restini & |
---|
395 | & (sto_restname_in, iim, jjm, lon, lat, llm, lev, & |
---|
396 | & sto_restname_out, itau_dep, date0, dt_files, rest_id_sto) |
---|
397 | !- |
---|
398 | IF ( dt_files .NE. dtradia ) THEN |
---|
399 | WRITE(*,*) 'dt_files',dt_files |
---|
400 | WRITE(*,*) 'dtradia',dtradia |
---|
401 | STOP 'PROBLEM with time steps.' |
---|
402 | ENDIF |
---|
403 | !- |
---|
404 | ! a time step for STOMATE corresponds to itau_step timesteps in SECHIBA |
---|
405 | !- |
---|
406 | itau_step = NINT(dt_force/dtradia) |
---|
407 | ! |
---|
408 | CALL ioconf_startdate(date0) |
---|
409 | CALL intsurf_time( itau_dep, date0, dtradia ) |
---|
410 | !- |
---|
411 | ! Length of integration |
---|
412 | !- |
---|
413 | WRITE(time_str,'(a)') '1Y' |
---|
414 | CALL getin ('TIME_LENGTH', time_str) |
---|
415 | ! transform into itau |
---|
416 | CALL tlen2itau(time_str, dt_files, date0, itau_len) |
---|
417 | ! itau_len*dtradia must be a multiple of dt_force |
---|
418 | itau_len = NINT( MAX(1.,FLOAT(NINT(itau_len*dtradia/dt_force))) & |
---|
419 | & *dt_force/dtradia) |
---|
420 | !- |
---|
421 | ! set up STOMATE history file |
---|
422 | !- |
---|
423 | !Config Key = STOMATE_OUTPUT_FILE |
---|
424 | !Config Desc = Name of file in which STOMATE's output is going |
---|
425 | !Config to be written |
---|
426 | !Config Def = stomate_history.nc |
---|
427 | !Config Help = This file is going to be created by the model |
---|
428 | !Config and will contain the output from the model. |
---|
429 | !Config This file is a truly COADS compliant netCDF file. |
---|
430 | !Config It will be generated by the hist software from |
---|
431 | !Config the IOIPSL package. |
---|
432 | !- |
---|
433 | stom_histname='stomate_history.nc' |
---|
434 | CALL getin ('STOMATE_OUTPUT_FILE', stom_histname) |
---|
435 | WRITE(*,*) 'STOMATE_OUTPUT_FILE', TRIM(stom_histname) |
---|
436 | !- |
---|
437 | !Config Key = STOMATE_HIST_DT |
---|
438 | !Config Desc = STOMATE history time step (d) |
---|
439 | !Config Def = 10. |
---|
440 | !Config Help = Time step of the STOMATE history file |
---|
441 | !- |
---|
442 | hist_days_stom = 10. |
---|
443 | CALL getin ('STOMATE_HIST_DT', hist_days_stom) |
---|
444 | hist_dt_stom = NINT( hist_days_stom ) * one_day |
---|
445 | WRITE(*,*) 'output frequency for STOMATE history file (d): ', & |
---|
446 | hist_dt_stom/one_day |
---|
447 | !- |
---|
448 | ! initialize |
---|
449 | CALL histbeg(stom_histname, iim, lon, jjm, lat, 1, iim, 1, jjm, & |
---|
450 | & itau_dep, date0, dt_files, hori_id, hist_id_sto) |
---|
451 | ! define PFT axis |
---|
452 | hist_PFTaxis = (/ ( REAL(i,r_std), i=1,nvm ) /) |
---|
453 | ! declare this axis |
---|
454 | CALL histvert (hist_id_sto, 'PFT', 'Plant functional type', & |
---|
455 | & '-', nvm, hist_PFTaxis, hist_PFTaxis_id) |
---|
456 | !- define Pool_10 axis |
---|
457 | hist_pool_10axis = (/ ( REAL(i,r_std), i=1,10 ) /) |
---|
458 | !- declare this axis |
---|
459 | CALL histvert (hist_id_sto, 'P10', 'Pool 10 years', & |
---|
460 | & '-', 10, hist_pool_10axis, hist_pool_10axis_id) |
---|
461 | !- define Pool_100 axis |
---|
462 | hist_pool_100axis = (/ ( REAL(i,r_std), i=1,100 ) /) |
---|
463 | !- declare this axis |
---|
464 | CALL histvert (hist_id_sto, 'P100', 'Pool 100 years', & |
---|
465 | & '-', 100, hist_pool_100axis, hist_pool_100axis_id) |
---|
466 | !- define Pool_11 axis |
---|
467 | hist_pool_11axis = (/ ( REAL(i,r_std), i=1,11 ) /) |
---|
468 | !- declare this axis |
---|
469 | CALL histvert (hist_id_sto, 'P11', 'Pool 10 years + 1', & |
---|
470 | & '-', 11, hist_pool_11axis, hist_pool_11axis_id) |
---|
471 | !- define Pool_101 axis |
---|
472 | hist_pool_101axis = (/ ( REAL(i,r_std), i=1,101 ) /) |
---|
473 | !- declare this axis |
---|
474 | CALL histvert (hist_id_sto, 'P101', 'Pool 100 years + 1', & |
---|
475 | & '-', 101, hist_pool_101axis, hist_pool_101axis_id) |
---|
476 | ! define STOMATE history file |
---|
477 | CALL stom_define_history (hist_id_sto, nvm, iim, jjm, & |
---|
478 | & dt_files, hist_dt_stom, hori_id, hist_PFTaxis_id, & |
---|
479 | & hist_pool_10axis_id, hist_pool_100axis_id, & |
---|
480 | & hist_pool_11axis_id, hist_pool_101axis_id) |
---|
481 | ! end definition |
---|
482 | CALL histend(hist_id_sto) |
---|
483 | hist_id_sec = -1 |
---|
484 | hist_id_sec2 = -1 |
---|
485 | hist_id_stom_IPCC = -1 |
---|
486 | !- |
---|
487 | ! read some variables we need from SECHIBA's restart file |
---|
488 | !- |
---|
489 | CALL ioget_expval(val_exp) |
---|
490 | !- |
---|
491 | ! first call of slowproc to initialize variables |
---|
492 | !- |
---|
493 | itau = itau_dep |
---|
494 | ldrestart_read = .TRUE. |
---|
495 | ldrestart_write = .FALSE. |
---|
496 | !- |
---|
497 | !MM Problem here with dpu which depends on soil type |
---|
498 | DO jv = 1, nbdl-1 |
---|
499 | ! first 2.0 is dpu |
---|
500 | ! second 2.0 is average |
---|
501 | diaglev(jv) = 2.0/(2**(nbdl-1) -1) * ( ( 2**(jv-1) -1) + ( 2**(jv) -1) ) / 2.0 |
---|
502 | ENDDO |
---|
503 | diaglev(nbdl) = 2.0 |
---|
504 | !- |
---|
505 | ! For sequential use only, we must initialize data_para : |
---|
506 | ! |
---|
507 | CALL init_para(.FALSE.) |
---|
508 | ! |
---|
509 | CALL init_data_para(iim,jjm,kjpindex,indices) |
---|
510 | ! |
---|
511 | !- global index index_g is the index_l of root proc |
---|
512 | IF (is_root_prc) index_g(:)=indices(1:kjpindex) |
---|
513 | !- |
---|
514 | CALL slowproc_main & |
---|
515 | & (itau, kjpij, kjpindex, dt_force, date0, & |
---|
516 | & ldrestart_read, ldrestart_write, .FALSE., .TRUE., & |
---|
517 | & indices, indexveg, lalo, neighbours, resolution, contfrac, soiltype, & |
---|
518 | & t2m, t2m_min, temp_sol, soiltemp, & |
---|
519 | & humrel_x, soilhum, litterhum, precip_rain, precip_snow, gpp_x, & |
---|
520 | & deadleaf_cover, assim_param_x, lai_x, height_x, veget_x, & |
---|
521 | & frac_nobio, veget_max_x, totfrac_nobio, qsintmax_x, & |
---|
522 | & rest_id_sec, hist_id_sec, hist_id_sec2, rest_id_sto, hist_id_sto, hist_id_stom_IPCC, co2_flux) |
---|
523 | ! correct date |
---|
524 | day_counter = one_day - dt_force |
---|
525 | date=1 |
---|
526 | !- |
---|
527 | ! time loop |
---|
528 | !- |
---|
529 | DO itau = itau_dep+itau_step,itau_dep+itau_len,itau_step |
---|
530 | !-- next forcing state |
---|
531 | iisf = iisf+1 |
---|
532 | !--- |
---|
533 | IF (iisf .GT. nsfm) THEN |
---|
534 | !----- |
---|
535 | !---- we have to read new forcing states from the file |
---|
536 | !----- |
---|
537 | !---- determine blocks of forcing states that are contiguous in memory |
---|
538 | !----- |
---|
539 | nblocks = 0 |
---|
540 | ifirst(:) = 1 |
---|
541 | ilast(:) = 1 |
---|
542 | !----- |
---|
543 | DO iisf=1,nsfm |
---|
544 | IF ( (nblocks .NE. 0) ) THEN |
---|
545 | IF ( (isf(iisf) .EQ. isf(ilast(nblocks))+1) ) THEN |
---|
546 | !-------- element is contiguous with last element found |
---|
547 | ilast(nblocks) = iisf |
---|
548 | ELSE |
---|
549 | !-------- found first element of new block |
---|
550 | nblocks = nblocks+1 |
---|
551 | IF (nblocks .GT. nsfm) THEN |
---|
552 | ! IF (nblocks .GT. 2) THEN |
---|
553 | STOP 'Problem in teststomate' |
---|
554 | ENDIF |
---|
555 | ifirst(nblocks) = iisf |
---|
556 | ilast(nblocks) = iisf |
---|
557 | ENDIF |
---|
558 | ELSE |
---|
559 | !-------- found first element of new block |
---|
560 | nblocks = nblocks+1 |
---|
561 | IF (nblocks .GT. nsfm) THEN |
---|
562 | ! IF (nblocks .GT. 2) THEN |
---|
563 | STOP 'Problem in teststomate' |
---|
564 | ENDIF |
---|
565 | ifirst(nblocks) = iisf |
---|
566 | ilast(nblocks) = iisf |
---|
567 | ENDIF |
---|
568 | ENDDO |
---|
569 | !----- |
---|
570 | DO iblocks=1,nblocks |
---|
571 | IF (ifirst(iblocks) .NE. ilast(iblocks)) THEN |
---|
572 | ndim = 2; |
---|
573 | start(:) = 1; start(ndim) = isf(ifirst(iblocks)); |
---|
574 | count(1:ndim) = SHAPE(clay_fm) |
---|
575 | count(ndim) = isf(ilast(iblocks))-isf(ifirst(iblocks))+1 |
---|
576 | ier = NF90_INQ_VARID (force_id,'clay',v_id) |
---|
577 | a_er = a_er.OR.(ier.NE.0) |
---|
578 | ier = NF90_GET_VAR (force_id,v_id, & |
---|
579 | & clay_fm(:,ifirst(iblocks):ilast(iblocks)), & |
---|
580 | & start=start(1:ndim), count=count(1:ndim)) |
---|
581 | a_er = a_er.OR.(ier.NE.0) |
---|
582 | !--------- |
---|
583 | ndim = 3; |
---|
584 | start(:) = 1; start(ndim) = isf(ifirst(iblocks)); |
---|
585 | count(1:ndim) = SHAPE(humrel_x_fm) |
---|
586 | count(ndim) = isf(ilast(iblocks))-isf(ifirst(iblocks))+1 |
---|
587 | ier = NF90_INQ_VARID (force_id,'humrel',v_id) |
---|
588 | a_er = a_er.OR.(ier.NE.0) |
---|
589 | ier = NF90_GET_VAR (force_id,v_id, & |
---|
590 | & humrel_x_fm(:,:,ifirst(iblocks):ilast(iblocks)), & |
---|
591 | & start=start(1:ndim), count=count(1:ndim)) |
---|
592 | a_er = a_er.OR.(ier.NE.0) |
---|
593 | !--------- |
---|
594 | ndim = 2; |
---|
595 | start(:) = 1; start(ndim) = isf(ifirst(iblocks)); |
---|
596 | count(1:ndim) = SHAPE(litterhum_fm) |
---|
597 | count(ndim) = isf(ilast(iblocks))-isf(ifirst(iblocks))+1 |
---|
598 | ier = NF90_INQ_VARID (force_id,'litterhum',v_id) |
---|
599 | a_er = a_er.OR.(ier.NE.0) |
---|
600 | ier = NF90_GET_VAR (force_id,v_id, & |
---|
601 | & litterhum_fm(:,ifirst(iblocks):ilast(iblocks)), & |
---|
602 | & start=start(1:ndim), count=count(1:ndim)) |
---|
603 | a_er = a_er.OR.(ier.NE.0) |
---|
604 | !--------- |
---|
605 | ndim = 2; |
---|
606 | start(:) = 1; start(ndim) = isf(ifirst(iblocks)); |
---|
607 | count(1:ndim) = SHAPE(t2m_fm) |
---|
608 | count(ndim) = isf(ilast(iblocks))-isf(ifirst(iblocks))+1 |
---|
609 | ier = NF90_INQ_VARID (force_id,'t2m',v_id) |
---|
610 | a_er = a_er.OR.(ier.NE.0) |
---|
611 | ier = NF90_GET_VAR (force_id,v_id, & |
---|
612 | & t2m_fm(:,ifirst(iblocks):ilast(iblocks)), & |
---|
613 | & start=start(1:ndim), count=count(1:ndim)) |
---|
614 | a_er = a_er.OR.(ier.NE.0) |
---|
615 | !--------- |
---|
616 | ndim = 2; |
---|
617 | start(:) = 1; start(ndim) = isf(ifirst(iblocks)); |
---|
618 | count(1:ndim) = SHAPE(t2m_min_fm) |
---|
619 | count(ndim) = isf(ilast(iblocks))-isf(ifirst(iblocks))+1 |
---|
620 | ier = NF90_INQ_VARID (force_id,'t2m_min',v_id) |
---|
621 | a_er = a_er.OR.(ier.NE.0) |
---|
622 | ier = NF90_GET_VAR (force_id,v_id, & |
---|
623 | & t2m_min_fm(:,ifirst(iblocks):ilast(iblocks)), & |
---|
624 | & start=start(1:ndim), count=count(1:ndim)) |
---|
625 | a_er = a_er.OR.(ier.NE.0) |
---|
626 | !--------- |
---|
627 | ndim = 2; |
---|
628 | start(:) = 1; start(ndim) = isf(ifirst(iblocks)); |
---|
629 | count(1:ndim) = SHAPE(temp_sol_fm) |
---|
630 | count(ndim) = isf(ilast(iblocks))-isf(ifirst(iblocks))+1 |
---|
631 | ier = NF90_INQ_VARID (force_id,'tsurf',v_id) |
---|
632 | a_er = a_er.OR.(ier.NE.0) |
---|
633 | ier = NF90_GET_VAR (force_id,v_id, & |
---|
634 | & temp_sol_fm(:,ifirst(iblocks):ilast(iblocks)), & |
---|
635 | & start=start(1:ndim), count=count(1:ndim)) |
---|
636 | a_er = a_er.OR.(ier.NE.0) |
---|
637 | !--------- |
---|
638 | ndim = 3; |
---|
639 | start(:) = 1; start(ndim) = isf(ifirst(iblocks)); |
---|
640 | count(1:ndim) = SHAPE(soiltemp_fm) |
---|
641 | count(ndim) = isf(ilast(iblocks))-isf(ifirst(iblocks))+1 |
---|
642 | ier = NF90_INQ_VARID (force_id,'tsoil',v_id) |
---|
643 | a_er = a_er.OR.(ier.NE.0) |
---|
644 | ier = NF90_GET_VAR (force_id,v_id, & |
---|
645 | & soiltemp_fm(:,:,ifirst(iblocks):ilast(iblocks)), & |
---|
646 | & start=start(1:ndim), count=count(1:ndim)) |
---|
647 | a_er = a_er.OR.(ier.NE.0) |
---|
648 | !--------- |
---|
649 | ndim = 3; |
---|
650 | start(:) = 1; start(ndim) = isf(ifirst(iblocks)); |
---|
651 | count(1:ndim) = SHAPE(soilhum_fm) |
---|
652 | count(ndim) = isf(ilast(iblocks))-isf(ifirst(iblocks))+1 |
---|
653 | ier = NF90_INQ_VARID (force_id,'soilhum',v_id) |
---|
654 | a_er = a_er.OR.(ier.NE.0) |
---|
655 | ier = NF90_GET_VAR (force_id,v_id, & |
---|
656 | & soilhum_fm(:,:,ifirst(iblocks):ilast(iblocks)), & |
---|
657 | & start=start(1:ndim), count=count(1:ndim)) |
---|
658 | a_er = a_er.OR.(ier.NE.0) |
---|
659 | !--------- |
---|
660 | ndim = 2; |
---|
661 | start(:) = 1; start(ndim) = isf(ifirst(iblocks)); |
---|
662 | count(1:ndim) = SHAPE(precip_fm) |
---|
663 | count(ndim) = isf(ilast(iblocks))-isf(ifirst(iblocks))+1 |
---|
664 | ier = NF90_INQ_VARID (force_id,'precip',v_id) |
---|
665 | a_er = a_er.OR.(ier.NE.0) |
---|
666 | ier = NF90_GET_VAR (force_id,v_id, & |
---|
667 | & precip_fm(:,ifirst(iblocks):ilast(iblocks)), & |
---|
668 | & start=start(1:ndim), count=count(1:ndim)) |
---|
669 | a_er = a_er.OR.(ier.NE.0) |
---|
670 | !--------- |
---|
671 | ndim = 3; |
---|
672 | start(:) = 1; start(ndim) = isf(ifirst(iblocks)); |
---|
673 | count(1:ndim) = SHAPE(gpp_x_fm) |
---|
674 | count(ndim) = isf(ilast(iblocks))-isf(ifirst(iblocks))+1 |
---|
675 | ier = NF90_INQ_VARID (force_id,'gpp',v_id) |
---|
676 | a_er = a_er.OR.(ier.NE.0) |
---|
677 | ier = NF90_GET_VAR (force_id,v_id, & |
---|
678 | & gpp_x_fm(:,:,ifirst(iblocks):ilast(iblocks)), & |
---|
679 | & start=start(1:ndim), count=count(1:ndim)) |
---|
680 | a_er = a_er.OR.(ier.NE.0) |
---|
681 | !--------- |
---|
682 | ndim = 3; |
---|
683 | start(:) = 1; start(ndim) = isf(ifirst(iblocks)); |
---|
684 | count(1:ndim) = SHAPE(veget_force_x_fm) |
---|
685 | count(ndim) = isf(ilast(iblocks))-isf(ifirst(iblocks))+1 |
---|
686 | ier = NF90_INQ_VARID (force_id,'veget',v_id) |
---|
687 | a_er = a_er.OR.(ier.NE.0) |
---|
688 | ier = NF90_GET_VAR (force_id,v_id, & |
---|
689 | & veget_force_x_fm(:,:,ifirst(iblocks):ilast(iblocks)), & |
---|
690 | & start=start(1:ndim), count=count(1:ndim)) |
---|
691 | a_er = a_er.OR.(ier.NE.0) |
---|
692 | !--------- |
---|
693 | ndim = 3; |
---|
694 | start(:) = 1; start(ndim) = isf(ifirst(iblocks)); |
---|
695 | count(1:ndim) = SHAPE(veget_max_force_x_fm) |
---|
696 | count(ndim) = isf(ilast(iblocks))-isf(ifirst(iblocks))+1 |
---|
697 | ier = NF90_INQ_VARID (force_id,'veget_max',v_id) |
---|
698 | a_er = a_er.OR.(ier.NE.0) |
---|
699 | ier = NF90_GET_VAR (force_id,v_id, & |
---|
700 | & veget_max_force_x_fm(:,:,ifirst(iblocks):ilast(iblocks)), & |
---|
701 | & start=start(1:ndim), count=count(1:ndim)) |
---|
702 | a_er = a_er.OR.(ier.NE.0) |
---|
703 | !--------- |
---|
704 | ndim = 3; |
---|
705 | start(:) = 1; start(ndim) = isf(ifirst(iblocks)); |
---|
706 | count(1:ndim) = SHAPE(lai_force_x_fm) |
---|
707 | count(ndim) = isf(ilast(iblocks))-isf(ifirst(iblocks))+1 |
---|
708 | ier = NF90_INQ_VARID (force_id,'lai',v_id) |
---|
709 | a_er = a_er.OR.(ier.NE.0) |
---|
710 | ier = NF90_GET_VAR (force_id,v_id, & |
---|
711 | & lai_force_x_fm(:,:,ifirst(iblocks):ilast(iblocks)), & |
---|
712 | & start=start(1:ndim), count=count(1:ndim)) |
---|
713 | a_er = a_er.OR.(ier.NE.0) |
---|
714 | ENDIF |
---|
715 | ENDDO |
---|
716 | !----- |
---|
717 | !---- determine which forcing states must be read next time |
---|
718 | !----- |
---|
719 | isf(1) = isf(nsfm)+1 |
---|
720 | IF ( isf(1) .GT. nsft ) isf(1) = 1 |
---|
721 | DO iisf = 2, nsfm |
---|
722 | isf(iisf) = isf(iisf-1)+1 |
---|
723 | IF ( isf(iisf) .GT. nsft ) isf(iisf) = 1 |
---|
724 | ENDDO |
---|
725 | !---- start again at first forcing state |
---|
726 | iisf = 1 |
---|
727 | ENDIF |
---|
728 | soiltype(:,3) = clay_fm(:,iisf) |
---|
729 | humrel_x(:,:) = humrel_x_fm(:,:,iisf) |
---|
730 | litterhum(:) = litterhum_fm(:,iisf) |
---|
731 | t2m(:) = t2m_fm(:,iisf) |
---|
732 | t2m_min(:) = t2m_min_fm(:,iisf) |
---|
733 | temp_sol(:) = temp_sol_fm(:,iisf) |
---|
734 | soiltemp(:,:) = soiltemp_fm(:,:,iisf) |
---|
735 | soilhum(:,:) = soilhum_fm(:,:,iisf) |
---|
736 | precip_rain(:) = precip_fm(:,iisf) |
---|
737 | gpp_x(:,:) = gpp_x_fm(:,:,iisf) |
---|
738 | veget_force_x(:,:) = veget_force_x_fm(:,:,iisf) |
---|
739 | veget_max_force_x(:,:) = veget_max_force_x_fm(:,:,iisf) |
---|
740 | lai_force_x(:,:) = lai_force_x_fm(:,:,iisf) |
---|
741 | WHERE ( t2m(:) .LT. ZeroCelsius ) |
---|
742 | precip_snow(:) = precip_rain(:) |
---|
743 | precip_rain(:) = 0. |
---|
744 | ELSEWHERE |
---|
745 | precip_snow(:) = 0. |
---|
746 | ENDWHERE |
---|
747 | !--- |
---|
748 | !-- scale GPP to new lai and veget_max |
---|
749 | !--- |
---|
750 | WHERE ( lai_x(:,:) .EQ. 0.0 ) gpp_x(:,:) = 0.0 |
---|
751 | !-- scale GPP to new LAI |
---|
752 | WHERE (lai_force_x(:,:) .GT. 0.0 ) |
---|
753 | gpp_x(:,:) = gpp_x(:,:)*atan(2.*lai_x(:,:)) & |
---|
754 | & /atan( 2.*MAX(lai_force_x(:,:),0.01)) |
---|
755 | ENDWHERE |
---|
756 | !-- scale GPP to new veget_max |
---|
757 | WHERE (veget_max_force_x(:,:) .GT. 0.0 ) |
---|
758 | gpp_x(:,:) = gpp_x(:,:)*veget_max_x(:,:)/veget_max_force_x(:,:) |
---|
759 | ENDWHERE |
---|
760 | !--- |
---|
761 | !-- number crunching |
---|
762 | !--- |
---|
763 | CALL intsurf_time( itau, date0, dtradia ) |
---|
764 | ldrestart_read = .FALSE. |
---|
765 | ldrestart_write = .FALSE. |
---|
766 | CALL slowproc_main & |
---|
767 | & (itau, kjpij, kjpindex, dt_force, date0, & |
---|
768 | & ldrestart_read, ldrestart_write, .FALSE., .TRUE., & |
---|
769 | & indices, indexveg, lalo, neighbours, resolution, contfrac, soiltype, & |
---|
770 | & t2m, t2m_min, temp_sol, soiltemp, & |
---|
771 | & humrel_x, soilhum, litterhum, precip_rain, precip_snow, gpp_x, & |
---|
772 | & deadleaf_cover, assim_param_x, lai_x, height_x, veget_x, & |
---|
773 | & frac_nobio, veget_max_x, totfrac_nobio, qsintmax_x, & |
---|
774 | & rest_id_sec, hist_id_sec, hist_id_sec2, rest_id_sto, hist_id_sto, hist_id_stom_IPCC, co2_flux) |
---|
775 | day_counter = one_day - dt_force |
---|
776 | ENDDO ! end of the time loop |
---|
777 | !- |
---|
778 | itau = itau -itau_step |
---|
779 | !- |
---|
780 | ! write restart files |
---|
781 | !- |
---|
782 | ! first, read and write variables that are not managed otherwise |
---|
783 | taboo_vars = & |
---|
784 | & '$lat$ $lon$ $lev$ $veget_year$ '// & |
---|
785 | & '$height$ $day_counter$ $veget$ $veget_max$ $frac_nobio$ '// & |
---|
786 | & '$lai$ $soiltype_frac$ $clay_frac$ '// & |
---|
787 | & '$nav_lon$ $nav_lat$ $nav_lev$ $time$ $time_steps$' |
---|
788 | !- |
---|
789 | CALL ioget_vname(rest_id_sec, nbvar, varnames) |
---|
790 | !!$!- |
---|
791 | !!$! read and write some special variables (1D or variables that we need) |
---|
792 | !!$!- |
---|
793 | !!$ var_name = 'day_counter' |
---|
794 | !!$ CALL restget (rest_id_sto, var_name, 1, 1, 1, itau_dep, .TRUE., xtmp) |
---|
795 | !!$ CALL restput (rest_id_sto, var_name, 1, 1, 1, itau_dep, xtmp) |
---|
796 | !!$!- |
---|
797 | !!$ var_name = 'dt_days' |
---|
798 | !!$ CALL restget (rest_id_sto, var_name, 1, 1, 1, itau_dep, .TRUE., xtmp) |
---|
799 | !!$ CALL restput (rest_id_sto, var_name, 1, 1, 1, itau_dep, xtmp) |
---|
800 | !!$!- |
---|
801 | !!$ var_name = 'date' |
---|
802 | !!$ CALL restget (rest_id_sto, var_name, 1, 1, 1, itau_dep, .TRUE., xtmp) |
---|
803 | !!$ CALL restput (rest_id_sto, var_name, 1, 1, 1, itau_dep, xtmp) |
---|
804 | !- |
---|
805 | DO iv = 1, nbvar |
---|
806 | !-- check if the variable is to be written here |
---|
807 | IF (INDEX( taboo_vars,'$'//TRIM(varnames(iv))//'$') .EQ. 0 ) THEN |
---|
808 | !---- get variable dimensions, especially 3rd dimension |
---|
809 | CALL ioget_vdim & |
---|
810 | & (rest_id_sec, varnames(iv), varnbdim_max, varnbdim, vardims) |
---|
811 | l1d = ALL(vardims(1:varnbdim) .EQ. 1) |
---|
812 | ALLOCATE(var_3d(kjpindex,vardims(3)),stat=ier) |
---|
813 | IF (ier .NE. 0) STOP 'ALLOCATION PROBLEM' |
---|
814 | !---- read it |
---|
815 | IF (l1d) THEN |
---|
816 | CALL restget (rest_id_sec,TRIM(varnames(iv)), & |
---|
817 | 1,vardims(3),1,itau_dep,.TRUE.,var_3d) |
---|
818 | ELSE |
---|
819 | CALL restget (rest_id_sec,TRIM(varnames(iv)), & |
---|
820 | kjpindex,vardims(3),1,itau_dep,.TRUE.,var_3d, & |
---|
821 | "gather",kjpindex,indices) |
---|
822 | ENDIF |
---|
823 | !---- write it |
---|
824 | IF (l1d) THEN |
---|
825 | CALL restput (rest_id_sec,TRIM(varnames(iv)), & |
---|
826 | 1,vardims(3),1,itau,var_3d) |
---|
827 | ELSE |
---|
828 | CALL restput (rest_id_sec,TRIM(varnames(iv)), & |
---|
829 | kjpindex,vardims(3),1,itau,var_3d, & |
---|
830 | 'scatter',kjpindex,indices) |
---|
831 | ENDIF |
---|
832 | DEALLOCATE (var_3d) |
---|
833 | ENDIF |
---|
834 | ENDDO |
---|
835 | ! call slowproc to write restart files |
---|
836 | ldrestart_read = .FALSE. |
---|
837 | ldrestart_write = .TRUE. |
---|
838 | !- |
---|
839 | CALL slowproc_main & |
---|
840 | & (itau, kjpij, kjpindex, dt_force, date0, & |
---|
841 | & ldrestart_read, ldrestart_write, .FALSE., .TRUE., & |
---|
842 | & indices, indexveg, lalo, neighbours, resolution, contfrac, soiltype, & |
---|
843 | & t2m, t2m_min, temp_sol, soiltemp, & |
---|
844 | & humrel_x, soilhum, litterhum, precip_rain, precip_snow, gpp_x, & |
---|
845 | & deadleaf_cover, assim_param_x, lai_x, height_x, veget_x, & |
---|
846 | & frac_nobio, veget_max_x, totfrac_nobio, qsintmax_x, & |
---|
847 | & rest_id_sec, hist_id_sec, hist_id_sec2, rest_id_sto, hist_id_sto, hist_id_stom_IPCC, co2_flux) |
---|
848 | !- |
---|
849 | ! close files |
---|
850 | !- |
---|
851 | CALL restclo |
---|
852 | CALL histclo |
---|
853 | ier = NF90_CLOSE (force_id) |
---|
854 | !- |
---|
855 | ! write a new driver restart file with correct time step |
---|
856 | !- |
---|
857 | write(*,*) 'teststomate: writing driver restart file with correct time step.' |
---|
858 | dri_restname_in = 'driver_start.nc' |
---|
859 | CALL getin ('RESTART_FILEIN',dri_restname_in) |
---|
860 | dri_restname_out = 'driver_restart.nc' |
---|
861 | CALL getin ('RESTART_FILEOUT',dri_restname_out) |
---|
862 | CALL SYSTEM & |
---|
863 | & ('cp '//TRIM(dri_restname_in)//' '//TRIM(dri_restname_out)) |
---|
864 | !- |
---|
865 | iret = NF90_OPEN (TRIM(sec_restname_out),NF90_NOWRITE,ncfid) |
---|
866 | iret = NF90_INQ_VARID (ncfid,'time',v_id) |
---|
867 | iret = NF90_GET_VAR (ncfid,v_id,r1d) |
---|
868 | time_sec = r1d(1) |
---|
869 | iret = NF90_INQ_VARID (ncfid,'time_steps',v_id) |
---|
870 | iret = NF90_GET_VAR (ncfid,v_id,time_step_sec) |
---|
871 | iret = NF90_CLOSE (ncfid) |
---|
872 | !- |
---|
873 | iret = NF90_OPEN (TRIM(dri_restname_out),NF90_WRITE,ncfid) |
---|
874 | iret = NF90_INQ_VARID (ncfid,'time',v_id) |
---|
875 | iret = NF90_GET_VAR (ncfid,v_id,r1d) |
---|
876 | time_dri = r1d(1) |
---|
877 | r1d(1) = time_sec |
---|
878 | iret = NF90_PUT_VAR (ncfid,v_id,r1d) |
---|
879 | iret = NF90_INQ_VARID (ncfid,'time_steps',v_id) |
---|
880 | iret = NF90_GET_VAR (ncfid,v_id,time_step_dri) |
---|
881 | iret = NF90_PUT_VAR (ncfid,v_id,time_step_sec) |
---|
882 | iret = NF90_INQ_VARID (ncfid,'julian',v_id) |
---|
883 | iret = NF90_GET_VAR (ncfid,v_id,r1d) |
---|
884 | julian = r1d(1) |
---|
885 | djulian = (time_step_sec-time_step_dri)*dtradia/one_day |
---|
886 | julian = julian & |
---|
887 | & +djulian-FLOAT(INT((julian+djulian)/one_year))*one_year |
---|
888 | r1d(1) = julian |
---|
889 | iret = NF90_PUT_VAR (ncfid,v_id,r1d) |
---|
890 | iret = NF90_CLOSE (ncfid) |
---|
891 | |
---|
892 | CALL getin_dump |
---|
893 | !--------------- |
---|
894 | END PROGRAM teststomate |
---|
895 | ! |
---|
896 | !=== |
---|
897 | ! |
---|