source: branches/publications/ORCHIDEE_CN_CAN_r5698/src_global/solar.f90 @ 7346

Last change on this file since 7346 was 4292, checked in by sebastiaan.luyssaert, 7 years ago

DEV: tested at 4 pixels for 2 years for all PFTs. Updated the code to r3851.

  • Property svn:keywords set to Revision Date HeadURL Date Author Revision
File size: 20.8 KB
Line 
1
2! =================================================================================================================================
3! MODULE       : solar
4!
5! CONTACT      : orchidee-help _at_ ipsl.jussieu.fr
6!
7! LICENCE      : IPSL (2011)
8! This software is governed by the CeCILL licence see ORCHIDEE/ORCHIDEE_CeCILL.LIC
9!
10!>\BRIEF        This module define solar
11!!
12!!\n DESCRIPTION:
13!!
14!! RECENT CHANGE(S): None
15!!
16!! REFERENCE(S) : None
17!!
18!! SVN          :
19!! $HeadURL$
20!! $Date$
21!! $Revision$
22!! \n
23!_ ================================================================================================================================
24
25MODULE solar
26
27  USE defprec
28  USE constantes
29  USE ioipsl_para 
30
31  IMPLICIT NONE
32
33
34
35CONTAINS
36
37
38!! ================================================================================================================================
39!! SUBROUTINE   : solarang
40!!
41!>\BRIEF         This subroutine computes the solar angle according to the method used by GSWP and developed by J.C. Morill.
42!!
43!! DESCRIPTION  : None
44!!
45!! RECENT CHANGE(S): None
46!!
47!! MAIN OUTPUT VARIABLE(S): ::csang
48!!
49!! REFERENCE(S) : See for further details :
50!! http://www.atmo.arizona.edu/~morrill/swrad.html
51!!
52!! FLOWCHART    : None
53!! \n
54!_ ================================================================================================================================
55
56  SUBROUTINE solarang (julian, julian0, iim, jjm, lon, lat, csang)
57
58    USE calendar
59
60    IMPLICIT NONE
61
62    !! 0. Parameters and variables declaration
63   
64    !! 0.1 Input variables
65
66    INTEGER, INTENT(in)                      :: iim, jjm        !!
67    REAL, INTENT(in)                         :: julian          !!
68    REAL, INTENT(in)                         :: julian0         !!
69    REAL, DIMENSION(iim,jjm), INTENT(in)     :: lon, lat        !!
70
71    !! 0.2 Output variables
72
73    REAL, DIMENSION(iim,jjm), INTENT(out)    :: csang           !!
74
75    !! 0.4 Local variables
76
77    REAL                                     :: gamma           !!
78    REAL                                     :: dec             !!
79    REAL                                     :: decd            !!
80    REAL                                     :: et              !!
81    REAL                                     :: gmt             !!
82    REAL                                     :: le              !!
83    REAL                                     :: ls              !!
84    REAL                                     :: lcorr           !!
85    REAL                                     :: latime          !!
86    REAL                                     :: omega           !!
87    REAL                                     :: omegad          !!
88    REAL                                     :: llatd, llat     !!
89    INTEGER                                  :: igmt            !!
90    INTEGER                                  :: ilon, ilat      !!
91    INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: zone            !!
92!$OMP THREADPRIVATE(zone)
93    REAL, SAVE, ALLOCATABLE, DIMENSION(:)    :: lhour           !!
94!$OMP THREADPRIVATE(lhour)
95    !
96    LOGICAL                                  :: check = .FALSE. !!
97!$OMP THREADPRIVATE(check)
98
99!_ ================================================================================================================================
100
101    IF (check) WRITE(numout,*) 'Calendar information', julian, julian0
102    !-
103    !- 1) Day angle gamma
104    !-
105    !   gamma = 2.*pi*MOD(julian,one_year)/one_year
106    gamma = 2.*pi*(julian-julian0)/one_year
107    !-
108    !- 2) Solar declination (assumed constant for a 24 hour period)  page 7
109    !-    in radians
110    !-
111    IF (check) WRITE(numout,*) 'Solar declination', gamma, one_year
112    !
113    dec = ( 0.006918-0.399912*COS(gamma)+0.070257*SIN(gamma) &
114         &       -0.006758*COS(2*gamma)+0.000907*SIN(2*gamma)      &
115         &       -0.002697*COS(3*gamma)+0.00148*SIN(3*gamma))
116    decd = dec*(180/pi)
117    !-
118    !- maximum error 0.0006 rad (<3'), leads to error
119    !- of less than 1/2 degree in zenith angle
120    !-
121    IF (check) WRITE(numout,*) 'Equation of time', decd
122    !- 3)  Equation of time  page 11
123    !-
124    et = ( 0.000075+0.001868*COS(gamma)-0.032077*SIN(gamma)&
125         &      -0.014615*COS(2*gamma)-0.04089*SIN(2*gamma))*229.18
126    !-
127    !- Get the time zones for the current time
128    !-
129    gmt = 24.*(julian-INT(julian))
130    IF (.NOT.ALLOCATED(zone))  ALLOCATE(zone(iim))
131    IF (.NOT.ALLOCATED(lhour)) ALLOCATE(lhour(iim))
132    !-
133    !igmt = NINT(gmt)
134    IF (check) WRITE(numout,*) 'Get time zone', et, gmt
135    CALL time_zone(gmt, iim, jjm, lon, zone, lhour)
136    !-
137    !- Start a loop through the grid
138    !-
139    IF (check) WRITE(numout,*) 'Start a loop through the grid'
140    DO ilon=1,iim
141       !---
142       !--- 4) Local apparent time  page 13
143       !---
144       !--- ls     standard longitude (nearest 15 degree meridian)
145       !--- le     local longtitude
146       !--- lhour  local standard time
147       !--- latime local apparent time
148       !--- lcorr  longitudunal correction (minutes)
149       !---
150       le = lon(ilon,1)
151       ls = ((zone(ilon)-1)*15)-180.
152       lcorr = 4.*(ls-le)*(-1)
153       latime = lhour(ilon)+lcorr/60.+et/60.
154       IF (latime <  0.) latime = latime+24
155       IF (latime > 24.) latime = latime-24
156       !---
157       !--- 5) Hour angle omega  page 15
158       !---
159       !--- hour angle is zero at noon, positive in the morning
160       !--- It ranges from 180 to -180
161       !---
162       !--- omegad is omega in degrees, omega is in radians
163       !---
164       omegad = (latime-12.)*(-15.)
165       omega  = omegad*pi/180.
166       !---
167       DO ilat=1,jjm
168          !-----
169          !----- 6)  Zenith angle  page 15
170          !-----
171          !----- csang cosine of zenith angle (radians)
172          !----- llatd =  local latitude in degrees
173          !----- llat  =  local latitude in radians
174          !-----
175          llatd = lat(1,ilat)
176          llat  = llatd*pi/180.
177          csang(ilon,ilat) = &
178               &  MAX(zero,SIN(dec)*SIN(llat)+COS(dec)*COS(llat)*COS(omega))
179       ENDDO
180    ENDDO
181    !----------------------
182  END SUBROUTINE solarang
183
184
185!! ================================================================================================================================
186!! SUBROUTINE   : time_zone
187!!
188!>\BRIEF         
189!!
190!! DESCRIPTION  : None
191!!
192!! RECENT CHANGE(S): None
193!!
194!! MAIN OUTPUT VARIABLE(S): ::zone, ::lhour
195!!
196!! REFERENCE(S) : None
197!!
198!! FLOWCHART    : None
199!! \n
200!_ ================================================================================================================================
201
202  SUBROUTINE time_zone (gmt, iim, jjm, lon, zone, lhour)
203
204    IMPLICIT NONE
205
206    !! 0. Parameters and variables declaration
207   
208    !! 0.1 Input variables
209
210    INTEGER, INTENT(in)                   :: iim, jjm        !!
211    REAL, DIMENSION(iim,jjm), INTENT(in)  :: lon             !!
212    REAL, INTENT(in)                      :: gmt             !!
213
214    !! 0.2 Output variables
215
216    INTEGER, DIMENSION(iim), INTENT(out)  :: zone            !!
217    REAL, DIMENSION(iim), INTENT(out)     :: lhour           !!
218
219    !! 0.4 Local variables
220
221    INTEGER                               :: deg             !!
222    !!??   REAL :: deg
223    INTEGER                               :: i, ilon, change !! Indices
224
225!_ ================================================================================================================================
226
227    DO ilon=1,iim
228       !---
229       !--- Convert longitude index to longtitude in degrees
230       !---
231       deg = lon(ilon,1)
232       !---
233       !--- Determine into which time zone (15 degree interval) the
234       !--- longitude falls.
235       !---
236       DO i=1,25
237          IF (deg < (-187.5+(15*i))) THEN
238             zone(ilon) = i
239             IF (zone(ilon) == 25)   zone(ilon) = 1
240             EXIT
241          ENDIF
242       ENDDO
243       !---
244       !--- Calculate change (in number of hours) from GMT time to
245       !--- local hour.  Change will be negative for zones < 13 and
246       !--- positive for zones > 13.
247       !---
248       !--- There is also a correction for lhour < 0 and lhour > 23
249       !--- to lhour between 0 and 23.
250       !---
251       IF (zone(ilon) < 13) THEN
252          change = zone(ilon)-13
253          lhour(ilon) = gmt+change
254       ENDIF
255       !---
256       IF (zone(ilon) == 13) THEN
257          lhour(ilon) = gmt
258       ENDIF
259       !---
260       IF (zone(ilon) > 13) THEN
261          change = zone(ilon)-13
262          lhour(ilon) = gmt+change
263       ENDIF
264       IF (lhour(ilon) <  0) lhour(ilon) = lhour(ilon)+24
265       IF (lhour(ilon) >= 24) lhour(ilon) = lhour(ilon)-24
266       !---
267    ENDDO
268    !-----------------------
269  END SUBROUTINE time_zone
270
271
272
273!! ================================================================================================================================
274!! SUBROUTINE   : downward_solar_flux
275!!
276!>\BRIEF         
277!!
278!! DESCRIPTION  : None
279!!
280!! RECENT CHANGE(S): None
281!!
282!! MAIN OUTPUT VARIABLE(S): ::solad, ::solai
283!!
284!! REFERENCE(S) :See for further details :
285!! http://www.atmo.arizona.edu/~morrill/swrad.html
286!!
287!! FLOWCHART    : None
288!! \n
289!_ ================================================================================================================================
290
291  SUBROUTINE downward_solar_flux (npoi,latitude,calendar_str,jday,rtime,cloud,nband,solad,solai)
292
293    IMPLICIT NONE
294
295    !! 0. Parameter and variables declaration
296   
297    !! 0.1 Input variables
298   
299    INTEGER, INTENT(in)                      :: npoi               !! number of grid points (unitless)
300    REAL, DIMENSION(npoi), INTENT(in)        :: latitude           !! latitude in degrees
301    CHARACTER(LEN=20),INTENT(in)             :: calendar_str       !! Calendar type
302    REAL, INTENT(in)                         :: jday
303    REAL,INTENT(in)                          :: rtime
304    REAL, DIMENSION(npoi), INTENT(in)        :: cloud              !! cloud fraction (0-1, unitless)
305    INTEGER,INTENT(in)                       :: nband              !! number of visible bands (unitless)
306
307    !! 0.2 Output variables
308
309    REAL, DIMENSION(npoi,nband), INTENT(out) :: solad              !! direct downward solar flux
310                                                                   !! @tex $(W.m^{-2})$ @endtex
311    REAL, DIMENSION(npoi,nband), INTENT(out) :: solai              !! diffuse downward solar flux
312                                                                   !! @tex $(W.m^{-2})$ @endtex
313
314    !! 0.4 Local variables
315
316    ! Parametres orbitaux:
317
318    REAL, SAVE                               :: ecc                !! Eccentricity
319!$OMP THREADPRIVATE(ecc)
320    REAL, SAVE                               :: perh               !! Longitude of perihelie
321!$OMP THREADPRIVATE(perh)
322    REAL, SAVE                               :: xob                !! obliquity
323!$OMP THREADPRIVATE(xob)
324    REAL, PARAMETER                          :: pir = pi/180.      !!
325    REAL, SAVE                               :: step               !!
326!$OMP THREADPRIVATE(step)
327    REAL                                     :: xl                 !!
328    REAL                                     :: so                 !!
329    REAL                                     :: xllp               !!
330    REAL                                     :: xee                !!
331    REAL                                     :: xse                !! 
332    REAL                                     :: xlam               !!
333    REAL                                     :: dlamm              !!
334    REAL                                     :: anm                !!
335    REAL                                     :: ranm               !!
336    REAL                                     :: ranv               !!
337    REAL                                     :: anv                !!
338    REAL                                     :: tls                !!
339    REAL                                     :: rlam               !!
340    REAL                                     :: sd                 !!
341    REAL                                     :: cd                 !!
342    REAL                                     :: deltar             !!
343    REAL                                     :: delta              !!
344    REAL                                     :: Dis_ST             !! 
345    REAL                                     :: ddt                !!
346    REAL                                     :: orbit              !!
347    REAL                                     :: angle              !!
348    REAL                                     :: xdecl              !!
349    REAL                                     :: xlat               !!
350    REAL, DIMENSION(npoi)                    :: trans              !!
351    REAL, DIMENSION(npoi)                    :: fdiffuse           !!
352    REAL, DIMENSION(npoi)                    :: coszen             !! cosine of solar zenith angle
353    REAL                                     :: sw                 !!
354    REAL                                     :: frac               !!
355    INTEGER                                  :: i,ib               !!
356    LOGICAL, SAVE                            :: lstep_init = .TRUE. !!
357!$OMP THREADPRIVATE(lstep_init)
358
359!_ ================================================================================================================================
360   
361    IF (lstep_init) THEN
362       IF ( TRIM(calendar_str) .EQ. 'gregorian' ) THEN 
363          step = 1.0
364       ELSE
365          step = one_year/365.2425
366       ENDIF
367       !-
368       ! Read Orbital Parameters
369       !-
370
371       !Config Key   = ECCENTRICITY
372       !Config Desc  = Use prescribed values
373       !Config If    = ALLOW_WEATHERGEN
374       !Config Def   = 0.016724
375       !Config Help  =
376       !Config Units = [-]
377       ecc = 0.016724
378       CALL getin_p ('ECCENTRICITY',ecc)
379       WRITE(numout,*) 'ECCENTRICITY: ',ecc
380       !
381       !Config Key  = PERIHELIE
382       !Config Desc = Use prescribed values
383       !Config If   = ALLOW_WEATHERGEN
384       !Config Def  = 102.04
385       !Config Help  =
386       !Config Units = [-]
387       perh = 102.04
388       CALL getin_p ('PERIHELIE',perh)
389       WRITE(numout,*) 'PERIHELIE: ',perh
390       !
391       !Config Key  = OBLIQUITY
392       !Config Desc = Use prescribed values
393       !Config If   = ALLOW_WEATHERGEN
394       !Config Def  = 23.446
395       !Config Help  =
396       !Config Units = [Degrees]
397       xob = 23.446
398       CALL getin_p ('OBLIQUITY',xob)
399       WRITE(numout,*) 'OBLIQUITY: ',xob
400
401       lstep_init = .FALSE.
402    ENDIF
403
404    !-
405    ! calendar and orbital calculations
406    !-
407
408    !-
409    ! calculate the earth's orbital angle (around the sun) in radians
410    orbit = 2.0*pi*jday/365.2425
411    !-
412    ! calculate the solar hour angle in radians
413    angle = 2.0*pi*(rtime-12.0)/24.0
414    !-
415    ! calculate the current solar declination angle
416    ! ref: global physical climatology, hartmann, appendix a
417    !
418    ! xdecl = 0.006918 &
419    !        -0.399912*cos(orbit) &
420    !        +0.070257*sin(orbit) &
421    !        -0.006758*cos(2.0*orbit) &
422    !        +0.000907*sin(2.0*orbit) &
423    !        -0.002697*cos(3.0*orbit) &
424    !        +0.001480*sin(3.0*orbit)
425    !
426    ! calculate the effective solar constant,
427    ! including effects of eccentricity
428    ! ref: global physical climatology, hartmann, appendix a
429    !
430    ! sw = 1370.*( 1.000110 &
431    !             +0.034221*cos(orbit) &
432    !             +0.001280*sin(orbit) &
433    !             +0.000719*cos(2.0*orbit) &
434    !             +0.000077*sin(2.0*orbit))
435    !
436    ! correction Marie-France Loutre
437    !
438    !    orbital parameters
439    !
440    !    ecc = 0.016724
441    !    perh = 102.04
442    !    xob = 23.446
443    !-
444    xl = perh+180.0
445    ! so : sinus of obliquity
446    so = sin(xob*pir)
447    !-
448    xllp = xl*pir
449    xee  = ecc*ecc
450    xse  = sqrt(1.0d0-xee)
451    ! xlam : true long. sun for mean long. = 0
452    xlam = (ecc/2.0+ecc*xee/8.0d0)*(1.0+xse)*sin(xllp)-xee/4.0 &
453         &      *(0.5+xse)*sin(2.0*xllp)+ecc*xee/8.0*(1.0/3.0+xse) &
454         &      *sin(3.0*xllp)
455    xlam  =2.0d0*xlam/pir
456    ! dlamm : mean long. sun for ma-ja
457    dlamm =xlam+(INT(jday)-79)*step
458    anm   = dlamm-xl
459    ranm  = anm*pir
460    xee   = xee*ecc
461    ! ranv : anomalie vraie    (radian)
462    ranv = ranm+(2.0*ecc-xee/4.0)*sin(ranm)+5.0/4.0*ecc*ecc &
463         &      *sin(2.0*ranm)+13.0/12.0*xee*sin(3.0*ranm)
464    ! anv  : anomalie vraie    (degrees)
465    anv  = ranv/pir
466    ! tls  : longitude vraie   (degrees)
467    tls  = anv+xl
468    ! rlam : longitude vraie   (radian)
469    rlam = tls*pir
470    ! sd and cd: cosinus and sinus of solar declination angle (delta)
471    ! sinus delta = sin (obl)*sin(lambda) with lambda = real longitude
472    ! (Phd. thesis of Marie-France Loutre, ASTR-UCL, Belgium, 1993)
473    sd    = so*sin(rlam)
474    cd    = sqrt(1.0d0-sd*sd)
475    ! delta : Solar Declination (degrees, angle sun at equator)
476    deltar = atan(sd/cd)
477    delta  = deltar/pir
478    !-
479    ! Eccentricity Effect
480    !-
481    Dis_ST  =(1.0-ecc*ecc)/(1.0+ecc*cos(ranv))
482    ! ddt :    1 / normalized earth's sun distance
483    ddt = 1.0/Dis_ST
484    !-
485    ! Insolation normal to the atmosphere (W/m2)
486    !-
487    sw  = ddt *ddt *1370.d0
488    !-
489    xdecl = deltar
490    !-
491    ! solar calculations
492    !-
493    do i=1,npoi
494       !---
495       !-- calculate the latitude in radians
496       !---
497       xlat = latitude(i)*pi/180.0
498       !---
499       !-- calculate the cosine of the solar zenith angle
500       !---
501       coszen(i) = MAX(0.0, (sin(xlat)*sin(xdecl) &
502            &                      + cos(xlat)*cos(xdecl)*cos(angle)))
503       !---
504       !-- calculate the solar transmission through the atmosphere
505       !-- using simple linear function of tranmission and cloud cover
506       !---
507       !-- note that the 'cloud cover' data is typically obtained from
508       !-- sunshine hours -- not direct cloud observations
509       !---
510       !-- where, cloud cover = 1 - sunshine fraction
511       !---
512       !-- different authors present different values for the slope and
513       !-- intercept terms of this equation
514       !---
515       !-- Friend, A: Parameterization of a global daily weather generator
516       !-- for terrestrial ecosystem and biogeochemical modelling,
517       !-- Ecological Modelling
518       !---
519       !-- Spitters et al., 1986: Separating the diffuse and direct component
520       !-- of global radiation and its implications for modeling canopy
521       !-- photosynthesis, Part I: Components of incoming radiation,
522       !-- Agricultural and Forest Meteorology, 38, 217-229.
523       !---
524       !-- A. Friend       : trans = 0.251+0.509*(1.0-cloud(i))
525       !-- Spitters et al. : trans = 0.200+0.560*(1.0-cloud(i))
526       !---
527       !-- we are using the values from A. Friend
528       !---
529       trans(i) = 0.251+0.509*(1.0-cloud(i))
530       !---
531       !-- calculate the fraction of indirect (diffuse) solar radiation
532       !-- based upon the cloud cover
533       !---
534       !-- note that these relationships typically are measured for either
535       !-- monthly or daily timescales, and may not be exactly appropriate
536       !-- for hourly calculations -- however, in ibis, cloud cover is fixed
537       !-- through the entire day so it may not make much difference
538       !---
539       !-- method i --
540       !---
541       !-- we use a simple empirical relationships
542       !-- from Nikolov and Zeller (1992)
543       !---
544       !-- Nikolov, N. and K.F. Zeller, 1992:  A solar radiation algorithm
545       !-- for ecosystem dynamics models, Ecological Modelling, 61, 149-168.
546       !---
547       fdiffuse(i) = 1.0045+trans(i)*( 0.0435+trans(i) &
548            &                                 *(-3.5227+trans(i)*2.6313))
549       !---
550       IF (trans(i) > 0.75) fdiffuse(i) = 0.166
551       !---
552       !-- method ii --
553       !---
554       !-- another method was suggested by Spitters et al. (1986) based on
555       !-- long-term data from the Netherlands
556       !--
557       !-- Spitters et al., 1986: Separating the diffuse and direct component
558       !-- of global radiation and its implications for modeling canopy
559       !-- photosynthesis, Part I: Components of incoming radiation,
560       !-- Agricultural and Forest Meteorology, 38, 217-229.
561       !---
562       !--       if ((trans == 0.00).and.(trans < 0.07)) then
563       !--         fdiffuse = 1.0
564       !--       else if ((trans >= 0.07).and.(trans < 0.35)) then
565       !--         fdiffuse = 1.0-2.3*(trans-0.07)**2
566       !--       else if ((trans >= 0.35).and.(trans < 0.75)) then
567       !--         fdiffuse = 1.33-1.46*trans
568       !--       ELSE
569       !--         fdiffuse = 0.23
570       !--       endif
571       !---
572    ENDDO
573    !-----------------------
574    !-
575    ! do for each waveband
576    !-
577    DO ib=1,nband
578       !---
579       !-- calculate the fraction in each waveband
580       !---
581       !-- GK010200
582       IF (nband == 2) then
583          frac = 0.46+0.08*REAL(ib-1)
584       ELSE
585          frac = 1./REAL(nband)
586       ENDIF
587       !---
588       DO i=1,npoi
589          !-----
590          !---- calculate the direct and indirect solar radiation
591          !-----
592          solad(i,ib) = sw*coszen(i)*frac*trans(i)*(1.-fdiffuse(i))
593          solai(i,ib) = sw*coszen(i)*frac*trans(i)*fdiffuse(i)
594       ENDDO
595    ENDDO
596
597
598  END SUBROUTINE downward_solar_flux
599
600 
601END MODULE solar
Note: See TracBrowser for help on using the repository browser.