source: branches/publications/ORCHIDEE_CAMEO_gmd_2022/src_sechiba/chemistry.f90 @ 7693

Last change on this file since 7693 was 5973, checked in by josefine.ghattas, 5 years ago

Use area instead of resolution. See ticket #563

  • Property svn:keywords set to Date Revision HeadURL
File size: 103.0 KB
Line 
1! ================================================================================================================================
2!  MODULE       : chemistry
3!
4!  CONTACT      : orchidee-help _at_ listes.ipsl.fr
5!
6!  LICENCE      : IPSL (2006)
7!  This software is governed by the CeCILL licence see ORCHIDEE/ORCHIDEE_CeCILL.LIC
8!
9!>\BRIEF   
10!!
11!!\n DESCRIPTION:
12!!
13!! RECENT CHANGE(S): The content of this module were previously part of the diffuco module.
14!!                   ok_inca changed name into ok_bvoc and DIFFUCO_OK_INCA changed into CHEMISTRY_BVOC
15!!                   LEAFAGE_OK_INCA changed name into CHEMISTRY_LEAFAGE
16!!
17!! $HeadURL$
18!! $Date$
19!! $Revision$
20!! \n
21!_ ================================================================================================================================
22
23MODULE chemistry
24
25  USE ioipsl
26  USE constantes
27  USE qsat_moisture
28  USE sechiba_io_p
29  USE ioipsl
30  USE pft_parameters
31  USE grid
32  USE time, ONLY : one_day, dt_sechiba, julian_diff, one_hour, one_year
33  USE ioipsl_para 
34  USE xios_orchidee
35
36  IMPLICIT NONE
37
38  PRIVATE
39  PUBLIC :: chemistry_xios_initialize, chemistry_initialize, chemistry_bvoc, chemistry_flux_interface, chemistry_clear
40
41
42  REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:,:) :: flx_iso            !! Isoprene emissions from vegetation (kgC.m^{-2}.s^{-1})
43!$OMP THREADPRIVATE(flx_iso)
44  REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:,:) :: flx_mono           !! Monoterpene emissions from vegetation (kgC.m^{-2}.s^{-1})
45!$OMP THREADPRIVATE(flx_mono)
46  REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:,:) :: flx_ORVOC          !! Other Volatile Organic Compound emissions from vegetation (kgC.m^{-2}.s^{-1})
47!$OMP THREADPRIVATE(flx_ORVOC)
48  REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:,:) :: flx_MBO            !! 2-methyl-3-buten-2-ol emissions from vegetation (mainly pines in America) (kgC.m^{-2}.s^{-1})
49!$OMP THREADPRIVATE(flx_MBO)
50  REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:,:) :: flx_methanol       !! Methanol emissions from vegetation (kgC.m^{-2}.s^{-1})
51!$OMP THREADPRIVATE(flx_methanol)
52  REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:,:) :: flx_acetone        !! Acetone emissions from vegetation (kgC.m^{-2}.s^{-1})
53!$OMP THREADPRIVATE(flx_acetone)
54  REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:,:) :: flx_acetal         !! Acetaldehyde emissions from vegetation (kgC.m^{-2}.s^{-1})
55!$OMP THREADPRIVATE(flx_acetal)
56  REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:,:) :: flx_formal         !! Formaldehyde emissions from vegetation (kgC.m^{-2}.s^{-1})
57!$OMP THREADPRIVATE(flx_formal)
58  REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:,:) :: flx_acetic         !! Acetic acid emissions from vegetation (kgC.m^{-2}.s^{-1})
59!$OMP THREADPRIVATE(flx_acetic)
60  REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:,:) :: flx_formic         !! Formic acid emissions from vegetation (kgC.m^{-2}.s^{-1})
61!$OMP THREADPRIVATE(flx_formic)
62  REAL(r_std),ALLOCATABLE, SAVE, DIMENSION(:,:)  :: flx_apinen         !! Alpha pinene emissions from vegetation (kgC.m^{-2}.s^{-1})
63!$OMP THREADPRIVATE(flx_apinen)
64  REAL(r_std),ALLOCATABLE, SAVE, DIMENSION(:,:)  :: flx_bpinen         !! Beta pinene emissions from vegetation (kgC.m^{-2}.s^{-1})
65!$OMP THREADPRIVATE(flx_bpinen)
66  REAL(r_std),ALLOCATABLE, SAVE, DIMENSION(:,:)  :: flx_limonen        !! Limonene emissions from vegetation (kgC.m^{-2}.s^{-1})
67!$OMP THREADPRIVATE(flx_limonen)
68  REAL(r_std),ALLOCATABLE, SAVE, DIMENSION(:,:)  :: flx_myrcen         !! Myrcene emissions from vegetation (kgC.m^{-2}.s^{-1})
69!$OMP THREADPRIVATE(flx_myrcen)
70  REAL(r_std),ALLOCATABLE, SAVE, DIMENSION(:,:)  :: flx_sabinen        !! Sabinene emissions from vegetation (kgC.m^{-2}.s^{-1})
71!$OMP THREADPRIVATE(flx_sabinen)
72  REAL(r_std),ALLOCATABLE, SAVE, DIMENSION(:,:)  :: flx_camphen        !! Camphene emissions from vegetation (kgC.m^{-2}.s^{-1})
73!$OMP THREADPRIVATE(flx_camphen)
74  REAL(r_std),ALLOCATABLE, SAVE, DIMENSION(:,:)  :: flx_3caren         !! 3-Carene emissions from vegetation (kgC.m^{-2}.s^{-1})
75!$OMP THREADPRIVATE(flx_3caren)
76  REAL(r_std),ALLOCATABLE, SAVE, DIMENSION(:,:)  :: flx_tbocimen       !! T-beta Ocimene emissions from vegetation (kgC.m^{-2}.s^{-1})
77!$OMP THREADPRIVATE(flx_tbocimen)
78  REAL(r_std),ALLOCATABLE, SAVE, DIMENSION(:,:)  :: flx_othermono      !! Emissions of other monoterpenes from vegetation (kgC.m^{-2}.s^{-1})
79!$OMP THREADPRIVATE(flx_othermono)
80  REAL(r_std),ALLOCATABLE, SAVE, DIMENSION(:,:)  :: flx_sesquiter      !! Sesquiterpene emissions from vegetation (kgC.m^{-2}.s^{-1})
81!$OMP THREADPRIVATE(flx_sesquiter)
82  REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:,:) :: flx_fertil_no      !! Biogenic nitrogen oxide soil emission due to N-fertilisation (ngN.m^{-2}.s^{-1})
83!$OMP THREADPRIVATE(flx_fertil_no)
84  REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:,:) :: flx_no_soil        !! Nitrogen Oxide emissions from soil, before deposition on canopy
85                                                                       !! (ngN.m^{-2}.s^{-1})
86!$OMP THREADPRIVATE(flx_no_soil)
87  REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:,:) :: flx_no             !! Net nitrogen Oxide emissions from soil, after deposition on canopy
88!$OMP THREADPRIVATE(flx_no)
89  REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:,:) :: CRF                !! Canopy reduction factor for net NO flux calculation (kjpindex,nvm)   
90!$OMP THREADPRIVATE(CRF)
91
92  ! variables used inside diffuco_inca module
93  LOGICAL, ALLOCATABLE, SAVE, DIMENSION(:)     :: ok_siesta        !! Flag for controlling post-pulse period (true/false)
94!$OMP THREADPRIVATE(ok_siesta)
95  LOGICAL, ALLOCATABLE, SAVE, DIMENSION(:)     :: allow_pulse      !! Flag for controlling pulse period (true/false)
96!$OMP THREADPRIVATE(allow_pulse)
97  REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: pulse            !! Pulse fonction
98!$OMP THREADPRIVATE(pulse)
99  REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: pulseday         !! Counter for pulse period
100!$OMP THREADPRIVATE(pulseday)
101  REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: siestaday        !! Counter for post-pulse period
102!$OMP THREADPRIVATE(siestaday)
103  REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: pulselim         !! Pulse period length
104!$OMP THREADPRIVATE(pulselim)
105  REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: siestalim        !! Post-pulse period length
106!$OMP THREADPRIVATE(siestalim)
107  REAL(r_std), SAVE                            :: nbre_precip 
108!$OMP THREADPRIVATE(nbre_precip)
109  REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: flx_co2_bbg_year !! CO2 emissions from biomass burning
110                                                                   !! Read in an input file (kgC.m^{-2}.year^{-1})
111!$OMP THREADPRIVATE(flx_co2_bbg_year)
112  REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: N_qt_WRICE_year  !! N fertilizers applied on wetland rice
113                                                                   !! Read in an input file (kgN.yr^{-1})
114!$OMP THREADPRIVATE(N_qt_WRICE_year)
115  REAL(r_std), ALLOCATABLE, SAVE, DIMENSION(:) :: N_qt_OTHER_year  !! N fertilizers applied on other crops and grasses
116                                                                   !! Read in an input file (kgN.yr^{-1})
117!$OMP THREADPRIVATE(N_qt_OTHER_year)
118  LOGICAL, SAVE   :: l_first_chemistry_inca=.TRUE.                 !! Initialisation for chemistry_flux_interface
119!$OMP THREADPRIVATE(l_first_chemistry_inca)
120     
121
122
123CONTAINS
124
125
126!! =============================================================================================================================
127!! SUBROUTINE:    chemistry_xios_initialize
128!!
129!>\BRIEF          Initialize xios dependant defintion before closing context defintion
130!!
131!! DESCRIPTION:   Initialize xios dependant defintion before closing context defintion
132!!
133!! \n
134!_ ==============================================================================================================================
135
136  SUBROUTINE chemistry_xios_initialize
137
138    CHARACTER(LEN=255) :: filename, name
139     
140
141    !! 1. Treat the file for fertilzation needed for option ok_cropsfertil_Nox
142
143    ! Read the input file name from run.def
144    filename = 'orchidee_fertilizer_1995.nc'
145    CALL getin_p('N_FERTIL_FILE',filename) 
146
147    ! Remove suffix .nc from filename
148    name = filename(1:LEN_TRIM(FILENAME)-3)
149    CALL xios_orchidee_set_file_attr("fertilizer_file",name=name)
150
151    ! Check if the file will be read by XIOS, by IOIPSL or not at all
152    IF (xios_interpolation .AND. ok_cropsfertil_Nox) THEN
153       IF (printlev>=2) WRITE (numout,*) 'The fertilizer file will be read later by XIOS'
154    ELSE
155       IF (ok_cropsfertil_Nox) THEN
156          IF (printlev>=2) WRITE (numout,*) 'The fertilizer file will be read later by IOIPSL'
157       ELSE
158          IF (printlev>=2) WRITE (numout,*) 'The fertilizer file will not be read'
159       END IF
160
161       ! The fertilizer file will not be read by XIOS. Now deactivate it for XIOS.
162       CALL xios_orchidee_set_file_attr("fertilizer_file",enabled=.FALSE.)
163       CALL xios_orchidee_set_field_attr("N_qt_WRICE_year_interp",enabled=.FALSE.)
164       CALL xios_orchidee_set_field_attr("N_qt_OTHER_year_interp",enabled=.FALSE.)
165    END IF
166
167
168    !! 2. Treat the file for bbg fertil needed for option ok_bbgfertil_Nox
169
170    ! Read the input file name from run.def
171    filename = 'orchidee_bbg_clim.nc'
172    CALL getin_p('CO2_BBG_FILE',filename)
173
174    ! Remove suffix .nc from filename
175    name = filename(1:LEN_TRIM(FILENAME)-3)
176    CALL xios_orchidee_set_file_attr("bbg_clim_file",name=name)
177
178    ! Check if the file will be read by XIOS, by IOIPSL or not at all
179    IF (xios_interpolation .AND. ok_bbgfertil_Nox) THEN
180       IF (printlev>=2) WRITE (numout,*) 'The bbgfertil file will be read later by XIOS'
181    ELSE
182       IF (ok_bbgfertil_Nox) THEN
183          IF (printlev>=2) WRITE (numout,*) 'The bbgfertil file will be read later by IOIPSL'
184       ELSE
185          IF (printlev>=2) WRITE (numout,*) 'The bbgfertil file will not be read'
186       END IF
187
188       ! This file will not be read by XIOS. Now deactivate it for XIOS.
189       CALL xios_orchidee_set_file_attr("bbg_clim_file",enabled=.FALSE.)
190       CALL xios_orchidee_set_field_attr("flx_co2_bbg_year_interp",enabled=.FALSE.)
191    END IF
192
193
194  END SUBROUTINE chemistry_xios_initialize
195
196!! ================================================================================================================================
197!! SUBROUTINE   : chemistry_initialize
198!!
199!>\BRIEF         This subroutine initializes the chemistry module
200!!
201!! DESCRIPTION  : Some of the variables and flags used chemistry_bvoc are allocated and initialised here.
202!!
203!! RECENT CHANGE(S): Changed name from diffuco_inca_init to chemistry_initialize
204!!
205!! MAIN OUTPUT VARIABLE(S): None
206!!
207!! REFERENCE(S) : None
208!!
209!! FLOWCHART    : None
210!_ ================================================================================================================================
211  SUBROUTINE chemistry_initialize(kjpindex, lalo, neighbours, resolution)
212
213    USE interpweight
214
215    IMPLICIT NONE
216   
217    !! 0. Variables and parameter declaration
218
219    !! 0.1 Input variables
220
221    INTEGER(i_std), INTENT(in)                         :: kjpindex         !! Domain size (unitless)
222    REAL(r_std), DIMENSION(kjpindex,2), INTENT (in)    :: lalo             !! Geographical coordinates
223    INTEGER(i_std), DIMENSION(kjpindex,8), INTENT (in) :: neighbours       !! Vector of neighbours for each
224                                                                           !! grid point (1=N, 2=E, 3=S, 4=W)
225    REAL(r_std),DIMENSION (kjpindex,2), INTENT(in)     :: resolution       !! The size in km of each grid-box in X and Y
226
227    !! 0.2 Output variables
228    REAL(r_std), DIMENSION(kjpindex)                   :: achem_wrice      !! Availability of data for the interpolation
229    REAL(r_std), DIMENSION(kjpindex)                   :: achem_other      !! Availability of data for the interpolation
230    REAL(r_std), DIMENSION(kjpindex)                   :: achem_co2        !! Availability of data for the interpolation
231
232    !! 0.3 Modified variables
233
234    !! 0.4 Local variables
235    LOGICAL                                            :: allow_weathergen
236    CHARACTER(LEN=80)                                  :: filename, fieldname
237    INTEGER(i_std)                                     :: iml, jml, lml, tml, force_id
238    INTEGER(i_std)                                     :: ier
239    REAL(r_std)                                        :: vmin, vmax       !! min/max values to use for the
240                                                                           !!   renormalization
241    CHARACTER(LEN=250)                                 :: maskvname        !! Variable to read the mask from
242                                                                           !! the file
243    CHARACTER(LEN=80)                                  :: lonname, latname !! lon, lat names in input file
244    REAL(r_std), DIMENSION(:), ALLOCATABLE             :: variabletypevals !! Values for all the types of the variable
245                                                                           !!   (variabletypevals(1) = -un, not used)
246    CHARACTER(LEN=50)                                  :: fractype         !! method of calculation of fraction
247                                                                           !!   'XYKindTime': Input values are kinds
248                                                                           !!     of something with a temporal
249                                                                           !!     evolution on the dx*dy matrix'
250    LOGICAL                                            :: nonegative       !! whether negative values should be removed
251    CHARACTER(LEN=50)                                  :: maskingtype      !! Type of masking
252                                                                           !!   'nomask': no-mask is applied
253                                                                           !!   'mbelow': take values below maskvals(1)
254                                                                           !!   'mabove': take values above maskvals(1)
255                                                                           !!   'msumrange': take values within 2 ranges;
256                                                                           !!      maskvals(2) <= SUM(vals(k)) <= maskvals(1)
257                                                                           !!      maskvals(1) < SUM(vals(k)) <= maskvals(3)
258                                                                           !!        (normalized by maskvals(3))
259                                                                           !!   'var': mask values are taken from a
260                                                                           !!     variable inside the file (>0)
261    REAL(r_std), DIMENSION(3)                          :: maskvals         !! values to use to mask (according to
262                                                                           !!   `maskingtype')
263    CHARACTER(LEN=250)                                 :: namemaskvar      !! name of the variable to use to mask
264    REAL(r_std)                                        :: chem_norefinf    !! No value
265    REAL(r_std)                                        :: chem_default     !! Default value
266
267!_ ================================================================================================================================
268
269    ALLOCATE (pulse(kjpindex),stat=ier)
270    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_initialize','Problem in allocate of variable pulse','','')
271    pulse(:) = un
272
273    ! If we acount for NOx pulse emissions
274    IF (ok_pulse_NOx) THEN
275
276       ALLOCATE (ok_siesta(kjpindex),stat=ier)
277       IF (ier /= 0) CALL ipslerr_p(3,'chemistry_initialize','Problem in allocate of variable ok_siesta','','')
278       ok_siesta(:) = .FALSE.
279
280       ALLOCATE (allow_pulse(kjpindex),stat=ier)
281       IF (ier /= 0) CALL ipslerr_p(3,'chemistry_initialize','Problem in allocate of variable allow_pulse','','')
282       allow_pulse(:) = .FALSE.
283
284       ALLOCATE (pulseday(kjpindex),stat=ier)
285       IF (ier /= 0) CALL ipslerr_p(3,'chemistry_initialize','Problem in allocate of variable pulseday','','')
286       pulseday(:) = zero
287
288       ALLOCATE (siestaday(kjpindex),stat=ier)
289       IF (ier /= 0) CALL ipslerr_p(3,'chemistry_initialize','Problem in allocate of variable siestaday','','')
290       siestaday(:) = zero
291
292       ALLOCATE (pulselim(kjpindex),stat=ier)
293       IF (ier /= 0) CALL ipslerr_p(3,'chemistry_initialize','Problem in allocate of variable pulselim','','')
294       pulselim(:) = zero
295
296       ALLOCATE (siestalim(kjpindex),stat=ier)
297       IF (ier /= 0) CALL ipslerr_p(3,'chemistry_initialize','Problem in allocate of variable siestalim','','')
298       siestalim(:) = zero
299
300    END IF ! (ok_pulse_NOx)
301
302    ! If we acount for NOx emissions by N-fertilizers
303    IF (ok_cropsfertil_NOx) THEN
304
305       ALLOCATE (N_qt_WRICE_year(kjpindex),stat=ier)  !! N fertilizers on wetland rice, read in file
306       IF (ier /= 0) CALL ipslerr_p(3,'chemistry_initialize','Problem in allocate of variable N_qt_WRICE_year','','')
307       N_qt_WRICE_year(:) = zero
308   
309       ALLOCATE (N_qt_OTHER_year(kjpindex),stat=ier)  !! N fertilizers on other crops and grasses, read in file
310       IF (ier /= 0) CALL ipslerr_p(3,'chemistry_initialize','Problem in allocate of variable N_qt_OTHER_year','','')
311       N_qt_OTHER_year(:) = zero
312
313       WRITE (numout,*) ' *********************** Interpolating N fertilizers files for NOx emissions... '
314
315       !Config Key   = N_FERTIL_FILE
316       !Config Desc  = File name
317       !Config If    = CHEMISTRY_BVOC and NOx_FERTILIZERS_USE
318       !Config Def   = orchidee_fertilizer_1995.nc
319       !Config Help  =
320       !Config Units = -
321       filename = 'orchidee_fertilizer_1995.nc'
322       CALL getin_p('N_FERTIL_FILE',filename) 
323
324       IF ( xios_interpolation ) THEN
325   
326         CALL xios_orchidee_recv_field('N_qt_WRICE_year_interp',N_qt_WRICE_year)
327         CALL xios_orchidee_recv_field('N_qt_OTHER_year_interp',N_qt_OTHER_year)   
328
329         achem_wrice(:)=1
330         achem_other(:)=1
331       ELSE
332       
333         !! Variables for interpweight
334         vmin = 0.
335         vmax = 0.
336         ! Type of calculation of cell fractions
337         fractype = 'default'
338         ! Name of the longitude and latitude in the input file
339         lonname = 'lon'
340         latname = 'lat'
341         ! Default value when no value is get from input file
342         chem_default = 0.
343         ! Reference value when no value is get from input file
344         chem_norefinf = 0.
345         ! Should negative values be set to zero from input file?
346         nonegative = .TRUE.
347         ! Type of mask to apply to the input data (see header for more details)
348         maskingtype = 'nomask'
349         ! Values to use for the masking
350         maskvals = (/ min_sechiba, undef_sechiba, undef_sechiba /)
351         ! Name of the variable with the values for the mask in the input file (only if maskkingtype='var') (here not used)
352         namemaskvar = ''
353
354         fieldname= 'N_qt_WRICE_year'
355         IF (printlev >= 1) WRITE(numout,*) "chemistry_initialize: Read and interpolate file " &
356              // TRIM(filename) //" for variable N_qt_WRICE_year"
357         CALL interpweight_2Dcont(kjpindex, 0, 0, lalo, resolution, neighbours,                        &
358              contfrac, filename, fieldname, lonname, latname, vmin, vmax, nonegative, maskingtype,       &
359              maskvals, namemaskvar, -1, fractype, chem_default, chem_norefinf,                           &
360              N_qt_WRICE_year, achem_wrice)
361       
362         fieldname= 'N_qt_OTHER_year'
363         IF (printlev >= 1) WRITE(numout,*) "chemistry_initialize: Read and interpolate file " &
364              // TRIM(filename) //" for variable N_qt_OTHER_year"
365         CALL interpweight_2Dcont(kjpindex, 0, 0, lalo, resolution, neighbours,                        &
366              contfrac, filename, fieldname, lonname, latname, vmin, vmax, nonegative, maskingtype,       &
367              maskvals, namemaskvar, -1, fractype, chem_default, chem_norefinf,                           &
368              N_qt_OTHER_year, achem_other)
369       END IF   
370    END IF
371
372    ALLOCATE (flx_iso(kjpindex,nvm), stat=ier) 
373    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable flx_iso','','')
374    flx_iso(:,:) = 0. 
375
376    ALLOCATE (flx_mono(kjpindex,nvm), stat=ier) 
377    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable flx_mono','','')
378    flx_mono(:,:) = 0. 
379
380    ALLOCATE (flx_ORVOC(kjpindex,nvm), stat=ier) 
381    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable flx_ORVOC ','','')
382    flx_ORVOC(:,:) = 0. 
383
384    ALLOCATE (flx_MBO(kjpindex,nvm), stat=ier) 
385    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable flx_MBO','','')
386    flx_MBO(:,:) = 0. 
387
388    ALLOCATE (flx_methanol(kjpindex,nvm), stat=ier) 
389    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable flx_methanol','','')
390    flx_methanol(:,:) = 0. 
391
392    ALLOCATE (flx_acetone(kjpindex,nvm), stat=ier) 
393    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable flx_acetone','','')
394    flx_acetone(:,:) = 0. 
395
396    ALLOCATE (flx_acetal(kjpindex,nvm), stat=ier) 
397    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable flx_acetal','','')
398    flx_acetal(:,:) = 0. 
399
400    ALLOCATE (flx_formal(kjpindex,nvm), stat=ier) 
401    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable flx_formal','','')
402    flx_formal(:,:) = 0. 
403
404    ALLOCATE (flx_acetic(kjpindex,nvm), stat=ier) 
405    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable flx_acetic','','')
406    flx_acetic(:,:) = 0. 
407
408    ALLOCATE (flx_formic(kjpindex,nvm), stat=ier) 
409    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable flx_formic','','')
410    flx_formic(:,:) = 0. 
411
412    ALLOCATE (flx_no_soil(kjpindex,nvm), stat=ier) 
413    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable flx_no_soil','','')
414    flx_no_soil(:,:) = 0. 
415
416    ALLOCATE (flx_no(kjpindex,nvm), stat=ier) 
417    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable flx_no','','')
418    flx_no(:,:) = 0. 
419       
420    ALLOCATE (flx_fertil_no(kjpindex,nvm), stat=ier) 
421    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable flx_fertil_no','','')
422    flx_fertil_no(:,:) = 0. 
423
424    ALLOCATE (flx_apinen(kjpindex,nvm), stat=ier) 
425    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable flx_apinen','','')
426    flx_apinen(:,:) = 0.       
427
428    ALLOCATE (flx_bpinen (kjpindex,nvm), stat=ier) 
429    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable flx_bpinen','','')
430    flx_bpinen(:,:) = 0.     
431
432    ALLOCATE (flx_limonen  (kjpindex,nvm), stat=ier) 
433    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable flx_limonen','','')
434    flx_limonen(:,:) = 0.   
435
436    ALLOCATE (flx_myrcen(kjpindex,nvm), stat=ier) 
437    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable flx_myrcen','','')
438    flx_myrcen(:,:) = 0.       
439
440    ALLOCATE (flx_sabinen(kjpindex,nvm), stat=ier) 
441    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable flx_sabinen','','')
442    flx_sabinen(:,:) = 0.     
443
444    ALLOCATE (flx_camphen(kjpindex,nvm), stat=ier) 
445    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable flx_camphen','','')
446    flx_camphen(:,:) = 0.     
447
448    ALLOCATE (flx_3caren(kjpindex,nvm), stat=ier) 
449    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable flx_3caren','','')
450    flx_3caren(:,:) = 0.       
451
452    ALLOCATE (flx_tbocimen(kjpindex,nvm), stat=ier) 
453    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable flx_tbocimen','','')
454    flx_tbocimen(:,:) = 0.     
455
456    ALLOCATE (flx_othermono(kjpindex,nvm), stat=ier) 
457    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable flx_othermono','','')
458    flx_othermono(:,:) = 0.   
459
460    ALLOCATE (flx_sesquiter(kjpindex,nvm), stat=ier) 
461    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable flx_sesquiter','','')
462    flx_sesquiter(:,:) = 0.   
463
464    ALLOCATE(CRF(kjpindex,nvm), stat=ier) 
465    IF (ier /= 0) CALL ipslerr_p(3,'chemistry_init','Problem in allocate of variable CRF ','','')
466    CRF(:,:) = 0. 
467
468
469    ! If we acount for NOx emissions due to Biomass Burning
470    IF (ok_bbgfertil_NOx) THEN
471
472       ALLOCATE (flx_co2_bbg_year(kjpindex),stat=ier) !! CO2 emissions from bbg, read in file
473       IF (ier /= 0) CALL ipslerr_p(3,'chemistry_initialize','Problem in allocate of variable flx_co2_bbg_year','','')
474       flx_co2_bbg_year(:) = zero   
475
476       WRITE (numout,*) ' *********************** Interpolating CO2 bbg files for NOx emissions... '
477       !Config Key   = N_FERTIL_FILE
478       !Config Desc  = File name
479       !Config If    = CHEMISTRY_BVOC and NOx_FERTILIZERS_USE
480       !Config Def   = orchidee_fertilizer_1995.nc
481       !Config Help  = ...
482       !Config Units = -
483       filename = 'orchidee_bbg_clim.nc'
484       CALL getin_p('CO2_BBG_FILE',filename)
485
486       IF (xios_interpolation) THEN
487   
488         CALL xios_orchidee_recv_field('flx_co2_bbg_year_interp',flx_co2_bbg_year)
489
490         achem_co2(:)=1
491       ELSE
492 
493         !! Variables for interpweight
494         vmin = 0.
495         vmax = 0.
496         ! Type of calculation of cell fractions
497         fractype = 'default'
498         ! Name of the longitude and latitude in the input file
499         lonname = 'lon'
500         latname = 'lat'
501         ! Default value when no value is get from input file
502         chem_default = 0.
503         ! Reference value when no value is get from input file
504         chem_norefinf = 0.
505         ! Should negative values be set to zero from input file?
506         nonegative = .TRUE.
507         ! Type of mask to apply to the input data (see header for more details)
508         maskingtype = 'nomask'
509         ! Values to use for the masking
510         maskvals = (/ min_sechiba, undef_sechiba, undef_sechiba /)
511         ! Name of the variable with the values for the mask in the input file (only if maskkingtype='var') (here not used)
512         namemaskvar = ''
513
514         fieldname = 'flx_co2_bbg_year'
515         IF (printlev >= 1) WRITE(numout,*) "chemistry_initialize: Read and interpolate file " &
516              // TRIM(filename) //" for variable flx_co2_bbg_year"
517         CALL interpweight_2Dcont(kjpindex, 0, 0, lalo, resolution, neighbours,                         &
518              contfrac, filename, fieldname, lonname, latname, vmin, vmax, nonegative, maskingtype,     &
519              maskvals, namemaskvar, -1, fractype, chem_default, chem_norefinf,                         &
520              flx_co2_bbg_year, achem_co2)
521       END IF
522    END IF
523
524    IF ( OFF_LINE_MODE ) THEN
525
526       !-
527       !- What are the alowed options for the temporal interpolation
528       !-
529       ! ALLOW_WEATHERGEN : Allow weather generator to create data
530       ! This parameter is already read in the driver
531       allow_weathergen = .FALSE.
532       CALL getin_p('ALLOW_WEATHERGEN',allow_weathergen)
533       
534       ! FORCING_FILE : Name of file containing the forcing data
535       ! This parameter is already read in the driver
536       filename='forcing_file.nc'
537       CALL getin_p('FORCING_FILE',filename)
538       CALL flininfo(filename,iml, jml, lml, tml, force_id)   
539       WRITE(numout,*) 'Number of data per year in forcing file :', tml 
540       CALL flinclo(force_id)
541       WRITE(numout,*) 'Forcing file closed in chemistry_initialize'
542       
543       
544       IF ( allow_weathergen ) THEN
545          WRITE(numout,*) '**chemistry_initialize: Using weather generator, careful to precip division for NOx '
546          nbre_precip = un
547          WRITE(numout,*) 'Division pour les precip, NOx:', nbre_precip
548       ELSE
549          WRITE(numout,*) 'DT_SECHIBA :', dt_sechiba
550          nbre_precip = (one_day/dt_sechiba)/(tml/one_year)
551          WRITE(numout,*) 'Division pour les precip, NOx:', nbre_precip
552       END IF
553
554    ELSE ! (in coupled mode)
555
556       nbre_precip = un
557       
558    END IF  ! (OFF_LINE_MODE)
559
560    ! Write diagnostics
561    IF (ok_cropsfertil_NOx) THEN
562      CALL xios_orchidee_send_field("achem_wrice",achem_wrice)
563      CALL xios_orchidee_send_field("achem_other",achem_other)
564      CALL xios_orchidee_send_field("interp_diag_N_qt_WRICE_year",N_qt_WRICE_year)
565      CALL xios_orchidee_send_field("interp_diag_N_qt_OTHER_year",N_qt_OTHER_year)
566
567    END IF
568    IF (ok_bbgfertil_NOx) THEN
569      CALL xios_orchidee_send_field("achem_co2",achem_co2)
570      CALL xios_orchidee_send_field("interp_diag_flx_co2_bbg_year",flx_co2_bbg_year)
571    END IF
572
573       
574  END SUBROUTINE chemistry_initialize
575
576
577!! ================================================================================================================================
578!! SUBROUTINE   : chemistry_bvoc
579!!
580!>\BRIEF         This subroutine computes biogenic emissions of reactive compounds, that is of
581!!               VOCs (volatile organic compounds) from vegetation and NOx (nitrogen oxides) from soils.
582!!               Calculation are mostly based on the works by Guenther et al. (1995) and Yienger and Levy (1995).\n
583!!
584!! DESCRIPTION  : Biogenic VOC emissions from vegetation are based on the parameterisations developped by
585!!                Guenther et al. (1995). Biogenic VOCs considered here are: isoprene, monoterpenes, OVOC and ORVOC
586!!                as bulked emissions, methanol, acetone, acetaldehyde, formaldehyde, acetic acid, formic acid
587!!                as single emissions.\n
588!!                For every biogenic VOCs an emission factor (EF), depending on the PFT considered, is used.\n
589!!                Isoprene emissions depend on temperature and radiation. A partition between sunlit and shaded
590!!                leaves is taken into account and either one (if ok_multilayer = FALSE) or several layers
591!!                (if ok_multilayer = TRUE) in the canopy can be used.\n
592!!                When radiation extinction is considered, the canopy radiative transfer model takes into
593!!                account light extinction through canopy, calculating first need diffuse and direct radiation
594!!                based on Andrew Friend 2001 radiative model and Spitters et al. 1986. The calculation of lai,
595!!                parscat, parsh and parsun, laisun and laishabsed based on Guenther et al.(JGR, 1995) and Norman (1982).\n
596!!                Emissions for other BVOCs (monoterpenes, OVOC, ORVOC and other single compounds such as
597!!                methanol, acetone...) depend only on temperature.\n   
598!!                The impact of leaf age, using an emission activity prescribed for each of the 4 leaf age
599!!                classes can also be considered for isoprene and methanol emissions when ok_leafage = TRUE.\n
600!!                NOx emissions from soils are based on Yienger and Levy (1995) and depend on soil moisture
601!!                and temperature and PFT. The pulse effect, related to significant rain occuring after severe
602!!                drought can also be considered (ok_pulse_NOx = TRUE), as well as the increase in emissions related to
603!!                biomass buring (ok_bbgfertil_NOx = TRUE) or use of fertilizers (ok_cropsfertil_NOx = TRUE).
604!!                A net NO flux is eventually calculated taking into account loss by deposition on the surface, using
605!!                a Canopy Reduction Factor (CRF) based on stomatal and leaf area.\n
606!!                This subroutine is called by diffuco_main only if biogenic emissions are activated
607!!                for sechiba (flag CHEMISTRY_BVOC=TRUE).\n
608!!
609!! RECENT CHANGE(S): Changed name from diffuco_inca to chemistry_bvoc
610!!
611!! MAIN OUTPUT VARIABLE(S): :: PAR, :: PARsun, :: PARsh, :: laisun, :: laish,
612!!                          :: flx_iso, :: flx_mono, :: flx_ORVOC, :: flx_MBO,
613!!                          :: flx_methanol, :: flx_acetone, :: flx_acetal, :: flx_formal,
614!!                          :: flx_acetic, :: flx_formic, :: flx_no_soil, :: flx_no,
615!!                          :: CRF, :: flx_fertil_no, :: Trans, :: Fdf,
616!!                          :: PARdf, :: PARdr, :: PARsuntab, :: PARshtab
617!!
618!! REFERENCE(S) :
619!! - Andrew Friend (2001), Modelling canopy CO2 fluxes: are 'big-leaf' simplifications justified?
620!! Global Ecology and Biogeography, 10, 6, 603-619, doi: 10.1046/j.1466-822x.2001.00268.x
621!! - Spitters, C.J.T, Toussaint, H.A.J.M, Groudriaan, J. (1986), Separating the diffuse and direct
622!! component of global radiation and its implications for modeling canopy photosynthesis, Agricultural
623!! and Forest Meteorology, 38, 1-3, 217-229, doi:10.1016/0168-1923(86)90060-2
624!! - Norman JM (1982) Simulation of microclimates. In: Hatfield JL, Thomason IJ (eds)
625!!  Biometeorology in integrated pest management. Academic, New York, pp 65–99
626!! - Guenther, A., Hewitt, C. N., Erickson, D., Fall, R., Geron, C., Graedel, T., Harley, P.,
627!! Klinger, L., Lerdau, M., McKay, W. A., Pierce, T., Scholes, B., Steinbrecher, R., Tallamraju,
628!! R., Taylor, J. et Zimmerman, P. (1995), A global model of natural volatile organic compound
629!! emissions, J. Geophys. Res., 100, 8873-8892.
630!! - MacDonald, R. et Fall, R. (1993), Detection of substantial emissions of methanol from
631!! plants to the atmosphere, Atmos. Environ., 27A, 1709-1713.
632!! - Guenther, A., Geron, C., Pierce, T., Lamb, B., Harley, P. et Fall, R. (2000), Natural emissions
633!! of non-methane volatile organic compounds, carbon monoxide, and oxides of nitrogen from
634!! North America, Atmos. Environ., 34, 2205-2230.
635!! - Yienger, J. J. et Levy II, H. (1995), Empirical model of global soil-biogenic NOx emissions,
636!! J. Geophys. Res., 100, 11,447-11,464.
637!! - Lathiere, J., D.A. Hauglustaine, A. Friend, N. De Noblet-Ducoudre, N. Viovy, and
638!!  G. Folberth (2006), Impact of climate variability and land use changes on global biogenic volatile
639!! organic compound emissions, Atmospheric Chemistry and Physics, 6, 2129-2146.
640!! - Lathiere, J., D.A. Hauglustaine, N. De Noblet-Ducoudre, G. Krinner et G.A. Folberth (2005),
641!! Past and future changes in biogenic volatile organic compound emissions simulated with a global
642!! dynamic vegetation model, Geophysical Research Letters, 32, doi: 10.1029/2005GL024164.
643!! - Lathiere, J. (2005), Evolution des emissions de composes organiques et azotes par la biosphere
644!!  continentale dans le modele LMDz-INCA-ORCHIDEE, These de doctorat, Universite Paris VI.
645!!
646!! FLOWCHART    : None
647!_ ================================================================================================================================
648
649  SUBROUTINE chemistry_bvoc (kjpindex, swdown, coszang, temp_air, &
650       temp_sol, ptnlev1, precip_rain, humrel, veget_max, lai, &
651       frac_age, lalo, ccanopy, cim,  wind, snow, &
652       veget, hist_id, hist2_id,kjit, index, &
653       indexlai, indexveg)
654
655    !! 0. Variables and parameter declaration
656
657    !! 0.1 Input variables
658
659    INTEGER(i_std), INTENT(in)                                 :: kjpindex         !! Domain size - terrestrial pixels only (unitless)
660    INTEGER(i_std), INTENT(in)                                 :: kjit             !! Time step number (-)
661    INTEGER(i_std),INTENT (in)                                 :: hist_id          !! History file identifier (-)
662    INTEGER(i_std),INTENT (in)                                 :: hist2_id         !! History file 2 identifier (-)
663    INTEGER(i_std),DIMENSION (kjpindex), INTENT (in)           :: index            !! Indeces of the points on the map (-)
664    INTEGER(i_std),DIMENSION (kjpindex*(nlai+1)), INTENT (in)  :: indexlai         !! Indeces of the points on the 3D map
665    INTEGER(i_std),DIMENSION (kjpindex*nvm), INTENT (in)       :: indexveg         !! Indeces of the points on the 3D map (-)
666    REAL(r_std), DIMENSION(kjpindex), INTENT(in)               :: swdown           !! Down-welling surface short-wave flux
667                                                                                   !! (W.m^{-2})
668    REAL(r_std), DIMENSION(kjpindex), INTENT(in)               :: coszang          !! Cosine of the solar zenith angle (unitless)
669    REAL(r_std), DIMENSION(kjpindex), INTENT(in)               :: temp_air         !! Air temperature (K)
670    REAL(r_std), DIMENSION(kjpindex), INTENT(in)               :: temp_sol         !! Skin temperature (K)
671    REAL(r_std), DIMENSION(kjpindex), INTENT(in)               :: ptnlev1          !! 1st level of soil temperature (K)
672    REAL(r_std), DIMENSION(kjpindex), INTENT(in)               :: precip_rain      !! Rain precipitation !!?? init
673    REAL(r_std), DIMENSION(kjpindex,nvm), INTENT(in)           :: humrel           !! Soil moisture stress (0-1, unitless)
674    REAL(r_std), DIMENSION(kjpindex,nvm), INTENT(in)           :: veget_max        !! Max. vegetation fraction (0-1, unitless)
675    REAL(r_std),DIMENSION (kjpindex), INTENT (in)              :: snow             !! Snow mass (kg)
676    REAL(r_std),DIMENSION (kjpindex,nvm), INTENT (in)          :: veget            !! Fraction of vegetation type (-)
677    REAL(r_std), DIMENSION(kjpindex,nvm), INTENT(in)           :: lai              !! Leaf area index (m^2.m^{-2})
678    REAL(r_std), DIMENSION(kjpindex,nvm,nleafages), INTENT(in) :: frac_age         !! Age efficacity from STOMATE for iso
679    REAL(r_std), DIMENSION(kjpindex,2), INTENT(in)             :: lalo             !! Geographical coordinates for pixels (degrees)
680    REAL(r_std),DIMENSION (kjpindex), INTENT (in)              :: ccanopy          !! CO2 concentration inside the canopy
681    REAL(r_std),DIMENSION (kjpindex,nvm), INTENT (in)          :: cim              !! Intercellular CO2 over nlai
682    REAL(r_std), DIMENSION (kjpindex), INTENT(in)              :: wind             !! Wind module (m s^{-1})
683
684    !! 0.2 Output variables
685
686    !! 0.3 Modified variables
687
688    !! 0.4 Local variables
689
690    INTEGER(i_std)                             :: ji, jv, jf, jl    !! Indices (unitless)
691    REAL(r_std), DIMENSION(kjpindex,nvm)       :: fol_dens          !! foliar density (gDM.m^{-2})
692    REAL(r_std), DIMENSION(kjpindex)           :: tleaf             !! Foliar temperature (K)
693    REAL(r_std), DIMENSION(kjpindex)           :: t_no              !! Temperature used for soil NO emissions (C)
694    REAL(r_std), DIMENSION(kjpindex)           :: exp_1             !! First exponential used in the calculation of
695                                                                    !! isoprene dependancy to Temperature
696    REAL(r_std), DIMENSION(kjpindex)           :: exp_2             !! Second exponential used in the calculation of
697                                                                    !! Isoprene dependancy to Temperature
698    REAL(r_std), DIMENSION(kjpindex)           :: Ct_iso            !! Isoprene dependancy to Temperature
699    REAL(r_std), DIMENSION(kjpindex)           :: Cl_iso            !! Isoprene dependancy to Light
700    REAL(r_std), DIMENSION(kjpindex)           :: Ct_mono           !! Monoterpene dependancy to Temperature
701    REAL(r_std), DIMENSION(kjpindex)           :: Ct_sesq           !! Sesquiterpenes dependancy to Temperature
702    REAL(r_std), DIMENSION(kjpindex)           :: Ct_meth           !! Methanol dependancy to Temperature
703    REAL(r_std), DIMENSION(kjpindex)           :: Ct_acet           !! Acetone dependancy to Temperature
704    REAL(r_std), DIMENSION(kjpindex)           :: Ct_oxyVOC         !! Other oxygenated BVOC dependancy to Temperature
705    REAL(r_std)                                :: GAMMA_iso         !! Temperature and light dependancy for isoprene and fo each PFT
706    REAL(r_std)                                :: GAMMA_iso_m       !! Temperature and light dependancy for isoprene and fo each PFT for multilayer
707    REAL(r_std), DIMENSION(kjpindex)           :: Ylt_mono          !! Total Temperature and light dependancy for monoterpenes
708    REAL(r_std), DIMENSION(kjpindex)           :: Ylt_sesq          !! Total Temperature and light dependancy for sesquiterpens
709    REAL(r_std), DIMENSION(kjpindex)           :: Ylt_meth          !! Total Temperature and light dependancy for methanol
710    REAL(r_std), DIMENSION(kjpindex)           :: Ylt_acet          !! Total Temperature and light dependancy for acetone
711    REAL(r_std), DIMENSION(kjpindex)           :: Ct_MBO            !! MBO dependance to Temperature
712    REAL(r_std), DIMENSION(kjpindex)           :: Cl_MBO            !! MBO dependance to Light
713    REAL(r_std), DIMENSION(kjpindex)           :: Xvar              !! Parameter used in the calculation
714                                                                    !! of MBO dependance to Temperature
715    REAL(r_std), DIMENSION(kjpindex,nvm)       :: flx_OVOC          !! Biogenic OVOC emission -
716                                                                    !! Other Volatil Organic Components (kgC.m^{-2}.s^{-1})
717    !!Canopy radiative transfer model
718    REAL(r_std), DIMENSION(kjpindex)           :: So                !! Maximum radiation at the Earth surface (W.m^{-2})
719    REAL(r_std), DIMENSION(kjpindex)           :: Rfrac             !! Parameter in the regression of diffuse
720                                                                    !! share on transmission
721    REAL(r_std), DIMENSION(kjpindex)           :: Kfrac             !! Parameter in the regression of diffuse
722                                                                    !! share on transmission
723    REAL(r_std), DIMENSION(kjpindex)           :: swdf              !! Sw diffuse radiation (W.m^{-2})
724    REAL(r_std), DIMENSION(kjpindex)           :: swdr              !! Sw direct radiation (W.m^{-2})
725    REAL(r_std), DIMENSION(kjpindex,nvm)       :: PARscat           !! Scatter PAR @tex ($\mu mol m^{-2} s^{-1}$) @endtex
726    REAL(r_std), DIMENSION(kjpindex,nvm)       :: Clsun_iso         !! Isoprene dependance to light for sun leaves
727    REAL(r_std), DIMENSION(kjpindex,nvm)       :: Clsh_iso          !! Isoprene dependance to light for shaded leaves
728    !! for multilayer canopy model for iso flux
729    REAL(r_std), DIMENSION(kjpindex,nlai+1)    :: PARscattab        !! Scatter PAR @tex ($\mu mol m^{-2} s^{-1}$) @endtex
730    REAL(r_std), DIMENSION(nlai+1)             :: laitab            !! LAI per layer (m^2.m^{-2})
731    REAL(r_std), DIMENSION(kjpindex,nlai)      :: laisuntabdep      !! LAI of sun leaves in each layer (m^2.m^{-2})
732    REAL(r_std), DIMENSION(kjpindex,nlai)      :: laishtabdep       !! LAI of shaded leaves in each layer
733                                                                    !! (m^2.m^{-2})
734    REAL(r_std)                                :: Clsun_iso_tab     !! Isoprene dependance to light
735                                                                    !! for sun leaves and per layer
736    REAL(r_std)                                :: Clsh_iso_tab      !! Isoprene dependance to light
737                                                                    !! for shaded leaves and per layer
738    !for multilayer canopy model Spitter et al. 1986
739    REAL(r_std), DIMENSION(kjpindex,nlai+1)    :: PARnotscat        !! Not-Scattered PAR
740    REAL(r_std), DIMENSION(kjpindex,nlai+1)    :: PARabsdir         !! Absorbed light of the PAR direct flux 
741    REAL(r_std), DIMENSION(kjpindex,nlai+1)    :: PARabsdiff        !! Absorbed light of the PAR diffuse flux 
742    REAL(r_std), PARAMETER                     :: sigma = 0.20      !! scattering coefficient of single leaves and for visible radiation
743    REAL(r_std), PARAMETER                     :: cluster = 0.85    !! clustering coefficient for leaves, the same that is setting for default in MEGAN V2.10
744    REAL(r_std)                                :: rho               !! reflection index of a green, closed vegetation
745    REAL(r_std)                                :: kbl               !! extinction coefficient of black leaves
746    REAL(r_std)                                :: kdf               !! extinction coefficient of diffuse flux
747    !!Leaf age
748    REAL(r_std), DIMENSION(kjpindex,nvm)       :: Eff_age_iso       !! Isoprene emission dependance to Leaf Age
749    REAL(r_std), DIMENSION(kjpindex,nvm)       :: Eff_age_meth      !! Methanol emission dependance to Leaf Age
750    REAL(r_std), DIMENSION(kjpindex,nvm)       :: Eff_age_VOC       !! Other VOC emission dependance to Leaf Age
751    !!BBG and Fertilizers for NOx soil emission
752    REAL(r_std), DIMENSION(kjpindex)           :: veget_max_nowoody !! sum of veget_max for non-woody PFT
753    REAL(r_std), DIMENSION(kjpindex,nvm)       :: N_qt_WRICE_pft    !! N fertiliser on RICE
754                                                                    !! (kgN per year per grid cell)
755    REAL(r_std), DIMENSION(kjpindex,nvm)       :: N_qt_OTHER_pft    !! N fertiliser on other veg
756                                                                    !! (kgN per year per grid cell)
757    !! CO2 inhibition effect on isoprene
758    REAL(r_std),DIMENSION (kjpindex,nvm)       :: fco2_wshort       !! Wilkinson short term function for CO2 impact on BVOC (isoprene)
759    REAL(r_std),DIMENSION (kjpindex)           :: fco2_wlong        !! Wilkinson long term function for CO2 impact on BVOC (isoprene)
760    REAL(r_std)                                :: fco2_ctrl
761    REAL(r_std),DIMENSION (kjpindex,nvm)       :: fco2              !! Function for CO2 impact on BVOC (isoprene)
762    REAL(r_std), DIMENSION(kjpindex)           :: Ismax_short
763    REAL(r_std), DIMENSION(kjpindex)           :: h_short
764    REAL(r_std), DIMENSION(kjpindex)           :: Cstar_short
765    REAL(r_std)                                :: Ismax_long
766    REAL(r_std)                                :: h_long
767    REAL(r_std)                                :: Cstar_long
768
769    !! 0.5 Parameters values
770
771    REAL(r_std), PARAMETER :: CT1 = 95000.0       !! Empirical coeffcient (see Guenther .et. al, 1995, eq(10)) (J.mol^{-1})
772    REAL(r_std), PARAMETER :: CT2 = 230000.0      !! Empirical coefficient (see Guenther .et. al, 1995, eq(10)) (J.mol^{-1})
773    REAL(r_std), PARAMETER :: TS = 303.0          !! Leaf temperature at standard condition
774                                                  !! (see Guenther .et. al, 1995, eq(10)) (K)
775    REAL(r_std), PARAMETER :: TM = 314.0          !! Leaf temperature (see Guenther .et. al, 1995, eq(10)) (K)
776
777    REAL(r_std), PARAMETER :: alpha_ = 0.0027     !! Empirical coeffcient (see Guenther .et. al, 1995, eq(9)) (unitless)
778    REAL(r_std), PARAMETER :: CL1 = 1.066         !! Empirical coeffcient (see Guenther .et. al, 1995, eq(9)) (unitless)
779    REAL(r_std), PARAMETER :: beta = 0.09         !! Empirical coeffcient (see Guenther .et. al, 1995, eq(11)) (K^{-1})
780    REAL(r_std), PARAMETER :: lai_threshold = 11. !! Lai threshold for the calculation of scattered radiation
781                                                  !! based on Guenther .et. al (1995) (m^2.m^{-2})
782
783                                             
784    ! Biogenic emissions
785    REAL(r_std),DIMENSION(kjpindex)          :: PAR                !! Photosynthetic active radiation, half of swdown
786                                                                   !! @tex ($\mu mol photons. m^{-2} s^{-1}$) @endtex
787    REAL(r_std),DIMENSION(kjpindex,nvm)      :: PARsun             !! PAR received by sun leaves
788                                                                   !! @tex ($\mu mol m^{-2} s^{-1}$) @endtex
789    REAL(r_std),DIMENSION(kjpindex,nvm)      :: PARsh              !! PAR received by shaded leaves
790                                                                   !! @tex ($\mu mol m^{-2} s^{-1}$) @endtex
791    REAL(r_std),DIMENSION(kjpindex,nvm)      :: laisun             !! Leaf area index of Sun leaves (m^2.m^{-2})
792    REAL(r_std),DIMENSION(kjpindex,nvm)      :: laish              !! Leaf area index of Shaded leaves (m^2.m^{-2})
793
794    CHARACTER(LEN=14)                        :: tleafsun_name      !! To store variables names for I/O
795    CHARACTER(LEN=13)                        :: tleafsh_name       !! To store variables names for I/O
796    REAL(r_std), DIMENSION(kjpindex,nlai+1)  :: Tleafsun_temp      !! temporary sunlit leaf temperature matrix for writing
797    REAL(r_std), DIMENSION(kjpindex,nlai+1)  :: Tleafsh_temp       !! temporary shade leaf temperature matrix for writing
798    REAL(r_std),DIMENSION(kjpindex)          :: Fdf                !! Fraction of the radiation which is diffused (0-1, unitless)
799    REAL(r_std),DIMENSION(kjpindex,nlai+1)   :: PARsuntab          !! PAR received by sun leaves over LAI layers
800                                                                   !! @tex ($\mu mol m^{-2} s^{-1}$) @endtex
801    REAL(r_std),DIMENSION(kjpindex,nlai+1)   :: PARshtab           !! PAR received by shaded leaves over LAI layers
802                                                                   !! @tex ($\mu mol m^{-2} s^{-1}$) @endtex
803    REAL(r_std),DIMENSION(kjpindex)          :: PARdf              !! Diffuse photosynthetic active radiation
804                                                                   !! @tex ($\mu mol m^{-2} s^{-1}$) @endtex
805    REAL(r_std),DIMENSION(kjpindex)          :: PARdr              !! Direct photosynthetic active radiation
806                                                                   !! @tex ($\mu mol m^{-2} s^{-1}$) @endtex
807    REAL(r_std),DIMENSION(kjpindex)          :: Trans              !! Atmospheric Transmissivity (unitless)
808
809!_ ================================================================================================================================
810        fco2 = 0.
811        fco2_wshort = 0.
812        fco2_wlong = 0.
813        Fdf(:) = 0.
814        PAR(:) = 0. 
815        PARsun(:,:) = 0. 
816        PARsh(:,:) = 0. 
817        laisun(:,:) = 0. 
818        laish(:,:) = 0. 
819        CRF(:,:) = 0.               
820        Trans(:) = 0.           
821        PARdf(:) = 0.           
822        PARdr(:) = 0.           
823        PARsuntab(:,:) = 0.       
824        PARshtab(:,:) = 0.       
825   
826
827    !! 1. Canopy radiative transfer model
828
829    !! Canopy radiative transfer model: takes into account light extinction through canopy
830    !! First need to calculate diffuse and direct radiation
831    !! Based on Andrew Friend radiative model (Global Ecology & Biogeography, 2001)
832    !! And Spitters et al. (Agricultural and Forest Meteorology, 1986)
833       
834    IF ( ok_radcanopy ) THEN
835
836       DO ji = 1, kjpindex
837          IF (coszang(ji) .GT. zero) THEN
838             !! 1.1 Extra-terrestrial solar irradiance at a plan parallel to Earh's surface
839             So(ji) = Sct*( un + 0.033*COS(360.*pi/180.*julian_diff/365.))*coszang(ji)
840             !! 1.2 Atmospheric transmissivity
841             Trans(ji) = swdown(ji)/So(ji)
842             !! 1.3 Empirical calculation of fraction diffuse from transmission based on Spitters et al. (1986)
843             Rfrac(ji) = 0.847 - 1.61*coszang(ji) + 1.04*(coszang(ji)**2.)
844             Kfrac(ji) = (1.47 - Rfrac(ji))/1.66     
845             IF (Trans(ji) .LE. 0.22) THEN
846                Fdf(ji) = un
847             ELSE IF (Trans(ji) .LE. 0.35) THEN
848                Fdf(ji) = un - 6.4*((Trans(ji) - 0.22)**2.) 
849             ELSE IF (Trans(ji) .LE. Kfrac(ji)) THEN
850                Fdf(ji) = 1.47 - 1.66*Trans(ji)
851             ELSE
852                Fdf(ji) = Rfrac(ji)
853             END IF
854             !! 1.4 Direct and diffuse sw radiation in W.m^{-2}
855             swdf(ji) = swdown(ji)*Fdf(ji)
856             swdr(ji) = swdown(ji)*(un-Fdf(ji))
857          ELSE
858             swdf(ji) = zero
859             swdr(ji) = zero
860          END IF
861
862          !! 1.5 PAR diffuse and direct in umol/m^2/s
863          PARdf(ji) = swdf(ji) * W_to_mol * RG_to_PAR
864          PARdr(ji) = swdr(ji) * W_to_mol * RG_to_PAR 
865       END DO
866
867       !! 1.6 Calculation of lai, parscat, parsh and parsun, laisun and laish !!?? define the terms
868       !! Based on Guenther et al. (JGR, 1995) and Norman (1982)
869       ! One-layer canopy model or multi-layer canopy model
870       IF (ok_multilayer) THEN 
871
872
873          ! Calculation PER LAYER
874          DO jl = nlai+1, 1, -1
875            laitab(jl) = laimax*(EXP(lai_level_depth*(jl-1) - un)/(EXP(lai_level_depth*nlai) - un))
876
877         !introduction of the Spitter way to calculate radiation over the levels !!!!!!!
878             DO ji = 1, kjpindex
879                kdf = cluster*0.8*SQRT((1 - sigma))
880                IF (ABS(ACOS(coszang(ji))) .LT. pi/2. .AND. coszang(ji) .NE. zero) THEN
881                   ! Coefficients calculation:
882                   kbl = cluster*0.5/ABS(coszang(ji))
883                   rho = ((1-SQRT((1 - sigma)))/(1+SQRT((1 - sigma))))*(2/(1+1.6*ABS(coszang(ji))))
884
885                   PARnotscat(ji,jl) = (1 - sigma)*PARdr(ji)*kbl*EXP(-SQRT((1 - sigma))*kbl*laitab(jl))
886                   PARabsdir(ji,jl) = (1 - rho)*SQRT((1 - sigma))*PARdr(ji)*kbl*EXP(-SQRT((1 - sigma))*kbl*laitab(jl))
887                   PARabsdiff(ji,jl) = (1 - rho)*PARdf(ji)*kdf*EXP(-kdf*laitab(jl))
888                   PARshtab(ji,jl) = (PARabsdiff(ji,jl) + (PARabsdir(ji,jl) - PARnotscat(ji,jl)))/(1 - sigma)
889                   PARsuntab(ji,jl) = PARshtab(ji,jl) + (1-sigma)*kbl*PARdr(ji)/(1 - sigma) 
890
891                   !correction using the equation (4) in Bodin et al 2012 and (7) or (8) in Spitter et al 1986
892                   !using the extinction coefficient kbl = 0.5/coszang and not only 0.5
893                   IF (jl .NE. nlai+1) THEN
894                      laisuntabdep(ji,jl) =(laitab(jl+1)-laitab(jl))*EXP(-kbl*laitab(jl))
895                      laishtabdep(ji,jl) =(laitab(jl+1)-laitab(jl))*(1.-EXP(-kbl*laitab(jl)))
896                   ENDIF
897
898                ELSE
899
900                   PARshtab(ji,jl) = PARdf(ji)*kdf*EXP(-kdf*laitab(jl))
901                   PARsuntab(ji,jl) = 0.
902
903                   IF (jl .NE. nlai+1) THEN
904                      laisuntabdep(ji,jl) = 0.
905                      laishtabdep(ji,jl) = laitab(jl+1)-laitab(jl)
906
907                   ENDIF 
908                END IF
909             END DO
910          END DO
911
912
913
914       ! introduction of the Spitter way to calculate radiation over the levels !!!!!!!
915       ELSE
916          ! Calculation FOR one layer
917          DO jv = 1, nvm
918             DO ji = 1, kjpindex
919                IF (lai(ji,jv) .LE. lai_threshold) THEN
920                   PARscat(ji,jv) = 0.07*PARdr(ji)*(1.1 - 0.1*lai(ji,jv))*exp(-coszang(ji))
921                ELSE
922                   PARscat(ji,jv) = zero
923                END IF
924
925                IF (coszang(ji) .NE. zero ) THEN
926                   PARsh(ji,jv) = PARdf(ji)*exp(-0.5*((lai(ji,jv))**0.7)) + PARscat(ji,jv)
927                   PARsun(ji,jv) = PARdr(ji)*COS(60.*pi/180.)/coszang(ji) + PARsh(ji,jv)
928                ELSE
929                   PARsh(ji,jv) = PARdf(ji)*exp(-0.5*(lai(ji,jv)**0.7)) + PARscat(ji,jv)
930                   PARsun(ji,jv) = zero 
931                END IF
932                IF (ABS(ACOS(coszang(ji))) .LT. pi/2. .AND. ABS(coszang(ji)) .GT. min_sechiba) THEN 
933                   ! calculation as in Lathiere (2005) = with correction removing lai in Guenther (1995)
934                   laisun(ji,jv) = (un - exp(-0.5*lai(ji,jv)/(coszang(ji))))*coszang(ji)/COS(60.*pi/180.)
935                   laish(ji,jv) = lai(ji,jv) - laisun(ji,jv)
936                ELSE
937                   laisun(ji,jv) = zero
938                   laish(ji,jv) = lai(ji,jv)
939                END IF
940             END DO
941          END DO
942       ENDIF
943    END IF
944
945
946    !! 2. Calculation of non-PFT dependant parameters used for VOC emissions
947    DO ji = 1, kjpindex ! (loop over # pixels)
948       !! 2.1 Calculation of Tleaf (based on Lathiere, 2005)
949
950
951       tleaf(ji) = temp_air(ji)
952
953       !! 2.2 Isoprene emission dependency - with no PARsun/PARshaded partitioning - Guenther et al. (1995) and Lathiere (2005)
954       !> @codeinc $$?? ecrire les equation en latex ?
955       exp_1(ji) = exp( (CT1 * ( tleaf(ji) - TS )) / (RR*TS*tleaf(ji)) )
956       exp_2(ji) = exp( (CT2 *( tleaf(ji) - TM )) / (RR*TS*tleaf(ji)) )
957       PAR(ji)   = swdown(ji) * W_to_mol * RG_to_PAR        ! from W/m^2 to umol photons/m^2/s and half of sw for PAR
958       Ct_iso(ji)    = exp_1(ji)/(un + exp_2(ji))            ! temperature dependance 
959       Cl_iso(ji)    = alpha_*CL1*PAR(ji)/sqrt(un + (alpha_**2) * (PAR(ji)**2) ) ! light dependance
960       !> @endcodeinc
961 
962       !! 2.3 Monoterpene and oxy VOB emission dependency to Temperature
963       !!     light independant fraction
964       !> @codeinc
965       !Ct_mono(ji) = exp(beta*(tleaf(ji) - TS))  ! Old method
966       Ct_mono(ji) = exp(beta_mono*(tleaf(ji) - TS))
967       Ct_sesq(ji) = exp(beta_sesq*(tleaf(ji) - TS))
968       Ct_meth(ji) = exp(beta_meth*(tleaf(ji) - TS))
969       Ct_acet(ji) = exp(beta_acet*(tleaf(ji) - TS))
970       Ct_oxyVOC(ji) = exp(beta_oxyVOC*(tleaf(ji) - TS))     
971       !> @endcodeinc
972       !! 2.4 MBO biogenic emissions dependency, only from PFT7 and PFT4 for location of vegetation emitter
973       ! but in fact MBO fluxes only in America (ponderosa and lodgepole pines only found in these areas)
974       !> @codeinc
975       Xvar(ji) = ((un/312.3) - (un/tleaf(ji)))/RR
976       !> @endcodeinc
977       !! 2.4.1 temperature dependency
978       !> @codeinc
979       Ct_MBO(ji)    = (1.52*209000.0*exp(67000.0*Xvar(ji)))/(209000.0 - 67000.0*(un - exp(209000.0*Xvar(ji))))
980       !> @endcodeinc
981       !! 2.4.2 light dependency
982       Cl_MBO(ji)    = (0.0011*1.44*PAR(ji))/(sqrt(un + (0.0011**2)*(PAR(ji)**2)))
983       !! 2.5 NO biogenic emissions given in ngN/m^2/s, emission factor in ngN/m^2/s too
984       !! calculation of temperature used for NO soil emissions
985       t_no(ji) = ptnlev1(ji) - ZeroCelsius  !!temp must be in celsius to calculate no emissions
986       !! 2.6 calculation of non-woody veget_max fraction
987       IF (ok_cropsfertil_NOx) THEN
988          veget_max_nowoody(ji) = zero
989          DO jv = 1,nvm
990             IF ( (jv /= ibare_sechiba) .AND. .NOT.(is_tree(jv)) ) THEN
991                veget_max_nowoody(ji) = veget_max_nowoody(ji) + veget_max(ji,jv)
992             ENDIF
993          ENDDO
994       END IF
995    END DO ! (loop over # pixels)
996
997    !! 2bis. Calculation of CO2 function for inhibition effect on isoprene
998    ! 2 approaches can be used: either Possell et al. (2005) or Wilkinson et al. (2006)
999
1000!! 19/04/2010 and then implemented in version revised by Nicolas Vuichard, 08042014
1001!! Impact of atmospheric CO2 on isoprene emissions
1002!! Can be activated or not
1003!! If considered, can use either Possell 2005 or Wilkinson 2009 parameterisation
1004!! This is used to rescale the emission factor, considered to be measured at 350 ppm of CO2
1005!! to the CO2 conditions of the run
1006
1007IF ( ok_co2bvoc_poss ) THEN
1008   WRITE(numout,*) 'CO2 impact on isoprene: Possell calculation'
1009
1010   !! Possell function needs to be normalized (experiments at 400 ppm and EF before 1995)
1011   !! Normalized at 350 ppm
1012   fco2_ctrl = (-0.0123+(441.4795/350.)+(-1282.65/(350.**2)))
1013
1014   !! 2 tests: using the canopy (atmospheric) CO2 'ccanopy'
1015   !! or the intercellular CO2 over nlai 'cim'
1016   !! with cim = ccanopy*0.667
1017   !! in the end I go for ccanopy for the Possell function since too much differences
1018   !! when using cim and also the function has been derived based on atmospheric CO2
1019   DO ji = 1, kjpindex
1020
1021      fco2(ji,:) = (-0.0123+(441.4795/ccanopy(ji))+(-1282.65/(ccanopy(ji)*ccanopy(ji))))/fco2_ctrl
1022
1023   END DO
1024ELSE IF ( ok_co2bvoc_wilk ) THEN
1025   WRITE(numout,*) 'CO2 impact on isoprene: Wilkinson calculation'
1026
1027   !! In the Wilkinson function, 2 impacts are considered:
1028   !! -short-term impact for CO2 variation during a single day (seconds/minutes)
1029   !! -long-term impact for CO2 variation during leaf-growth (weeks/month)
1030
1031
1032   !! Long-term parameters
1033   !! Constant
1034   Ismax_long = 1.344
1035   h_long = 1.4614
1036   Cstar_long = 585.
1037   !! Short-term parameters
1038   !! They have to be calculated based on atmospheric CO2
1039   !! 10/05/2010
1040   !! For atmospheric CO2 lower than 400 ppm or higher than 1200 ppm
1041   !! (min and max CO2 level tested for short-term effect in Wilkinson et al. 2009)
1042   !! we use the parameters calculated at 400/1200 ppm. For intermediate CO2 concentration,
1043   !! parameters are calculated.
1044   !! Linear interpolation
1045
1046   DO ji = 1, kjpindex
1047
1048      IF (ccanopy(ji) .LE. 400.) THEN
1049
1050         Ismax_short(ji) = 1.072
1051         h_short(ji) = 1.7
1052         Cstar_short(ji) = 1218.
1053
1054      ELSE IF (ccanopy(ji) .EQ. 600.) THEN
1055
1056         Ismax_short(ji) = 1.036
1057         h_short(ji) = 2.0125
1058         Cstar_short(ji) = 1150.
1059
1060      ELSE IF (ccanopy(ji) .EQ. 800.) THEN
1061
1062         Ismax_short(ji) = 1.046
1063         h_short(ji) = 1.5380
1064         Cstar_short(ji) = 2025.
1065
1066      ELSE IF (ccanopy(ji) .GE. 1200.) THEN
1067
1068         Ismax_short(ji) = 1.014
1069         h_short(ji) = 2.8610
1070         Cstar_short(ji) = 1525.
1071
1072
1073      ELSE IF ((ccanopy(ji) .GT. 400.) .AND. (ccanopy(ji) .LT. 600.)) THEN
1074
1075         Ismax_short(ji) = 1.036 + (ccanopy(ji)-600.)*(1.036-1.072)/(600.-400.)
1076         h_short(ji) = 2.0125 + (ccanopy(ji)-600.)*(2.0125-1.7)/(600.-400.)
1077         Cstar_short(ji) =  1150. + (ccanopy(ji)-600.)*(1150.-1218.)/(600.-400.)
1078
1079      ELSE IF ((ccanopy(ji) .GT. 600.) .AND. (ccanopy(ji) .LT. 800.)) THEN
1080
1081         Ismax_short(ji) = 1.046 + (ccanopy(ji)-800.)*(1.046-1.036)/(800.-600.)
1082         h_short(ji) = 1.5380 + (ccanopy(ji)-800.)*(1.5380-2.0125)/(800.-600.)
1083         Cstar_short(ji) = 2025. + (ccanopy(ji)-800.)*(2025.-1150.)/(800.-600.)
1084
1085      ELSE IF ((ccanopy(ji) .GT. 800.) .AND. (ccanopy(ji) .LT. 1200.)) THEN
1086
1087        Ismax_short(ji) = 1.014 + (ccanopy(ji)-1200.)*(1.014-1.046)/(1200.-800.)
1088        h_short(ji) = 2.8610 + (ccanopy(ji)-1200.)*(2.8610-1.5380)/(1200.-800.)
1089        Cstar_short(ji) = 1525. + (ccanopy(ji)-1200.)*(1525.-2025.)/(1200.-800.)
1090
1091
1092      END IF
1093
1094   END DO
1095
1096   WRITE(numout,*) '***Wilkinson BVOC-CO2 function: parameters***'
1097   WRITE(numout,*) 'Ismax_long: ', Ismax_long
1098   WRITE(numout,*) 'h_long: ', h_long
1099   WRITE(numout,*) 'Cstar_long: ', Cstar_long
1100   WRITE(numout,*) 'Ismax_short: ', MAXVAL(Ismax_short(:)) , MINVAL(Ismax_short(:))
1101   WRITE(numout,*) 'h_short: ', MAXVAL(h_short(:)) , MINVAL(h_short(:))
1102   WRITE(numout,*) 'Cstar_short: ', MAXVAL(Cstar_short(:)) , MINVAL(Cstar_short(:))
1103   WRITE(numout,*) '******'
1104
1105   DO ji = 1, kjpindex
1106      fco2_wlong(ji) = (Ismax_long-((Ismax_long*((0.667*ccanopy(ji))**h_long))/&
1107                     & ((Cstar_long**h_long)+(0.667*ccanopy(ji))**h_long)))/1.06566
1108      DO jv = 1, nvm
1109         fco2_wshort(ji,jv) = (Ismax_short(ji)-((Ismax_short(ji)*((cim(ji,jv))**h_short(ji)))/&
1110                            & ((Cstar_short(ji)**h_short(ji))+(cim(ji,jv))**h_short(ji))))/1.010803
1111      END DO
1112   END DO
1113
1114   DO ji = 1, kjpindex
1115      DO jv = 1, nvm
1116         fco2(ji,jv) = fco2_wshort(ji,jv)*fco2_wlong(ji)
1117      END DO
1118   END DO
1119
1120ELSE
1121      WRITE(numout,*) 'CO2 impact on isoprene not considered'
1122      fco2(:,:) = 1.
1123END IF
1124
1125
1126    !! 3. Calculation of PFT dependant parameters and
1127    ! Calculation of VOC emissions flux
1128
1129    Eff_age_iso(:,:) = zero
1130    Eff_age_meth(:,:) = zero
1131
1132
1133    DO jv = 1, nvm ! loop over the PDFs
1134       DO ji = 1, kjpindex ! loop over the grid cell
1135          ! 6-Calculation of Leaf Age Function (Lathiere 2005)
1136          IF ( ok_leafage ) THEN
1137             DO jf = 1, nleafages
1138                !> @codeinc
1139                Eff_age_iso(ji,jv) = Eff_age_iso(ji,jv) + frac_age(ji,jv,jf)*iso_activity(jf)
1140                Eff_age_meth(ji,jv) = Eff_age_meth(ji,jv) + frac_age(ji,jv,jf)*methanol_activity(jf)
1141                !> @endcodeinc
1142             END DO
1143             !> @codeinc
1144             Eff_age_VOC(ji,jv) = un
1145             !> @endcodeinc
1146          ELSE
1147             Eff_age_iso(ji,jv) = un
1148             Eff_age_meth(ji,jv) = un
1149             Eff_age_VOC(ji,jv) = un
1150          END IF
1151          !! 5. Calculation of foliar density
1152          IF ( sla(jv) .eq. zero ) THEN
1153             fol_dens(ji,jv) = zero
1154          ELSE
1155             ! 2 factor for conversion from gC to gDM
1156             fol_dens(ji,jv) = 2 * lai(ji,jv)/sla(jv)
1157          ENDIF
1158          !! 6. Calculation of VOC emissions from vegetation
1159          IF ( ok_radcanopy ) THEN
1160             ! if multi-layer canopy model
1161             IF (ok_multilayer) THEN
1162
1163                laisun(ji,jv) = zero
1164                laish(ji,jv) = zero
1165                GAMMA_iso_m  = zero
1166                flx_iso(ji,jv) = zero
1167                flx_mono(ji,jv) = zero
1168                flx_apinen(ji,jv) = zero
1169                flx_bpinen(ji,jv) = zero
1170                flx_limonen(ji,jv) = zero
1171                flx_myrcen(ji,jv) =  zero
1172                flx_sabinen(ji,jv) =  zero
1173                flx_camphen(ji,jv) = zero
1174                flx_3caren(ji,jv) = zero
1175                flx_tbocimen(ji,jv) = zero
1176                flx_othermono(ji,jv) = zero
1177                flx_sesquiter(ji,jv) =  zero
1178                flx_methanol(ji,jv) = zero
1179                flx_acetone(ji,jv) =  zero
1180                flx_acetal(ji,jv) = zero
1181                flx_formal(ji,jv) = zero
1182                flx_acetic(ji,jv) = zero
1183                flx_formic(ji,jv) = zero
1184                ! loop over the NLAI canopy layers
1185                DO jl = 1, nlai
1186                   IF ((laitab(jl) .LE. lai(ji,jv)).AND.(lai(ji,jv).NE.zero)) THEN
1187                      !sunlit vegetation
1188                      Clsun_iso_tab   = alpha_*CL1*PARsuntab(ji,jl)/sqrt(un + (alpha_**2) * (PARsuntab(ji,jl)**2) )
1189                      ! shaded vegetation
1190                      Clsh_iso_tab    = alpha_*CL1*PARshtab(ji,jl)/sqrt(un + (alpha_**2) * (PARshtab(ji,jl)**2) ) 
1191                      flx_iso(ji,jv) = flx_iso(ji,jv) + (laisuntabdep(ji,jl)*Clsun_iso_tab+ &
1192                           & laishtabdep(ji,jl)*Clsh_iso_tab)* &
1193                           & fol_dens(ji,jv)/lai(ji,jv)*Ct_iso(ji)*em_factor_isoprene(jv)* &
1194                           & Eff_age_iso(ji,jv)*fco2(ji,jv)*1e-9/one_hour
1195
1196                      GAMMA_iso_m = GAMMA_iso_m + (laisuntabdep(ji,jl)*Clsun_iso_tab+ &
1197                           & laishtabdep(ji,jl)*Clsh_iso_tab)* &
1198                           & fol_dens(ji,jv)/lai(ji,jv)*Ct_iso(ji)*1e-9/one_hour
1199
1200                      laisun(ji,jv) = laisun(ji,jv) + laisuntabdep(ji,jl)
1201                      laish(ji,jv)  = laish(ji,jv) + laishtabdep(ji,jl)
1202                   END IF
1203                END DO
1204
1205                !! 6.1 Calculation of monoterpene biogenic emissions
1206                flx_mono(ji,jv) = ((1-LDF_mono)*Ct_mono(ji)*1e-9/one_hour*fol_dens(ji,jv) + LDF_mono*GAMMA_iso_m)* &
1207                     & em_factor_monoterpene(jv)*Eff_age_VOC(ji,jv) 
1208                !! 6.12 Calculation of sesquiterpenes biogenic emission
1209                flx_sesquiter(ji,jv) = ((1-LDF_sesq)*Ct_sesq(ji)*1e-9/one_hour*fol_dens(ji,jv) +LDF_sesq*GAMMA_iso_m)* &
1210                     & em_factor_sesquiterp(jv)*Eff_age_VOC(ji,jv)
1211                !! 6.13 Calculation of methanol biogenic emissions
1212                flx_methanol(ji,jv) = ((1-LDF_meth)*Ct_meth(ji)*1e-9/one_hour*fol_dens(ji,jv) +LDF_meth*GAMMA_iso_m)* &
1213                     & em_factor_methanol(jv)*Eff_age_meth(ji,jv)
1214                !! 6.14 Calculation of acetone biogenic emissions
1215                flx_acetone(ji,jv) = ((1-LDF_acet)*Ct_acet(ji)*1e-9/one_hour*fol_dens(ji,jv) +LDF_acet*GAMMA_iso_m)* &
1216                     & em_factor_acetone(jv)*Eff_age_VOC(ji,jv)
1217                !! 6.14 Calculation of acetaldehyde biogenic emissions
1218                flx_acetal(ji,jv) = ((1-LDF_meth)*Ct_meth(ji)*1e-9/one_hour*fol_dens(ji,jv) +LDF_meth*GAMMA_iso_m)* &
1219                     & em_factor_acetal(jv)*Eff_age_VOC(ji,jv)
1220                !! 6.16 Calculation of formaldehyde biogenic emissions
1221                flx_formal(ji,jv) = ((1-LDF_meth)*Ct_meth(ji)*1e-9/one_hour*fol_dens(ji,jv) +LDF_meth*GAMMA_iso_m)* &
1222                     & em_factor_formal(jv)*Eff_age_VOC(ji,jv)
1223                !! 6.17 Calculation of acetic acid biogenic emissions
1224                flx_acetic(ji,jv) = ((1-LDF_meth)*Ct_meth(ji)*1e-9/one_hour*fol_dens(ji,jv) +LDF_meth*GAMMA_iso_m)* &
1225                     & em_factor_acetic(jv)*Eff_age_VOC(ji,jv)
1226                !! 6.18 Calculation of formic acid biogenic emissions
1227                flx_formic(ji,jv) = ((1-LDF_meth)*Ct_meth(ji)*1e-9/one_hour*fol_dens(ji,jv) +LDF_meth*GAMMA_iso_m)* &
1228                     & em_factor_formic(jv)*Eff_age_VOC(ji,jv)
1229
1230
1231                !! 6.3 Calculation of alfa pinene biogenic emission
1232                flx_apinen(ji,jv) = em_factor_apinene(jv)*flx_mono(ji,jv) 
1233                !! 6.4 Calculation of beta pinene biogenic emission
1234                flx_bpinen(ji,jv) = em_factor_bpinene(jv)*flx_mono(ji,jv) 
1235                !! 6.5 Calculation of limonene biogenic emission
1236                flx_limonen(ji,jv) = em_factor_limonene(jv)*flx_mono(ji,jv) 
1237                !! 6.6 Calculation of myrcene biogenic emission !!
1238                flx_myrcen(ji,jv) = em_factor_myrcene(jv)*flx_mono(ji,jv) 
1239                !! 6.7 Calculation of sabinene biogenic emission
1240                flx_sabinen(ji,jv) = em_factor_sabinene(jv)*flx_mono(ji,jv) 
1241                !! 6.8 Calculation of camphene biogenic emission
1242                flx_camphen(ji,jv) = em_factor_camphene(jv)*flx_mono(ji,jv) 
1243                !! 6.9 Calculation of 3-carene biogenic emission
1244                flx_3caren(ji,jv) = em_factor_3carene(jv)*flx_mono(ji,jv) 
1245                !! 6.10 Calculation of t-beta-ocimene biogenic emission
1246                flx_tbocimen(ji,jv) = em_factor_tbocimene(jv)*flx_mono(ji,jv) 
1247                !! 6.11 Calculation of other monoterpenes biogenic emission
1248                flx_othermono(ji,jv) = em_factor_othermonot(jv)*flx_mono(ji,jv) 
1249
1250                ! if mono-layer canopy model
1251             ELSE
1252                !sunlit vegetation
1253                Clsun_iso(ji,jv)   = alpha_*CL1*PARsun(ji,jv)/sqrt(un + (alpha_**2) * (PARsun(ji,jv)**2) )
1254                ! shaded vegetation     
1255                Clsh_iso(ji,jv)    = alpha_*CL1*PARsh(ji,jv)/sqrt(un + (alpha_**2) * (PARsh(ji,jv)**2) )       
1256                IF (lai(ji,jv) .NE. zero) THEN
1257                   !! 6.1 Calculation of isoprene biogenic emissions
1258                   GAMMA_iso = (laisun(ji,jv)*Clsun_iso(ji,jv) + laish(ji,jv)*Clsh_iso(ji,jv))/lai(ji,jv)*Ct_iso(ji)
1259                   flx_iso(ji,jv) = GAMMA_iso*fol_dens(ji,jv)*em_factor_isoprene(jv)*Eff_age_iso(ji,jv)*fco2(ji,jv)*1e-9/one_hour
1260                   !! 6.2 Calculation of monoterpene biogenic emissions
1261                   flx_mono(ji,jv) = ((1-LDF_mono)*Ct_mono(ji)+LDF_mono*GAMMA_iso)*fol_dens(ji,jv)* &
1262                        & em_factor_monoterpene(jv)*Eff_age_VOC(ji,jv)*1e-9/one_hour
1263                   !! 6.3 Calculation of alfa pinene biogenic emission
1264                   flx_apinen(ji,jv) = em_factor_apinene(jv)*flx_mono(ji,jv)
1265                   !! 6.4 Calculation of beta pinene biogenic emission
1266                   flx_bpinen(ji,jv) = em_factor_bpinene(jv)*flx_mono(ji,jv)
1267                   !! 6.5 Calculation of limonene biogenic emission
1268                   flx_limonen(ji,jv) = em_factor_limonene(jv)*flx_mono(ji,jv)
1269                   !! 6.6 Calculation of myrcene biogenic emission
1270                   flx_myrcen(ji,jv) = em_factor_myrcene(jv)*flx_mono(ji,jv)
1271                   !! 6.7 Calculation of sabinene biogenic emission
1272                   flx_sabinen(ji,jv) = em_factor_sabinene(jv)*flx_mono(ji,jv)
1273                   !! 6.8 Calculation of camphene biogenic emission
1274                   flx_camphen(ji,jv) = em_factor_camphene(jv)*flx_mono(ji,jv)
1275                   !! 6.9 Calculation of 3-carene biogenic emission
1276                   flx_3caren(ji,jv) = em_factor_3carene(jv)*flx_mono(ji,jv)
1277                   !! 6.10 Calculation of t-beta-ocimene biogenic emission
1278                   flx_tbocimen(ji,jv) = em_factor_tbocimene(jv)*flx_mono(ji,jv)
1279                   !! 6.11 Calculation of other monoterpenes biogenic emission
1280                   flx_othermono(ji,jv) = em_factor_othermonot(jv)*flx_mono(ji,jv)
1281                   !! 6.12 Calculation of sesquiterpenes biogenic emission
1282                   flx_sesquiter(ji,jv) = ((1-LDF_sesq)*Ct_sesq(ji)+LDF_sesq*GAMMA_iso)*fol_dens(ji,jv)* &
1283                        & em_factor_sesquiterp(jv)*Eff_age_VOC(ji,jv)*1e-9/one_hour
1284                   !! 6.13 Calculation of methanol biogenic emissions
1285                   flx_methanol(ji,jv) = ((1-LDF_meth)*Ct_meth(ji)+LDF_meth*GAMMA_iso)*fol_dens(ji,jv)* &
1286                        & em_factor_methanol(jv)*Eff_age_meth(ji,jv)*1e-9/one_hour
1287                   !! 6.14 Calculation of acetone biogenic emissions
1288                   flx_acetone(ji,jv) = ((1-LDF_acet)*Ct_acet(ji)+LDF_acet*GAMMA_iso)*fol_dens(ji,jv)* &
1289                        & em_factor_acetone(jv)*Eff_age_VOC(ji,jv)*1e-9/one_hour
1290                   !! 6.15 Calculation of acetaldehyde biogenic emissions
1291                   flx_acetal(ji,jv) = ((1-LDF_meth)*Ct_meth(ji)+LDF_meth*GAMMA_iso)*fol_dens(ji,jv)* &
1292                        & em_factor_acetal(jv)*Eff_age_VOC(ji,jv)*1e-9/one_hour
1293                   !! 6.16 Calculation of formaldehyde biogenic emissions
1294                   flx_formal(ji,jv) = ((1-LDF_meth)*Ct_meth(ji)+LDF_meth*GAMMA_iso)*fol_dens(ji,jv)* &
1295                        & em_factor_formal(jv)*Eff_age_VOC(ji,jv)*1e-9/one_hour
1296                   !! 6.17 Calculation of acetic acid biogenic emissions
1297                   flx_acetic(ji,jv) = ((1-LDF_meth)*Ct_meth(ji)+LDF_meth*GAMMA_iso)*fol_dens(ji,jv)* &
1298                        & em_factor_acetic(jv)*Eff_age_VOC(ji,jv)*1e-9/one_hour
1299                   !! 6.18 Calculation of formic acid biogenic emissions
1300                   flx_formic(ji,jv) = ((1-LDF_meth)*Ct_meth(ji)+LDF_meth*GAMMA_iso)*fol_dens(ji,jv)* &
1301                        & em_factor_formic(jv)*Eff_age_VOC(ji,jv)*1e-9/one_hour
1302 
1303                ELSE
1304                   !
1305                   flx_iso(ji,jv) = zero
1306                   flx_mono(ji,jv) = zero
1307                   flx_apinen(ji,jv) = zero 
1308                   flx_bpinen(ji,jv) = zero 
1309                   flx_limonen(ji,jv) = zero 
1310                   flx_myrcen(ji,jv) =  zero
1311                   flx_sabinen(ji,jv) =  zero 
1312                   flx_camphen(ji,jv) = zero 
1313                   flx_3caren(ji,jv) = zero 
1314                   flx_tbocimen(ji,jv) = zero
1315                   flx_othermono(ji,jv) = zero 
1316                   flx_sesquiter(ji,jv) =  zero 
1317                   flx_methanol(ji,jv) = zero
1318                   flx_acetone(ji,jv) =  zero 
1319                   flx_acetal(ji,jv) = zero
1320                   flx_formal(ji,jv) = zero 
1321                   flx_acetic(ji,jv) = zero 
1322                   flx_formic(ji,jv) = zero 
1323                END IF
1324             END IF
1325             ! if no light extinction due to vegetation 
1326          ELSE
1327             !! Isoprene emissions - general equation
1328             flx_iso(ji,jv) = fol_dens(ji,jv)*Ct_iso(ji)*Cl_iso(ji)*Eff_age_iso(ji,jv)*fco2(ji,jv)* &
1329                  em_factor_isoprene(jv)*1e-9/one_hour
1330             !! 6.2 Calculation of monoterpene biogenic emissions
1331             Ylt_mono(ji) = ((1-LDF_mono)*Ct_mono(ji)+LDF_mono*Ct_iso(ji)*Cl_iso(ji)) 
1332             flx_mono(ji,jv) = fol_dens(ji,jv)*em_factor_monoterpene(jv)*Ylt_mono(ji)*Eff_age_VOC(ji,jv)*&
1333                  1e-9/one_hour
1334             !! 6.3 Calculation of alfa pinene biogenic emission
1335             flx_apinen(ji,jv) = em_factor_apinene(jv)*flx_mono(ji,jv) 
1336             !! 6.4 Calculation of beta pinene biogenic emission
1337             flx_bpinen(ji,jv) = em_factor_bpinene(jv)*flx_mono(ji,jv)                       
1338             !! 6.5 Calculation of limonene biogenic emission
1339             flx_limonen(ji,jv) = em_factor_limonene(jv)*flx_mono(ji,jv)                     
1340             !! 6.6 Calculation of myrcene biogenic emission
1341             flx_myrcen(ji,jv) = em_factor_myrcene(jv)*flx_mono(ji,jv)                       
1342             !! 6.7 Calculation of sabinene biogenic emission
1343             flx_sabinen(ji,jv) = em_factor_sabinene(jv)*flx_mono(ji,jv)           
1344             !! 6.8 Calculation of camphene biogenic emission
1345             flx_camphen(ji,jv) = em_factor_camphene(jv)*flx_mono(ji,jv)
1346             !! 6.9 Calculation of 3-carene biogenic emission
1347             flx_3caren(ji,jv) = em_factor_3carene(jv)*flx_mono(ji,jv)                       
1348             !! 6.10 Calculation of t-beta-ocimene biogenic emission
1349             flx_tbocimen(ji,jv) = em_factor_tbocimene(jv)*flx_mono(ji,jv)                     
1350             !! 6.11 Calculation of other monoterpenes biogenic emission
1351             flx_othermono(ji,jv) = em_factor_othermonot(jv)*flx_mono(ji,jv)                   
1352             !! 6.12 Calculation of sesquiterpenes biogenic emission
1353             Ylt_sesq(ji) = ((1-LDF_sesq)*Ct_sesq(ji)+LDF_sesq*Ct_iso(ji)*Cl_iso(ji))
1354             flx_sesquiter(ji,jv) = fol_dens(ji,jv)*em_factor_sesquiterp(jv)*Ylt_sesq(ji)*Eff_age_VOC(ji,jv)*1e-9/one_hour   
1355             !! 6.16 Calculation of methanol biogenic emissions
1356             Ylt_meth(ji) = ((1-LDF_meth)*Ct_meth(ji)+LDF_meth*Ct_iso(ji)*Cl_iso(ji))
1357             flx_methanol(ji,jv) = fol_dens(ji,jv)*em_factor_methanol(jv)*Ylt_meth(ji)*Eff_age_meth(ji,jv)*1e-9/one_hour
1358             !! 6.17 Calculation of acetone biogenic emissions
1359             Ylt_acet(ji) = ((1-LDF_acet)*Ct_acet(ji)+LDF_acet*Ct_iso(ji)*Cl_iso(ji))
1360             flx_acetone(ji,jv) = fol_dens(ji,jv)*em_factor_acetone(jv)*Ylt_acet(ji)*Eff_age_VOC(ji,jv)*1e-9/one_hour
1361             !! 6.18 Calculation of acetaldehyde biogenic emissions
1362             flx_acetal(ji,jv) = fol_dens(ji,jv)*em_factor_acetal(jv)*Ylt_meth(ji)*Eff_age_VOC(ji,jv)*1e-9/one_hour
1363             !! 6.19 Calculation of formaldehyde biogenic emissions
1364             flx_formal(ji,jv) = fol_dens(ji,jv)*em_factor_formal(jv)*Ylt_meth(ji)*Eff_age_VOC(ji,jv)*1e-9/one_hour
1365             !! 6.20 Calculation of acetic acid biogenic emissions
1366             flx_acetic(ji,jv) = fol_dens(ji,jv)*em_factor_acetic(jv)*Ylt_meth(ji)*Eff_age_VOC(ji,jv)*1e-9/one_hour
1367             !! 6.21 Calculation of formic acid biogenic emissions
1368             flx_formic(ji,jv) = fol_dens(ji,jv)*em_factor_formic(jv)*Ylt_meth(ji)*Eff_age_VOC(ji,jv)*1e-9/one_hour
1369
1370          END IF
1371
1372          !! 6.22 Calculation of ORVOC biogenic emissions
1373          !! Other Reactive Volatile Organic Compounds
1374          !> @codeinc
1375          flx_ORVOC(ji,jv) = fol_dens(ji,jv)*em_factor_ORVOC(jv)*Ct_mono(ji)*Eff_age_VOC(ji,jv)*1e-9/one_hour
1376          !> @endcodeinc
1377          !! 6.4 Calculation of OVOC biogenic emissions
1378          !! Other Volatile Organic Compounds
1379          flx_OVOC(ji,jv) = fol_dens(ji,jv)*em_factor_OVOC(jv)*Ct_mono(ji)*Eff_age_VOC(ji,jv)*1e-9/one_hour
1380          !! 6.5 Calculation of MBO biogenic emissions
1381          !! 2-Methyl-3-Buten-2-ol
1382          IF(lalo(ji,1) .GE. 20. .AND. lalo(ji,2) .LE. -100) THEN
1383             flx_MBO(ji,jv) = fol_dens(ji,jv)*em_factor_MBO(jv)*Ct_MBO(ji)*Cl_MBO(ji)*Eff_age_VOC(ji,jv)*1e-9/one_hour
1384          ELSE
1385             flx_MBO(ji,jv) = zero
1386          END IF
1387       END DO
1388
1389    END DO
1390
1391
1392    !! 7. Calculation of NOx emissions from soils
1393    ! Based on Yienger & Levy (1995) and Lathiere (2005, chapter 3)
1394    DO ji = 1, kjpindex
1395       !! 7.1 Precipitation-related pulse function
1396       IF (ok_pulse_NOx) THEN
1397          ! if we are during a period where pulses are not allowed
1398          IF (ok_siesta(ji)) THEN
1399             ! if this period is not over
1400             IF (FLOOR(siestaday(ji)) .LE. siestalim(ji)) THEN
1401                siestaday(ji) = siestaday(ji) + (dt_sechiba/one_day)
1402                ! if this period is over
1403             ELSE
1404                ok_siesta(ji) = .FALSE.
1405                siestaday(ji) = zero
1406             END IF
1407          END IF
1408          ! if we are during a period where pulses are allowed
1409          IF ((.NOT. ok_siesta(ji)) .AND. (.NOT. allow_pulse(ji))) THEN
1410             IF (humrel(ji,1) .LT. 0.15) THEN
1411                ! if precip exceeds 1 mm/day over one time step => a pulse occurs
1412                IF(precip_rain(ji)/nbre_precip .GE. un/(one_day/dt_sechiba)) THEN
1413                   ! if precip is up to 5 mm/day => pulse length is 3 days
1414                   IF (precip_rain(ji)/nbre_precip .LT. 5./(one_day/dt_sechiba)) THEN
1415                      pulselim(ji) = 3.
1416                      ! if precip is up to 15 mm/day => pulse length is 7 days
1417                   ELSE IF (precip_rain(ji)/nbre_precip .LT. 15./(one_day/dt_sechiba)) THEN
1418                      pulselim(ji) = 7.
1419                      ! if precip is upper than 15 mm/day => pulse length is 14 days
1420                   ELSE IF (precip_rain(ji)/nbre_precip .GE. 15./(one_day/dt_sechiba)) THEN
1421                      pulselim(ji) = 14.
1422                   END IF
1423                   allow_pulse(ji)=.TRUE.
1424                   pulseday(ji) = un
1425                END IF
1426             END IF
1427          END IF
1428          ! if we were during a pulse period
1429          IF (allow_pulse(ji)) THEN
1430             ! if we are still during the pulse period
1431             ! 16/06/2010 We assume a (pulselim-1) days for the pulse length (NVui+Jlath)
1432             IF(FLOOR(pulseday(ji)) .LT. pulselim(ji)) THEN
1433                ! calculation of the pulse function
1434                IF (pulselim(ji).EQ.3) THEN
1435                   pulse(ji) = 11.19*exp(-0.805*pulseday(ji))
1436                ELSE IF (pulselim(ji).EQ.7) THEN
1437                   pulse(ji) = 14.68*exp(-0.384*pulseday(ji))
1438                ELSE IF (pulselim(ji).EQ.14) THEN
1439                   pulse(ji) = 18.46*exp(-0.208*pulseday(ji))
1440                END IF
1441                pulseday(ji) = pulseday(ji) + (dt_sechiba/one_day)
1442                ! if the pulse period is over
1443             ELSE
1444                ! pulse function is set to 1
1445                pulse(ji) = un
1446                allow_pulse(ji) = .FALSE.
1447                siestaday(ji) = un
1448                siestalim(ji) = pulselim(ji)
1449                ok_siesta(ji) = .TRUE. 
1450             END IF
1451          END IF
1452          ! no precipitation-related pulse function
1453       ELSE
1454          pulse(ji) = un
1455       END IF
1456    END DO
1457
1458    !! 7.2 Calculation of NO basal emissions including pulse effect
1459    DO jv = 1, nvm
1460       DO ji = 1, kjpindex
1461          !Tropical forests
1462          IF ( is_tropical(jv) .AND. is_evergreen(jv) ) THEN
1463             ! Wet soils
1464             IF (humrel(ji,1) .GT. 0.3) THEN
1465                flx_no_soil(ji,jv) = 2.6*pulse(ji)
1466                ! Dry soils
1467             ELSE
1468                flx_no_soil(ji,jv) = 8.6*pulse(ji)
1469             END IF
1470             !Else If agricultural lands OR Wet soils
1471          ELSE IF ( ( .NOT.(natural(jv)) ) .OR. ( humrel(ji,1) .GT. 0.3 ) ) THEN
1472             ! Calculation of NO emissions depending of Temperature
1473             IF (t_no(ji) .LT. zero) THEN
1474                flx_no_soil(ji,jv) = zero
1475             ELSE IF (t_no(ji) .LE. 10.) THEN
1476                flx_no_soil(ji,jv) = 0.28*em_factor_no_wet(jv)*t_no(ji)*pulse(ji)
1477             ELSE IF (t_no(ji) .LE. 30.) THEN
1478                flx_no_soil(ji,jv) = em_factor_no_wet(jv)*exp(0.103*t_no(ji))*pulse(ji)
1479             ELSE
1480                flx_no_soil(ji,jv) = 21.97*em_factor_no_wet(jv)*pulse(ji)
1481             END IF
1482             !Else if Temp negative
1483          ELSE IF (t_no(ji) .LT. zero) THEN
1484             flx_no_soil(ji,jv) = zero
1485             !Else if Temp <= 30
1486          ELSE IF (t_no(ji) .LE. 30.) THEN
1487             flx_no_soil(ji,jv) = (em_factor_no_dry(jv)*t_no(ji))/30.*pulse(ji)
1488          ELSE
1489             flx_no_soil(ji,jv) = em_factor_no_dry(jv)*pulse(ji)
1490          END IF
1491
1492          !! 7.3 IF ACTIVATED (ok_bbgfertil_NOx = TRUE) calculation of NOx soil emission increase due to biomass burning
1493          ! Calculation of Biomass-Burning-induced NOx emissions (Lathiere, 2005)
1494          ! => NOx emissions 3-fold increase
1495          IF (ok_bbgfertil_NOx) THEN
1496             IF ( natural(jv) ) THEN
1497                ! North Tropical zone from May to June
1498                IF ((lalo(ji,1) .LE. 30. .AND. lalo(ji,1) .GE. zero).AND. &
1499                     (julian_diff .GE. 121. .AND. julian_diff .LE. 181).AND.(flx_co2_bbg_year(ji) .GT. 0.1)) THEN
1500                   flx_no_soil(ji,jv) = flx_no_soil(ji,jv)*3.
1501                   ! South Tropical zone from November to December
1502                ELSE IF ((lalo(ji,1) .GE. -30. .AND. lalo(ji,1) .LT. zero).AND.(julian_diff .GE. 305.).AND. & 
1503                        (flx_co2_bbg_year(ji) .GT. 0.1)) THEN
1504                   flx_no_soil(ji,jv) = flx_no_soil(ji,jv)*3.
1505                END IF
1506             END IF
1507          END IF
1508
1509          !! 7.4 IF ACTIVATED (ok_cropsfertil_NOx = TRUE) calculation of NOx soil emission increase due to fertilizer use
1510          ! Calculation of N-fertiliser-induced NOx emissions
1511          flx_fertil_no(ji,jv) = zero
1512          IF (ok_cropsfertil_NOx) THEN
1513             IF (veget_max_nowoody(ji) .NE. zero) THEN
1514                ! Non-agricultural lands
1515                IF ( (jv == ibare_sechiba) .OR. is_tree(jv) ) THEN
1516                   N_qt_WRICE_pft(ji,jv) = zero
1517                   N_qt_OTHER_pft(ji,jv) = zero
1518                ! Grasslands or Croplands
1519                ELSE
1520                   N_qt_WRICE_pft(ji,jv) = N_qt_WRICE_year(ji)*veget_max(ji,jv)/veget_max_nowoody(ji)
1521                   N_qt_OTHER_pft(ji,jv) = N_qt_OTHER_year(ji)*veget_max(ji,jv)/veget_max_nowoody(ji)
1522                END IF
1523             ELSE
1524                N_qt_WRICE_pft(ji,jv) = zero
1525                N_qt_OTHER_pft(ji,jv) = zero
1526             END IF
1527
1528             ! North temperate regions from May to August
1529             ! OR South Temperate regions from November to February
1530             IF (((lalo(ji,1) .GT. 30.) .AND. (julian_diff .GE. 121. .AND. julian_diff .LE. 243.).AND. &
1531                  (veget_max(ji,jv) .GT. min_sechiba)) .OR. &
1532                  ((lalo(ji,1) .LT. -30.) .AND. (julian_diff .GE. 305. .OR. julian_diff .LE. 59.) .AND.&
1533                  (veget_max(ji,jv) .GT. min_sechiba))) THEN
1534                ! 1e12 for conversion from kg to Ng
1535                ! 1/(365/12*24*60*60*4) for conversion from year to seconds corrected for 4 months of emissions
1536                flx_fertil_no(ji,jv) = (N_qt_WRICE_pft(ji,jv)*(1/30.)+N_qt_OTHER_pft(ji,jv))*(2.5/100.) &
1537                     & *1e12/(365*24*60*60*4/12)/(area(ji)*veget_max(ji,jv))
1538                ! OR Tropical regions all the year
1539             ELSE IF ((lalo(ji,1) .GE. -30.).AND.(lalo(ji,1) .LE. 30.).AND.(veget_max(ji,jv) .GT. min_sechiba)) THEN
1540                flx_fertil_no(ji,jv) = (N_qt_WRICE_pft(ji,jv)*(1/30.)+N_qt_OTHER_pft(ji,jv))*(2.5/100.) &
1541                     & *1e12/(365*24*60*60)/(area(ji)*veget_max(ji,jv))
1542             END IF
1543             flx_no_soil(ji,jv) = flx_no_soil(ji,jv) + flx_fertil_no(ji,jv)
1544          END IF
1545
1546          !! 7.5 Calculation of net NO flux above soil accounting for surface deposition,
1547          !! based on the Canopy Reduction Factor (CRF), calculated using Leaf Area and Stomatal Area
1548          !kc=cuticle absorptivity = 0.24m^2/m^2
1549          !ks=stomatal absorptivity = 8.75m^2/m^2
1550          !Larch=Larcher SAI/LAI ratio given in Larcher 1991
1551          !> @codeinc
1552          CRF(ji,jv) = (exp(-8.75*Larch(jv)*lai(ji,jv)) + exp(-0.24*lai(ji,jv)))/2.
1553          flx_no(ji,jv) = flx_no_soil(ji,jv)*CRF(ji,jv)
1554          !> @endcodeinc
1555       END DO
1556    END DO
1557
1558
1559    ! Write output with XIOS
1560    CALL xios_orchidee_send_field("PAR",PAR)
1561    CALL xios_orchidee_send_field("flx_fertil_no",flx_fertil_no)
1562    CALL xios_orchidee_send_field("flx_iso",flx_iso)
1563    CALL xios_orchidee_send_field("flx_mono",flx_mono)
1564    CALL xios_orchidee_send_field("flx_ORVOC",flx_ORVOC)
1565    CALL xios_orchidee_send_field("flx_MBO",flx_MBO)
1566    CALL xios_orchidee_send_field("flx_methanol",flx_methanol)
1567    CALL xios_orchidee_send_field("flx_acetone",flx_acetone)
1568    CALL xios_orchidee_send_field("flx_acetal",flx_acetal)
1569    CALL xios_orchidee_send_field("flx_formal",flx_formal)
1570    CALL xios_orchidee_send_field("flx_acetic",flx_acetic)
1571    CALL xios_orchidee_send_field("flx_formic",flx_formic)
1572    CALL xios_orchidee_send_field("flx_no_soil",flx_no_soil)
1573    CALL xios_orchidee_send_field("flx_no",flx_no)
1574    CALL xios_orchidee_send_field('flx_apinen'   , flx_apinen)
1575    CALL xios_orchidee_send_field('flx_bpinen'   , flx_bpinen)
1576    CALL xios_orchidee_send_field('flx_limonen'  ,flx_limonen)
1577    CALL xios_orchidee_send_field('flx_myrcen'   , flx_myrcen)
1578    CALL xios_orchidee_send_field('flx_sabinen'  ,flx_sabinen)
1579    CALL xios_orchidee_send_field('flx_camphen'  ,flx_camphen)
1580    CALL xios_orchidee_send_field('flx_3caren'   , flx_3caren)
1581    CALL xios_orchidee_send_field('flx_tbocimen' ,flx_tbocimen)
1582    CALL xios_orchidee_send_field('flx_othermono',flx_othermono)
1583    CALL xios_orchidee_send_field('flx_sesquiter',flx_sesquiter)
1584    CALL xios_orchidee_send_field("CRF",CRF)
1585    CALL xios_orchidee_send_field('fco2', fco2)
1586
1587    IF ( ok_radcanopy ) THEN
1588       CALL xios_orchidee_send_field("PARdf",PARdf)
1589       CALL xios_orchidee_send_field("PARdr",PARdr)
1590       
1591       IF (ok_multilayer) THEN
1592          CALL xios_orchidee_send_field("PARsuntab",PARsuntab)
1593          CALL xios_orchidee_send_field("PARshtab",PARshtab)
1594       ELSE
1595          CALL xios_orchidee_send_field("PARsun",PARsun)
1596          CALL xios_orchidee_send_field("PARsh",PARsh)
1597          CALL xios_orchidee_send_field("laisun",laisun)
1598          CALL xios_orchidee_send_field("laish",laish)
1599       ENDIF
1600    ENDIF
1601
1602    IF ( ok_bbgfertil_Nox ) THEN
1603       CALL xios_orchidee_send_field("flx_co2_bbg_year",flx_co2_bbg_year)
1604    END IF
1605
1606    IF ( ok_cropsfertil_Nox ) THEN
1607       CALL xios_orchidee_send_field("N_qt_WRICE_year",N_qt_WRICE_year)
1608       CALL xios_orchidee_send_field("N_qt_OTHER_year",N_qt_OTHER_year)
1609    END IF
1610   
1611
1612    ! Write output with IOIPSL
1613    IF ( .NOT. almaoutput ) THEN
1614
1615       CALL histwrite_p(hist_id, 'PAR', kjit, PAR, kjpindex, index)
1616       IF ( ok_radcanopy ) THEN
1617          CALL histwrite_p(hist_id, 'laisun', kjit, laisun, kjpindex*nvm, indexveg)
1618          CALL histwrite_p(hist_id, 'laish', kjit, laish, kjpindex*nvm, indexveg)
1619          CALL histwrite_p(hist_id, 'Fdf', kjit, Fdf, kjpindex, index)
1620          IF (ok_multilayer) THEN
1621             CALL histwrite_p(hist_id, 'PARsuntab', kjit, PARsuntab, kjpindex*(nlai+1), indexlai)
1622             CALL histwrite_p(hist_id, 'PARshtab', kjit, PARshtab, kjpindex*(nlai+1), indexlai)
1623          ELSE
1624             CALL histwrite_p(hist_id, 'PARsun', kjit, PARsun, kjpindex*nvm, indexveg)
1625             CALL histwrite_p(hist_id, 'PARsh', kjit, PARsh, kjpindex*nvm, indexveg)
1626          END IF
1627          CALL histwrite_p(hist_id, 'coszang', kjit, coszang, kjpindex, index)
1628          CALL histwrite_p(hist_id, 'PARdf', kjit, PARdf, kjpindex, index)
1629          CALL histwrite_p(hist_id, 'PARdr', kjit, PARdr, kjpindex, index)
1630          CALL histwrite_p(hist_id, 'Trans', kjit, Trans, kjpindex, index)
1631       END IF
1632       CALL histwrite_p(hist_id, 'flx_fertil_no', kjit, flx_fertil_no, kjpindex*nvm, indexveg)
1633       CALL histwrite_p(hist_id, 'CRF', kjit, CRF, kjpindex*nvm, indexveg)
1634       CALL histwrite_p(hist_id, 'fco2', kjit, fco2, kjpindex*nvm, indexveg)
1635
1636       IF ( ok_bbgfertil_Nox ) THEN
1637          CALL histwrite_p(hist_id, 'flx_co2_bbg_year', 1, flx_co2_bbg_year, kjpindex, index)
1638       ENDIF
1639       IF ( ok_cropsfertil_Nox ) THEN
1640          CALL histwrite_p(hist_id, 'N_qt_WRICE_year', 1, N_qt_WRICE_year, kjpindex, index)
1641          CALL histwrite_p(hist_id, 'N_qt_OTHER_year', 1, N_qt_OTHER_year, kjpindex, index)
1642       ENDIF
1643       CALL histwrite_p(hist_id, 'flx_iso', kjit, flx_iso, kjpindex*nvm, indexveg)
1644       CALL histwrite_p(hist_id, 'flx_mono', kjit, flx_mono, kjpindex*nvm, indexveg)
1645       CALL histwrite_p(hist_id, 'flx_apinen', kjit, flx_apinen, kjpindex*nvm, indexveg)
1646       CALL histwrite_p(hist_id, 'flx_bpinen', kjit, flx_bpinen, kjpindex*nvm, indexveg)
1647       CALL histwrite_p(hist_id, 'flx_limonen', kjit, flx_limonen, kjpindex*nvm, indexveg)
1648       CALL histwrite_p(hist_id, 'flx_myrcen', kjit, flx_myrcen, kjpindex*nvm, indexveg)
1649       CALL histwrite_p(hist_id, 'flx_sabinen', kjit, flx_sabinen, kjpindex*nvm, indexveg)
1650       CALL histwrite_p(hist_id, 'flx_camphen', kjit, flx_camphen, kjpindex*nvm, indexveg)
1651       CALL histwrite_p(hist_id, 'flx_3caren', kjit, flx_3caren, kjpindex*nvm, indexveg)
1652       CALL histwrite_p(hist_id, 'flx_tbocimen', kjit, flx_tbocimen, kjpindex*nvm, indexveg)
1653       CALL histwrite_p(hist_id, 'flx_othermono', kjit, flx_othermono, kjpindex*nvm, indexveg)
1654       CALL histwrite_p(hist_id, 'flx_sesquiter', kjit, flx_sesquiter, kjpindex*nvm, indexveg)
1655       CALL histwrite_p(hist_id, 'flx_ORVOC', kjit, flx_ORVOC, kjpindex*nvm, indexveg)
1656       CALL histwrite_p(hist_id, 'flx_MBO', kjit, flx_MBO, kjpindex*nvm, indexveg)
1657       CALL histwrite_p(hist_id, 'flx_methanol', kjit, flx_methanol, kjpindex*nvm, indexveg)
1658       CALL histwrite_p(hist_id, 'flx_acetone', kjit, flx_acetone, kjpindex*nvm, indexveg)
1659       CALL histwrite_p(hist_id, 'flx_acetal', kjit, flx_acetal, kjpindex*nvm, indexveg)
1660       CALL histwrite_p(hist_id, 'flx_formal', kjit, flx_formal, kjpindex*nvm, indexveg)
1661       CALL histwrite_p(hist_id, 'flx_acetic', kjit, flx_acetic, kjpindex*nvm, indexveg)
1662       CALL histwrite_p(hist_id, 'flx_formic', kjit, flx_formic, kjpindex*nvm, indexveg)
1663       CALL histwrite_p(hist_id, 'flx_no_soil', kjit, flx_no_soil, kjpindex*nvm, indexveg)
1664       CALL histwrite_p(hist_id, 'flx_no', kjit, flx_no, kjpindex*nvm, indexveg)
1665       
1666       IF ( hist2_id > 0 ) THEN
1667          CALL histwrite_p(hist2_id, 'PAR', kjit, PAR, kjpindex, index)
1668          IF ( ok_radcanopy ) THEN
1669             CALL histwrite_p(hist2_id, 'PARsun', kjit, PARsun, kjpindex*nvm, indexveg)
1670             CALL histwrite_p(hist2_id, 'PARsh', kjit, PARsh, kjpindex*nvm, indexveg)
1671             CALL histwrite_p(hist2_id, 'laisun', kjit, laisun, kjpindex*nvm, indexveg)
1672             CALL histwrite_p(hist2_id, 'laish', kjit, laish, kjpindex*nvm, indexveg)
1673          ENDIF
1674          CALL histwrite_p(hist2_id, 'flx_fertil_no', kjit, flx_fertil_no, kjpindex*nvm, indexveg)
1675          CALL histwrite_p(hist2_id, 'CRF', kjit, CRF, kjpindex*nvm, indexveg)
1676          IF ( ok_bbgfertil_Nox ) THEN
1677             CALL histwrite_p(hist2_id, 'flx_co2_bbg_year', 1, flx_co2_bbg_year, kjpindex, index)
1678          ENDIF
1679          IF ( ok_cropsfertil_Nox ) THEN
1680             CALL histwrite_p(hist2_id, 'N_qt_WRICE_year', 1, N_qt_WRICE_year, kjpindex, index)
1681             CALL histwrite_p(hist2_id, 'N_qt_OTHER_year', 1, N_qt_OTHER_year, kjpindex, index)
1682          ENDIF
1683          CALL histwrite_p(hist2_id, 'flx_iso', kjit, flx_iso, kjpindex*nvm, indexveg)
1684          CALL histwrite_p(hist2_id, 'flx_mono', kjit, flx_mono, kjpindex*nvm, indexveg)
1685          CALL histwrite_p(hist2_id, 'flx_apinen', kjit, flx_apinen, kjpindex*nvm, indexveg)
1686          CALL histwrite_p(hist2_id, 'flx_bpinen', kjit, flx_bpinen, kjpindex*nvm, indexveg)
1687          CALL histwrite_p(hist2_id, 'flx_limonen', kjit, flx_limonen, kjpindex*nvm, indexveg)
1688          CALL histwrite_p(hist2_id, 'flx_myrcen', kjit, flx_myrcen, kjpindex*nvm, indexveg)
1689          CALL histwrite_p(hist2_id, 'flx_sabinen', kjit, flx_sabinen, kjpindex*nvm, indexveg)
1690          CALL histwrite_p(hist2_id, 'flx_camphen', kjit, flx_camphen, kjpindex*nvm, indexveg)
1691          CALL histwrite_p(hist2_id, 'flx_3caren', kjit, flx_3caren, kjpindex*nvm, indexveg)
1692          CALL histwrite_p(hist2_id, 'flx_tbocimen', kjit, flx_tbocimen, kjpindex*nvm, indexveg)
1693          CALL histwrite_p(hist2_id, 'flx_othermono', kjit, flx_othermono, kjpindex*nvm, indexveg)
1694          CALL histwrite_p(hist2_id, 'flx_sesquiter', kjit, flx_sesquiter, kjpindex*nvm, indexveg)
1695          CALL histwrite_p(hist2_id, 'flx_ORVOC', kjit, flx_ORVOC, kjpindex*nvm, indexveg)
1696          CALL histwrite_p(hist2_id, 'flx_MBO', kjit, flx_MBO, kjpindex*nvm, indexveg)
1697          CALL histwrite_p(hist2_id, 'flx_methanol', kjit, flx_methanol, kjpindex*nvm, indexveg)
1698          CALL histwrite_p(hist2_id, 'flx_acetone', kjit, flx_acetone, kjpindex*nvm, indexveg)
1699          CALL histwrite_p(hist2_id, 'flx_acetal', kjit, flx_acetal, kjpindex*nvm, indexveg)
1700          CALL histwrite_p(hist2_id, 'flx_formal', kjit, flx_formal, kjpindex*nvm, indexveg)
1701          CALL histwrite_p(hist2_id, 'flx_acetic', kjit, flx_acetic, kjpindex*nvm, indexveg)
1702          CALL histwrite_p(hist2_id, 'flx_formic', kjit, flx_formic, kjpindex*nvm, indexveg)
1703          CALL histwrite_p(hist2_id, 'flx_no_soil', kjit, flx_no_soil, kjpindex*nvm, indexveg)
1704          CALL histwrite_p(hist2_id, 'flx_no', kjit, flx_no, kjpindex*nvm, indexveg)
1705       ENDIF
1706    ENDIF
1707
1708    IF (printlev>=3) WRITE(numout,*) 'OK chemistry_bvoc'
1709
1710  END SUBROUTINE chemistry_bvoc
1711
1712!! ================================================================================================================================
1713!! SUBROUTINE   : chemistry_flux_interface
1714!!
1715!>\BRIEF         This subroutine will give to the interface between inca model and orchidee model in sechiba all the flux ask by inca
1716!!
1717!! DESCRIPTION  :  This subroutine allow the transfer of fluxes between surface and atmospheric chemistry. It is called from sechiba
1718!!
1719!! RECENT CHANGE(S): None
1720!!
1721!! MAIN OUTPUT VARIABLE(S): ::
1722!!
1723!! REFERENCE(S) : None
1724!!
1725!! FLOWCHART    : None
1726!_ ================================================================================================================================
1727
1728 SUBROUTINE chemistry_flux_interface( field_out_names, fields_out, field_in_names, fields_in)
1729
1730     !
1731     ! Optional arguments
1732     !
1733     ! Names and fields for emission variables : to be transport by Orchidee to Inca
1734     CHARACTER(LEN=*),DIMENSION(:), OPTIONAL, INTENT(IN) :: field_out_names
1735     REAL(r_std),DIMENSION(:,:,:), OPTIONAL, INTENT(OUT) :: fields_out
1736     !
1737     ! Names and fields for deposit variables : to be transport from chemistry model by INCA to ORCHIDEE.
1738     CHARACTER(LEN=*),DIMENSION(:), OPTIONAL, INTENT(IN) :: field_in_names
1739     REAL(r_std),DIMENSION(:,:), OPTIONAL, INTENT(IN)    :: fields_in
1740     !
1741     ! Number of fields to give (nb_fields_out) or get from (nb_fields_in) INCA :
1742     INTEGER(i_std), SAVE   :: nb_fields_out, nb_fields_in
1743!$OMP THREADPRIVATE(nb_fields_out, nb_fields_in)
1744
1745     ! Id of fields to give (nb_fields_out) or get from (nb_fields_in) INCA :
1746     INTEGER(i_std)  :: i_fields_out, i_fields_in
1747
1748     IF (l_first_chemistry_inca) THEN
1749
1750        ! il faut verifier que ok_bvoc (chemistry_ok_bvoc) est bien active sinon on arrete tout
1751        if (.not.ok_bvoc) then
1752          CALL ipslerr_p (3,'chemistry_inca', &
1753            &          'you need to activate chemistry_ok_bvoc in orchidee.def',&
1754            &          'This model won''t be able to continue.', &
1755            &          '')
1756        endif
1757
1758        ! Prepare fieds out/in for interface with INCA.
1759        IF (PRESENT(field_out_names)) THEN
1760           nb_fields_out=SIZE(field_out_names)
1761        ELSE
1762           nb_fields_out=0
1763        ENDIF
1764
1765        IF (PRESENT(field_in_names)) THEN
1766           nb_fields_in=SIZE(field_in_names)
1767        ELSE
1768           nb_fields_in=0
1769        ENDIF
1770        l_first_chemistry_inca = .FALSE.
1771
1772     ENDIF
1773
1774
1775
1776     IF (ok_bvoc) THEN 
1777        ! Fields for deposit variables : to be transport from chemistry model by INCA to ORCHIDEE.
1778        DO i_fields_in=1,nb_fields_in
1779           SELECT CASE(TRIM(field_in_names(i_fields_in)))
1780           CASE DEFAULT 
1781              CALL ipslerr_p (3,'chemistry_inca', &
1782                   &          'You ask in INCA an unknown field '//TRIM(field_in_names(i_fields_in))//&
1783                   &          ' to give to ORCHIDEE for this specific version.',&
1784                   &          'This model won''t be able to continue.', &
1785                   &          '(check your tracer parameters in INCA)')
1786           END SELECT
1787        ENDDO
1788       
1789        ! Fields for Biogenic emissions : to be transport by Orchidee to Inca
1790        DO i_fields_out=1,nb_fields_out
1791           SELECT CASE(TRIM(field_out_names(i_fields_out)))
1792           CASE("flx_iso") 
1793              fields_out(:,:,i_fields_out)=flx_iso(:,:)
1794           CASE("flx_mono") 
1795              fields_out(:,:,i_fields_out)=flx_mono(:,:)
1796           CASE("flx_ORVOC") 
1797              fields_out(:,:,i_fields_out)=flx_ORVOC(:,:)
1798           CASE("flx_MBO") 
1799              fields_out(:,:,i_fields_out)=flx_MBO(:,:)
1800           CASE("flx_methanol") 
1801              fields_out(:,:,i_fields_out)=flx_methanol(:,:)
1802           CASE("flx_acetone") 
1803              fields_out(:,:,i_fields_out)=flx_acetone(:,:)
1804           CASE("flx_acetal") 
1805              fields_out(:,:,i_fields_out)=flx_acetal(:,:)
1806           CASE("flx_formal") 
1807              fields_out(:,:,i_fields_out)=flx_formal(:,:)
1808           CASE("flx_acetic") 
1809              fields_out(:,:,i_fields_out)=flx_acetic(:,:)
1810           CASE("flx_formic") 
1811              fields_out(:,:,i_fields_out)=flx_formic(:,:)
1812           CASE("flx_no_soil") 
1813              fields_out(:,:,i_fields_out)=flx_no_soil(:,:)
1814           CASE("flx_nox") 
1815              fields_out(:,:,i_fields_out)=flx_no(:,:)
1816           CASE("flx_fertil_no") 
1817              fields_out(:,:,i_fields_out)=flx_fertil_no(:,:)
1818           CASE("flx_apinen")
1819              fields_out(:,:,i_fields_out)=flx_apinen(:,:)
1820           CASE("flx_bpinen")
1821              fields_out(:,:,i_fields_out)=flx_bpinen(:,:)
1822           CASE("flx_limonen")
1823              fields_out(:,:,i_fields_out)=flx_limonen(:,:)
1824           CASE("flx_myrcen")
1825              fields_out(:,:,i_fields_out)=flx_myrcen(:,:)
1826           CASE("flx_sabinen")
1827              fields_out(:,:,i_fields_out)=flx_sabinen(:,:)
1828           CASE("flx_camphen")
1829              fields_out(:,:,i_fields_out)=flx_camphen(:,:)
1830           CASE("flx_3caren")
1831              fields_out(:,:,i_fields_out)=flx_3caren(:,:)
1832           CASE("flx_tbocimen")
1833              fields_out(:,:,i_fields_out)=flx_tbocimen(:,:)
1834           CASE("flx_othermono")
1835              fields_out(:,:,i_fields_out)=flx_othermono(:,:)
1836           CASE("flx_sesquiter")
1837              fields_out(:,:,i_fields_out)=flx_sesquiter(:,:)
1838             
1839           CASE DEFAULT 
1840              CALL ipslerr_p (3,'chemistry_inca', &
1841                   &          'You ask from INCA an unknown field '//TRIM(field_out_names(i_fields_out))//&
1842                   &          ' to ORCHIDEE for this specific version.',&
1843                   &          'This model won''t be able to continue.', &
1844                   &          '(check your tracer parameters in INCA)')
1845           END SELECT
1846        ENDDO
1847       
1848     ENDIF
1849   END SUBROUTINE chemistry_flux_interface
1850
1851!! ================================================================================================================================
1852!! SUBROUTINE   : chemistry_clear
1853!!
1854!>\BRIEF         Clear chemistry module
1855!!
1856!! DESCRIPTION  :  Deallocate memory and reset initialization variables to there original values
1857!!
1858!_ ================================================================================================================================
1859   SUBROUTINE chemistry_clear
1860
1861     !! 1. Initialize as for first run
1862     l_first_chemistry_inca = .TRUE.
1863     
1864     !! 2. Deallocate dynamic variables
1865     IF (ALLOCATED(pulse)) DEALLOCATE (pulse)
1866     IF (ALLOCATED (ok_siesta)) DEALLOCATE (ok_siesta) 
1867     IF (ALLOCATED (allow_pulse)) DEALLOCATE (allow_pulse)
1868     IF (ALLOCATED (pulseday)) DEALLOCATE (pulseday)
1869     IF (ALLOCATED (siestaday)) DEALLOCATE (siestaday)
1870     IF (ALLOCATED (pulselim)) DEALLOCATE (pulselim)
1871     IF (ALLOCATED (siestalim)) DEALLOCATE (siestalim)
1872     IF (ALLOCATED (N_qt_WRICE_year)) DEALLOCATE (N_qt_WRICE_year)
1873     IF (ALLOCATED (N_qt_OTHER_year)) DEALLOCATE (N_qt_OTHER_year)
1874     IF (ALLOCATED (flx_iso)) DEALLOCATE (flx_iso)
1875     IF (ALLOCATED (flx_mono)) DEALLOCATE (flx_mono)
1876     IF (ALLOCATED (flx_ORVOC)) DEALLOCATE (flx_ORVOC)
1877     IF (ALLOCATED (flx_MBO)) DEALLOCATE (flx_MBO)
1878     IF (ALLOCATED (flx_methanol)) DEALLOCATE (flx_methanol)
1879     IF (ALLOCATED (flx_acetone)) DEALLOCATE (flx_acetone)
1880     IF (ALLOCATED (flx_acetal)) DEALLOCATE (flx_acetal)
1881     IF (ALLOCATED (flx_formal)) DEALLOCATE (flx_formal)
1882     IF (ALLOCATED (flx_acetic)) DEALLOCATE (flx_acetic)
1883     IF (ALLOCATED (flx_formic)) DEALLOCATE (flx_formic)
1884     IF (ALLOCATED (flx_no_soil)) DEALLOCATE (flx_no_soil)
1885     IF (ALLOCATED (flx_no)) DEALLOCATE (flx_no)
1886     IF (ALLOCATED (flx_fertil_no)) DEALLOCATE (flx_fertil_no)
1887     IF (ALLOCATED (flx_apinen)) DEALLOCATE (flx_apinen)
1888     IF (ALLOCATED (flx_bpinen)) DEALLOCATE (flx_bpinen)
1889     IF (ALLOCATED (flx_limonen)) DEALLOCATE (flx_limonen)
1890     IF (ALLOCATED (flx_myrcen)) DEALLOCATE (flx_myrcen)
1891     IF (ALLOCATED (flx_sabinen)) DEALLOCATE (flx_sabinen)
1892     IF (ALLOCATED (flx_camphen)) DEALLOCATE (flx_camphen)
1893     IF (ALLOCATED (flx_3caren)) DEALLOCATE (flx_3caren)
1894     IF (ALLOCATED (flx_tbocimen)) DEALLOCATE (flx_tbocimen)
1895     IF (ALLOCATED (flx_othermono)) DEALLOCATE (flx_othermono)
1896     IF (ALLOCATED (flx_sesquiter)) DEALLOCATE (flx_sesquiter)
1897     IF (ALLOCATED (CRF)) DEALLOCATE (CRF)
1898     IF (ALLOCATED (flx_co2_bbg_year)) DEALLOCATE (flx_co2_bbg_year)
1899
1900   END SUBROUTINE chemistry_clear
1901   
1902END MODULE chemistry
Note: See TracBrowser for help on using the repository browser.