SURFEX  V8_0
Surfex V8_0 release
 All Classes Files Functions Variables
fapair.F90
Go to the documentation of this file.
1 !SFX_LIC Copyright 1994-2014 CNRS, Meteo-France and Universite Paul Sabatier
2 !SFX_LIC This is part of the SURFEX software governed by the CeCILL-C licence
3 !SFX_LIC version 1. See LICENSE, CeCILL-C_V1-en.txt and CeCILL-C_V1-fr.txt
4 !SFX_LIC for details. version 1.
5 ! ######
6 SUBROUTINE fapair(PABC, PFD_SKY, PIA, PLAI, PXMUS, PSSA_SUP, PSSA_INF, &
7  pb_sup, pb_inf, palb_veg, palb_soil, oshade, &
8  pfapr, pfapr_bs, plai_eff, piacan, &
9  piacan_shade, piacan_sunlit, pfrac_sun )
10 ! #########################################################################
11 !
12 !!**** *FAPAIR*
13 !!
14 !! PURPOSE
15 !! -------
16 !! Calculates fraction of absorbed photosynthetic active radiation (FAPAR) and
17 !! fraction of absorbed near infrared (FAPIR) of vegetation and bare soil.
18 !!
19 !!** METHOD
20 !! ------
21 !! Carrer et al.
22 !!
23 !! EXTERNAL
24 !! --------
25 !! none
26 !!
27 !! IMPLICIT ARGUMENTS
28 !! ------------------
29 !!
30 !! USE MODD_SURF_PAR
31 !! USE MODD_CSTS
32 !! USE MODD_CO2V_PAR
33 !! USE MODI_CCETR_PAIR
34 !!
35 !! REFERENCE
36 !! ---------
37 !! Carrer et al, 2013 (doi:10.1002/jgrg20070)
38 !!
39 !! AUTHOR
40 !! ------
41 !! D. Carrer * Meteo-France *
42 !!
43 !! MODIFICATIONS
44 !! -------------
45 !! Original 01/04/2011
46 !! Commented by C. Delire 07/13
47 !! C. Delire 08/13 : moved calculation of diffuse fraction from here to radiative_transfert.F90
48 !!
49 !!-------------------------------------------------------------------------------
50 USE modd_surf_par, ONLY : xundef
51 USE modd_co2v_par, ONLY : xk_sup, xk_inf, xxsi_sup, xxsi_inf ! clumping index parameters (Carrer et al 2.1.3)
52 !
53 USE modi_ccetr_pair
54 !
55 USE yomhook ,ONLY : lhook, dr_hook
56 USE parkind1 ,ONLY : jprb
57 !
58 !* 0. DECLARATIONS
59 ! ------------
60 !
61 IMPLICIT NONE
62 !
63 !* 0.1 declarations of arguments
64 !
65 REAL, DIMENSION(:), INTENT(IN) :: pabc ! abscissa needed for integration
66 ! ! of net assimilation and stomatal
67 ! ! conductance over canopy depth
68 REAL, DIMENSION(:), INTENT(IN) :: pfd_sky ! fraction of diffused radiation in sky
69 REAL, DIMENSION(:), INTENT(IN) :: pia ! incident PAR or NIR
70 REAL, DIMENSION(:), INTENT(IN) :: plai ! leaf area index
71 REAL, DIMENSION(:), INTENT(IN) :: pxmus ! cosine of solar zenith angle
72 REAL, INTENT(IN) :: pssa_sup, pssa_inf !single scatering albedo (PAR or NIR)
73 REAL, DIMENSION(:), INTENT(IN) :: pb_sup, pb_inf !sigma parameter in clumping (Table 1, eq4)
74 REAL, DIMENSION(:), INTENT(IN) :: palb_veg, palb_soil
75 LOGICAL, DIMENSION(:), INTENT(IN) :: oshade ! OSHADE = if 1 shading activated
76 !
77 REAL, DIMENSION(:), INTENT(OUT) :: pfapr
78 REAL, DIMENSION(:), INTENT(OUT) :: pfapr_bs
79 REAL, DIMENSION(:), OPTIONAL, INTENT(OUT) :: plai_eff
80 !
81 REAL, DIMENSION(:,:), OPTIONAL, INTENT(OUT) :: piacan ! PAR in the canopy at different gauss level
82 REAL, DIMENSION(:,:), OPTIONAL, INTENT(OUT) :: piacan_shade ! PAR in the canopy at different gauss level
83 REAL, DIMENSION(:,:), OPTIONAL, INTENT(OUT) :: piacan_sunlit ! PAR in the canopy at different gauss level
84 REAL, DIMENSION(:,:), OPTIONAL, INTENT(OUT) :: pfrac_sun ! fraction of sunlit leaves
85 !
86 !* 0.2 declarations of local variables
87 !
88 !
89 REAL, DIMENSION(SIZE(PLAI)) :: zxia, zxia_sup, zkmusp_sup, zkmusp_inf
90 REAL, DIMENSION(SIZE(PLAI)) :: zb_dr_sup, zb_dr_inf, zomega_dr_sup, zomega_dr_inf, &
91  zomega_df_sup, zomega_df_inf
92 ! ZXIA = abs. radiation of vegetation
93 REAL, DIMENSION(SIZE(PLAI)) :: ztr, zfd_veg, zfd_sup, zlai_eff0, zlai_eff
94 ! ZTR = transmittance
95 ! ZFD_VEG, ZFD_SUP = fraction of radiation diffused by the considered medium (vegetation)
96 !REAL, DIMENSION(SIZE(PLAI)) :: ZXIA_SUNLIT, ZXIA_SHADE, ZLAI_SUNLIT, ZLAI_SHADE
97 ! ZXIA_SUNLIT = absorbed PAR of sunlit leaves
98 ! ZXIA_SHADE = absorbed PAR of shaded leaves
99 ! ZLAI_SUNLIT = LAI of sunlit leaves
100 ! ZLAI_SHADE = LAI of shaded leaves
101 !REAL, DIMENSION(SIZE(PLAI)) :: ZRN_SUNLIT, ZRN_SHADE
102 REAL, DIMENSION(SIZE(PLAI),SIZE(PABC)) :: ziacan, ziacan_sunlit, ziacan_shade, zfrac_sun
103 REAL :: zabc, zweight, zcoef, zsup, zinf, &
104  zssa_sup, zssa_inf, zb_df_sup, zb_df_inf
105 ! ZABC = abscissa needed for integration
106 ! of net assimilation and stomatal
107 ! conductance over canopy depth (working scalar)
108 ! ZSUP, ZINF = d_sup and d_inf from Table 1. Carrer et al.
109 !
110 INTEGER :: jint, i ! index for loops
111 REAL(KIND=JPRB) :: zhook_handle
112 !-------------------------------------------------------------------------------
113 IF (lhook) CALL dr_hook('FAPAIR',0,zhook_handle)
114 !
115 ! initialisation
116 !
117 ztr(:) = 1.0
118 !
119 zxia_sup(:) = 0.
120 !
121 zfd_veg(:) = 0.
122 zfd_sup(:) = 0.
123 !
124 !ZXIA_SUNLIT(:) = 0.
125 !ZXIA_SHADE(:) = 0.
126 !ZLAI_SUNLIT(:) = 0.
127 !ZLAI_SHADE(:) = 0.
128 !
129 zlai_eff(:) = 0.
130 !
131 ziacan(:,:) = 0.
132 ziacan_sunlit(:,:) = 0.
133 ziacan_shade(:,:) = 0.
134 zfrac_sun(:,:) = 0.
135 !
136 pfapr(:) = 0.
137 pfapr_bs(:) = 0.
138 !PRN_SHADE(:) = 0.
139 !PRN_SUNLIT(:) = 0.
140 !
141 !
142 IF (pabc(SIZE(pabc)).GT.0.8) zfd_veg(:) = min(pfd_sky(:),1.)
143 ! set param sup / inf
144 !
145 zssa_sup = sqrt(1.-pssa_sup)
146 zssa_inf = sqrt(1.-pssa_inf)
147 !
148 zsup = - 0.461 * xxsi_sup + 3.8 ! d_sup Carrer et al Table 1
149 zinf = - 0.461 * xxsi_inf + 3.8 ! d_inf Carrer et al Table 1
150 !
151 DO i=1,SIZE(pia)
152  IF (pia(i).NE.0.) THEN
153  zkmusp_sup(i) = exp(-xk_sup*(acos(pxmus(i)))**zsup)
154  zkmusp_inf(i) = exp(-xk_inf*(acos(pxmus(i)))**zinf)
155  ! direct case
156  ! Directional albedo of upper/lower layer
157  zb_dr_sup(i) = 1.-(1.-zssa_sup)/(1.+2.*pxmus(i)*zssa_sup) ! Carrer et al. b_dr Table 1
158  zb_dr_inf(i) = 1.-(1.-zssa_sup)/(1.+2.*pxmus(i)*zssa_inf)
159  ! CLUMPING INDEX
160  zomega_dr_sup(i) = 1. / (1.+ pb_sup(i)*zkmusp_sup(i)) ! Carrer et al. eq. 4a
161  zomega_dr_inf(i) = 1. / (1.+ pb_inf(i)*zkmusp_inf(i))
162  ! diffuse case
163  ! CLUMPING INDEX
164  zomega_df_sup(i) = (1.+pb_sup(i)/2.)/(1.+pb_sup(i)) ! Carrer et al. eq. 4b
165  zomega_df_inf(i) = (1.+pb_inf(i)/2.)/(1.+pb_inf(i))
166  ENDIF
167 ENDDO
168 !
169 ! Non-directional albedo of diffuse
170 zb_df_sup = 1.-(1.-zssa_sup)/(1.+ zssa_sup)
171 zb_df_inf = 1.-(1.-zssa_inf)/(1.+ zssa_inf)
172 !
173 ! Integration over the canopy: SIZE(PABC) increments
174 ! are used to approximate the integral. And to calculate
175 ! absorded fluxes within the canopy and in the bare soil
176 DO jint = SIZE(pabc),1,-1
177 !
178  zabc = 1. ! normalized height unit of the layer above
179  IF (jint.LT.SIZE(pabc)) zabc = pabc(jint+1)
180  zweight = zabc - pabc(jint)
181  !
182  IF (pabc(jint).GT.0.8) THEN
183  ! Compute transmittance of each level
184  CALL ccetr_pair(jint, pabc(jint), zabc, pia, pxmus, zb_dr_sup, &
185  zomega_dr_sup, zomega_df_sup, zb_df_sup, plai, &
186  palb_veg, palb_soil, pfd_sky, zfd_veg, ztr, &
187  zxia, zlai_eff0 )
188  ELSE
189  CALL ccetr_pair(jint, pabc(jint), zabc, pia, pxmus, zb_dr_inf, &
190  zomega_dr_inf, zomega_df_inf, zb_df_inf, plai, &
191  palb_veg, palb_soil, pfd_sky, zfd_veg, ztr, &
192  zxia, zlai_eff0 )
193  ENDIF
194  !
195  DO i=1,SIZE(pia)
196  !
197  zxia(i) = max(0.,zxia(i))
198  ziacan(i,jint) = max(0.,zxia(i)-zxia_sup(i))
199  zxia_sup(i) = zxia(i)
200  !
201  zlai_eff0(i) = max(0.,zlai_eff0(i))
202  zlai_eff(i) = zlai_eff(i) + zlai_eff0(i)
203  !
204  !calculate a FAPAR/FAPIR of the entire canopy
205  pfapr(i) = pfapr(i) + ziacan(i,jint)
206  !
207  !------------------------------------------------------
208  ! If LSHADE=0 no shading, only sunlit leaves
209  ! If LSHADE=1 shading
210  ! PIACAN is used to calculate An of each level within the canopy in cotwores
211  ! ZIACAN_SUNLIT used for net assimilation of a sunlit leave in COTWO
212  ! ZIACAN_SHADE used in A-gs for net assimilation of a shaded leave in COTWO
213  IF (oshade(i)) THEN
214  !
215  !sunlit leaves
216  !absorbed PAR of an equivalent canopy representative of the layer of leaves eq. (8)
217  zcoef = (1.0-zfd_sup(i))/ztr(i)+ zfd_sup(i)
218  ziacan_sunlit(i,jint) = zcoef/(zweight*max(0.0001,plai(i)))*ziacan(i,jint)
219  !not sunlit leaves
220  ziacan_shade(i,jint) = max(0.,zfd_sup(i)/(zweight*max(0.0001,plai(i)))*ziacan(i,jint))
221  !
222  !ZXIA_SUNLIT(I) = ZXIA_SUNLIT(I) + ZWEIGHT*ZTR(I) *ZIACAN_SUNLIT(I,JINT)
223  !ZLAI_SUNLIT(I) = ZLAI_SUNLIT(I) + ZWEIGHT*ZTR(I)*ZCOEF*PLAI(I)
224  !
225  !ZXIA_SHADE(I) = ZXIA_SHADE(I) + ZWEIGHT*(1-ZTR(I)) *ZIACAN_SHADE(I,JINT)
226  !ZLAI_SHADE(I) = ZLAI_SHADE(I) + ZWEIGHT*(1-ZTR(I))*ZFD_SUP(I)*PLAI(I)
227  !
228  zfrac_sun(i,jint) = ztr(i) !fraction of sunlit leaves
229  !
230  ELSE
231  !
232  ziacan_sunlit(i,jint) = max(0.,ziacan(i,jint)/(zweight*max(0.0001,plai(i))))
233  !ZLAI_SUNLIT(I) = ZLAI_SUNLIT(I) + ZWEIGHT*PLAI(I)
234  !
235  ENDIF
236  !
237  zfd_sup(i) = zfd_veg(i)
238  !
239  ENDDO
240  !
241 END DO
242 !
243 !
244 WHERE (pia(:).NE.0.)
245  pfapr(:) = pfapr(:) / pia(:)
246  pfapr_bs(:)=(1.-palb_veg(:))*(1-palb_soil(:))*(1.+palb_veg(:)*palb_soil(:))*ztr(:)
247  WHERE (plai(:).EQ.0) pfapr_bs(:) = 1-palb_soil(:)
248 END WHERE
249 !
250 !WHERE (ZLAI_SHADE(:) .NE.0.) ZRN_SHADE(:) = ZXIA_SHADE(:) / ZLAI_SHADE(:)
251 !WHERE (ZLAI_SUNLIT(:).NE.0.) ZRN_SUNLIT(:) = ZXIA_SUNLIT(:)/ ZLAI_SUNLIT(:)
252 !
253 IF (present(plai_eff)) plai_eff = zlai_eff
254 IF (present(piacan)) piacan = ziacan
255 IF (present(piacan_sunlit)) piacan_sunlit = ziacan_sunlit
256 IF (present(piacan_shade)) piacan_shade = ziacan_shade
257 IF (present(pfrac_sun)) pfrac_sun = zfrac_sun
258 !
259 IF (lhook) CALL dr_hook('FAPAIR',1,zhook_handle)
260 !
261 END SUBROUTINE fapair
subroutine ccetr_pair(KNIV, PABC, PABC_SUP, PIA, PXMUS, PB_DR, POMEGA_DR, POMEGA_DF, PB_DF, PLAI, PALB_VEG, PALB_SOIL, PFD_SKY, PFD_VEG, PTR, PXIA, PLAI_EFF)
Definition: ccetr_pair.F90:7
subroutine fapair(PABC, PFD_SKY, PIA, PLAI, PXMUS, PSSA_SUP, PSSA_INF, PB_SUP, PB_INF, PALB_VEG, PALB_SOIL, OSHADE, PFAPR, PFAPR_BS, PLAI_EFF, PIACAN, PIACAN_SHADE, PIACAN_SUNLIT, PFRAC_SUN)
Definition: fapair.F90:6