source: branches/ORCHIDEE_EXT/ORCHIDEE_OL/forcesoil.f90 @ 313

Last change on this file since 313 was 313, checked in by didier.solyga, 13 years ago

Update teststomate.f90 and forcesoil.f90 in order to take into account the modifications of the revisions 311 and 312

File size: 14.4 KB
Line 
1PROGRAM forcesoil
2  !---------------------------------------------------------------------
3  ! This program reads a forcing file for STOMATE's soil carbon routine
4  ! which was created with STOMATE.
5  ! It then integrates STOMATE's soil carbon parameterizations
6  ! for a given number of years using this forcing.
7  !---------------------------------------------------------------------
8  !- IPSL (2006)
9  !-  This software is governed by the CeCILL licence see ORCHIDEE/ORCHIDEE_CeCILL.LIC
10  USE netcdf
11  !-
12  USE defprec
13  USE constantes
14  USE pft_parameters 
15  USE stomate_data
16  USE ioipsl
17  USE stomate_soilcarbon
18  USE parallel
19  !-
20  IMPLICIT NONE
21  !-
22  CHARACTER(LEN=80) :: sto_restname_in,sto_restname_out,var_name
23  INTEGER(i_std)                             :: iim,jjm
24  INTEGER(i_std),PARAMETER                   :: llm = 1
25  INTEGER(i_std)                             :: kjpindex
26  INTEGER(i_std)                             :: itau_dep,itau_len
27  CHARACTER(LEN=30)                         :: time_str
28  INTEGER(i_std)                             :: ier,iret
29  REAL(r_std)                                :: dt_files
30  REAL(r_std)                                :: date0
31  INTEGER(i_std)                             :: rest_id_sto
32  INTEGER(i_std)                             :: ncfid
33  REAL(r_std)                                :: dt_force,dt_forcesoil
34  INTEGER                                   :: nparan
35  INTEGER,PARAMETER                         :: nparanmax=36
36  REAL(r_std)                                :: xbid1,xbid2
37  INTEGER(i_std)                             :: ibid
38  INTEGER(i_std),DIMENSION(:),ALLOCATABLE    :: indices
39  REAL(r_std),DIMENSION(llm)                 :: lev
40  REAL(r_std),DIMENSION(:,:,:,:),ALLOCATABLE :: soilcarbon_input
41  REAL(r_std),DIMENSION(:,:,:),ALLOCATABLE   :: &
42       &  carbon,control_moist,control_temp
43  REAL(r_std),DIMENSION(:,:),ALLOCATABLE     :: &
44       &  lon,lat,resp_hetero_soil,var_3d
45  REAL(r_std),DIMENSION(:),ALLOCATABLE       :: &
46       &  x_indices
47  REAL(r_std)                                :: time
48  INTEGER                                   :: i,j,m,iatt,iv
49  CHARACTER(LEN=400)                        :: taboo_vars
50  REAL(r_std),DIMENSION(1)                   :: xtmp
51  INTEGER(i_std),PARAMETER                   :: nbvarmax=300
52  INTEGER(i_std)                             :: nbvar
53  CHARACTER(LEN=50),DIMENSION(nbvarmax)     :: varnames
54  INTEGER(i_std)                             :: varnbdim
55  INTEGER(i_std),PARAMETER                   :: varnbdim_max=20
56  INTEGER,DIMENSION(varnbdim_max)           :: vardims
57  LOGICAL                                   :: l1d
58  REAL(r_std)                                :: x_tmp
59  ! clay fraction
60  REAL(r_std),ALLOCATABLE,SAVE,DIMENSION(:)  :: clay
61  !-
62  ! string suffix indicating an index
63  CHARACTER(LEN=10)  :: part_str
64  !
65  CHARACTER(LEN=100) :: Cforcing_name
66  INTEGER            :: Cforcing_id
67  INTEGER            :: v_id
68
69  REAL(r_std),ALLOCATABLE :: clay_loc(:)
70  REAL(r_std),ALLOCATABLE :: soilcarbon_input_loc(:,:,:,:)
71  REAL(r_std),ALLOCATABLE :: control_temp_loc(:,:,:)
72  REAL(r_std),ALLOCATABLE :: control_moist_loc(:,:,:)
73  REAL(r_std),ALLOCATABLE :: carbon_loc(:,:,:)
74  INTEGER :: ierr
75  !>> DS add for externalization
76  LOGICAL  :: l_error
77  ! >> DS
78
79  CALL Init_para(.FALSE.) 
80
81  !
82  ! DS : For externalization cause we decoupled forcesoil from ORCHIDEE
83  !
84 
85  ! 1. Read the number of PFTs
86  CALL getin('NVM',nvm)
87  ! 2. Allocation
88  l_error = .FALSE.
89  ALLOCATE(pft_to_mtc(nvm),stat=ier)
90  l_error = l_error .OR. (ier .NE. 0)
91  IF (l_error) THEN
92     STOP 'pft_to_mtc (forcesoil only) : error in memory allocation'
93  ENDIF
94
95  ! 3. Initialisation of the correspondance table
96  pft_to_mtc (:) = undef_int
97 
98  ! 4.Reading of the conrrespondance table in the .def file
99  CALL getin('PFT_TO_MTC',pft_to_mtc)
100
101  ! 4.1 if nothing is found, we use the standard configuration
102  IF(nvm .EQ. 13 ) THEN
103     IF(pft_to_mtc(1) .EQ. undef_int) THEN
104        WRITE(numout,*) 'Note to the user : we will use ORCHIDEE to its standard configuration'
105        pft_to_mtc(:) = (/ 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13 /)
106     ENDIF
107  ELSE   
108     IF(pft_to_mtc(1) .EQ. undef_int) THEN
109        WRITE(numout,*)' The array PFT_TO_MTC is empty : we stop'
110     ENDIF
111  ENDIF
112 
113  ! 5. What happened if pft_to_mtc(j) > nvmc (if the mtc doesn't exist)?
114  DO i = 1, nvm
115     IF(pft_to_mtc(i) .GT. nvmc) THEN
116        WRITE(numout,*) "the MTC you chose doesn't exist"
117        STOP 'we stop reading pft_to_mtc'
118     ENDIF
119  ENDDO
120 
121 
122  ! 6. Check if pft_to_mtc(1) = 1
123  IF(pft_to_mtc(1) .NE. 1) THEN
124     WRITE(numout,*) 'the first pft has to be the bare soil'
125     STOP 'we stop reading next values of pft_to_mtc'
126  ELSE
127     DO i = 2,nvm
128        IF(pft_to_mtc(i) .EQ.1) THEN
129           WRITE(numout,*) 'only pft_to_mtc(1) has to be the bare soil'
130           STOP 'we stop reading pft_to_mtc'
131        ENDIF
132     ENDDO
133  ENDIF
134 
135
136  ! 7. Allocate and initialize natural ans is_c4
137 
138  ! 7.1 Memory allocation
139  l_error = .FALSE.
140  ALLOCATE(natural(nvm),stat=ier)
141  l_error = l_error .OR. (ier .NE. 0)
142  ALLOCATE(is_c4(nvm),stat=ier)
143
144  IF (l_error) THEN
145     STOP 'natural or is_c4 (forcesoil only) : error in memory allocation'
146  ENDIF
147
148  ! 7.2 Initialisation
149  DO j= 1, nvm
150     natural(j) = natural_mtc(pft_to_mtc(j))
151     is_c4(j) = is_c4_mtc(pft_to_mtc(j))
152  ENDDO
153
154  !-
155  ! Stomate's restart files
156  !-
157  IF (is_root_prc) THEN
158     sto_restname_in = 'stomate_start.nc'
159     CALL getin ('STOMATE_RESTART_FILEIN',sto_restname_in)
160     WRITE(*,*) 'STOMATE INPUT RESTART_FILE: ',TRIM(sto_restname_in)
161     sto_restname_out = 'stomate_restart.nc'
162     CALL getin ('STOMATE_RESTART_FILEOUT',sto_restname_out)
163     WRITE(*,*) 'STOMATE OUTPUT RESTART_FILE: ',TRIM(sto_restname_out)
164     !-
165     ! We need to know iim, jjm.
166     ! Get them from the restart files themselves.
167     !-
168     iret = NF90_OPEN (sto_restname_in, NF90_NOWRITE, ncfid)
169     iret = NF90_INQUIRE_DIMENSION (ncfid,1,len=iim)
170     iret = NF90_INQUIRE_DIMENSION (ncfid,2,len=jjm)
171     iret = NF90_CLOSE (ncfid)
172     !-
173     ! Allocate longitudes and latitudes
174     !-
175     ALLOCATE (lon(iim,jjm))
176     ALLOCATE (lat(iim,jjm))
177     lon(:,:) = 0.0
178     lat(:,:) = 0.0
179     lev(1)   = 0.0
180     !-
181     CALL restini &
182          & (sto_restname_in, iim, jjm, lon, lat, llm, lev, &
183          &  sto_restname_out, itau_dep, date0, dt_files, rest_id_sto)
184  ENDIF
185  !-
186  ! calendar
187  !-
188  CALL bcast(date0)
189!!! MM : à revoir : choix du calendrier dans forcesoil ?? Il est dans le restart de stomate !
190  !  CALL ioconf_calendar ('noleap')
191  CALL ioget_calendar  (one_year,one_day)
192
193  CALL ioconf_startdate(date0)
194
195  IF (is_root_prc) THEN
196     !-
197     ! open FORCESOIL's forcing file to read some basic info
198     !-
199     Cforcing_name = 'stomate_Cforcing.nc'
200     CALL getin ('STOMATE_CFORCING_NAME',Cforcing_name)
201     !-
202     ier = NF90_OPEN (TRIM(Cforcing_name),NF90_NOWRITE,Cforcing_id)
203     !-
204     ier = NF90_GET_ATT (Cforcing_id,NF90_GLOBAL,'kjpindex',x_tmp)
205     kjpindex = NINT(x_tmp)
206     ier = NF90_GET_ATT (Cforcing_id,NF90_GLOBAL,'nparan',x_tmp)
207     nparan = NINT(x_tmp)
208     !-
209     ALLOCATE (indices(kjpindex))
210     ALLOCATE (clay(kjpindex))
211     !-
212     ALLOCATE (x_indices(kjpindex),stat=ier)
213     ier = NF90_INQ_VARID (Cforcing_id,'index',v_id)
214     ier = NF90_GET_VAR   (Cforcing_id,v_id,x_indices)
215     indices(:) = NINT(x_indices(:))
216     DEALLOCATE (x_indices)
217     !-
218     ier = NF90_INQ_VARID (Cforcing_id,'clay',v_id)
219     ier = NF90_GET_VAR   (Cforcing_id,v_id,clay)
220     !-
221     ! time step of forcesoil
222     !-
223     dt_forcesoil = one_year / FLOAT(nparan)
224     WRITE(*,*) 'time step (d): ',dt_forcesoil
225     !-
226     ! read (and partially write) the restart file
227     !-
228     ! read and write the variables we do not need
229     !-
230     taboo_vars = '$lon$ $lat$ $lev$ $nav_lon$ $nav_lat$ $nav_lev$ $time$ $time_steps$ '// &
231          &             '$day_counter$ $dt_days$ $date$ '// &
232          &             '$carbon_01$ $carbon_02$ $carbon_03$ $carbon_04$ $carbon_05$'// &
233          &             '$carbon_06$ $carbon_07$ $carbon_08$ $carbon_09$ $carbon_10$'// &
234          &             '$carbon_11$ $carbon_12$ $carbon_13$'
235     !-
236     CALL ioget_vname(rest_id_sto, nbvar, varnames)
237     !-
238     ! read and write some special variables (1D or variables that we need)
239     !-
240     var_name = 'day_counter'
241     CALL restget (rest_id_sto, var_name, 1, 1, 1, itau_dep, .TRUE., xtmp)
242     CALL restput (rest_id_sto, var_name, 1, 1, 1, itau_dep, xtmp)
243     !-
244     var_name = 'dt_days'
245     CALL restget (rest_id_sto, var_name, 1, 1, 1, itau_dep, .TRUE., xtmp)
246     CALL restput (rest_id_sto, var_name, 1, 1, 1, itau_dep, xtmp)
247     !-
248     var_name = 'date'
249     CALL restget (rest_id_sto, var_name, 1, 1, 1, itau_dep, .TRUE., xtmp)
250     CALL restput (rest_id_sto, var_name, 1, 1, 1, itau_dep, xtmp)
251     !-
252     DO iv=1,nbvar
253        !-- check if the variable is to be written here
254        IF (INDEX(taboo_vars,'$'//TRIM(varnames(iv))//'$') == 0 ) THEN
255           !---- get variable dimensions, especially 3rd dimension
256           CALL ioget_vdim &
257                &      (rest_id_sto, varnames(iv), varnbdim_max, varnbdim, vardims)
258           l1d = ALL(vardims(1:varnbdim) == 1)
259           !----
260           ALLOCATE( var_3d(kjpindex,vardims(3)), stat=ier)
261           IF (ier /= 0) STOP 'ALLOCATION PROBLEM'
262           !---- read it
263           IF (l1d) THEN
264              CALL restget &
265                   &        (rest_id_sto, TRIM(varnames(iv)), 1, vardims(3), &
266                   &         1, itau_dep, .TRUE., xtmp)
267           ELSE
268              CALL restget &
269                   &        (rest_id_sto, TRIM(varnames(iv)), kjpindex, vardims(3), &
270                   &         1, itau_dep, .TRUE., var_3d, "gather", kjpindex, indices)
271           ENDIF
272           !---- write it
273           IF (l1d) THEN
274              CALL restput &
275                   &        (rest_id_sto, TRIM(varnames(iv)), 1, vardims(3), &
276                   &         1, itau_dep, xtmp)
277           ELSE
278              CALL restput &
279                   &        (rest_id_sto, TRIM(varnames(iv)), kjpindex, vardims(3), &
280                   &         1, itau_dep, var_3d, 'scatter',  kjpindex, indices)
281           ENDIF
282           !----
283           DEALLOCATE(var_3d)
284        ENDIF
285     ENDDO
286     !-
287     ! read soil carbon
288     !-
289     ALLOCATE(carbon(kjpindex,ncarb,nvm))
290     carbon(:,:,:) = val_exp
291     DO m = 1, nvm
292        WRITE (part_str, '(I2)') m
293        IF (m<10) part_str(1:1)='0'
294        var_name = 'carbon_'//part_str(1:LEN_TRIM(part_str))
295        CALL restget_p &
296             &    (rest_id_sto, var_name, kjpindex, ncarb , 1, itau_dep, &
297             &     .TRUE., carbon(:,:,m), 'gather', kjpindex, indices)
298        IF (ALL(carbon(:,:,m) == val_exp)) carbon(:,:,m) = zero
299        !-- do not write this variable: it will be modified.
300     ENDDO
301     !-
302     ! Length of run
303     !-
304     WRITE(time_str,'(a)') '10000Y'
305     CALL getin('TIME_LENGTH', time_str)
306     ! transform into itau
307     CALL tlen2itau(time_str, dt_forcesoil*one_year, date0, itau_len)
308     write(*,*) 'Number of time steps to do: ',itau_len
309     !-
310     ! read the rest of the forcing file and store forcing in an array.
311     ! We read an average year.
312     !-
313     ALLOCATE(soilcarbon_input(kjpindex,ncarb,nvm,nparan))
314     ALLOCATE(control_temp(kjpindex,nlevs,nparan))
315     ALLOCATE(control_moist(kjpindex,nlevs,nparan))
316     !-
317     ier = NF90_INQ_VARID (Cforcing_id,'soilcarbon_input',v_id)
318     ier = NF90_GET_VAR   (Cforcing_id,v_id,soilcarbon_input)
319     ier = NF90_INQ_VARID (Cforcing_id,   'control_moist',v_id)
320     ier = NF90_GET_VAR   (Cforcing_id,v_id,control_moist)
321     ier = NF90_INQ_VARID (Cforcing_id,    'control_temp',v_id)
322     ier = NF90_GET_VAR   (Cforcing_id,v_id,control_temp)
323     !-
324     ier = NF90_CLOSE (Cforcing_id)
325     !-
326     !MM Problem here with dpu which depends on soil type           
327     DO iv = 1, nbdl-1
328        ! first 2.0 is dpu
329        ! second 2.0 is average
330        diaglev(iv) = 2.0/(2**(nbdl-1) -1) * ( ( 2**(iv-1) -1) + ( 2**(iv) -1) ) / 2.0
331     ENDDO
332     diaglev(nbdl) = 2.0
333     !-
334     ! For sequential use only, we must initialize data_para :
335     !
336     !
337  ENDIF
338
339  CALL bcast(iim)
340  CALL bcast(jjm)
341  call bcast(kjpindex)
342  CALL init_data_para(iim,jjm,kjpindex,indices)
343
344  !-
345  !-
346  ! there we go: time loop
347  !-
348  CALL bcast(nparan)
349  ALLOCATE(clay_loc(nbp_loc))
350  ALLOCATE(soilcarbon_input_loc(nbp_loc,ncarb,nvm,nparan))
351  ALLOCATE(control_temp_loc(nbp_loc,nlevs,nparan))
352  ALLOCATE(control_moist_loc(nbp_loc,nlevs,nparan))
353  ALLOCATE(carbon_loc(nbp_loc,ncarb,nvm))
354  ALLOCATE(resp_hetero_soil(nbp_loc,nvm))
355  iatt = 0
356
357  CALL bcast(itau_len)
358  CALL bcast(nparan)
359  CALL bcast(dt_forcesoil)
360  CALL Scatter(clay,clay_loc)
361  CALL Scatter(soilcarbon_input,soilcarbon_input_loc)
362  CALL Scatter(control_temp,control_temp_loc)
363  CALL Scatter(control_moist,control_moist_loc)
364  CALL Scatter(carbon,carbon_loc)
365
366!!$ DS 16/06/2011 : calling the new_values of soilcarbon parameters before loop
367     !
368     CALL getin_p('FRAC_CARB_AA',frac_carb_aa)
369     CALL getin_p('FRAC_CARB_AP',frac_carb_ap)   
370     CALL getin_p('FRAC_CARB_SS',frac_carb_ss)
371     CALL getin_p('FRAC_CARB_SA',frac_carb_sa)
372     CALL getin_p('FRAC_CARB_SP',frac_carb_sp)
373     CALL getin_p('FRAC_CARB_PP',frac_carb_pp)
374     CALL getin_p('FRAC_CARB_PA',frac_carb_pa)
375     CALL getin_p('FRAC_CARB_PS',frac_carb_ps)
376     !
377     CALL getin_p('ACTIVE_TO_PASS_CLAY_FRAC',active_to_pass_clay_frac)
378     CALL getin_p('CARBON_TAU_IACTIVE',carbon_tau_iactive)
379     CALL getin_p('CARBON_TAU_ISLOW',carbon_tau_islow)
380     CALL getin_p('CARBON_TAU_IPASSIVE',carbon_tau_ipassive)
381     CALL getin_p('FLUX_TOT_COEFF',flux_tot_coeff)
382
383  DO i=1,itau_len
384     iatt = iatt+1
385     IF (iatt > nparan) iatt = 1
386     CALL soilcarbon &
387          &    (nbp_loc, dt_forcesoil, clay_loc, &
388          &     soilcarbon_input_loc(:,:,:,iatt), &
389          &     control_temp_loc(:,:,iatt), control_moist_loc(:,:,iatt), &
390          &     carbon_loc, resp_hetero_soil)
391  ENDDO
392
393  CALL Gather(carbon_loc,carbon)
394  !-
395  ! write new carbon into restart file
396  !-
397  IF (is_root_prc) THEN
398     DO m=1,nvm
399        WRITE (part_str, '(I2)') m
400        IF (m<10) part_str(1:1)='0'
401        var_name = 'carbon_'//part_str(1:LEN_TRIM(part_str))
402        CALL restput_p &
403             &    (rest_id_sto, var_name, kjpindex, ncarb , 1, itau_dep, &
404             &     carbon(:,:,m), 'scatter', kjpindex, indices)
405     ENDDO
406     !-
407     CALL getin_dump
408     CALL restclo
409  ENDIF
410#ifdef CPP_PARA
411  CALL MPI_FINALIZE(ierr)
412#endif
413  !--------------------
414END PROGRAM forcesoil
Note: See TracBrowser for help on using the repository browser.