source: tags/ORCHIDEE_1_9_5/ORCHIDEE/src_global/solar.f90 @ 8

Last change on this file since 8 was 8, checked in by orchidee, 14 years ago

import first tag equivalent to CVS orchidee_1_9_5 + OOL_1_9_5

File size: 6.0 KB
Line 
1!! This module define solar
2!!
3!! @call sechiba_main
4!! @Version : $Revision: 1.1 $, $Date: 2010/04/06 14:26:07 $
5!!
6!! $Header: /home/ssipsl/CVSREP/ORCHIDEE/src_global/solar.f90,v 1.1 2010/04/06 14:26:07 ssipsl Exp $
7!!
8!! @author Marie-Alice Foujols, Jan Polcher and Martial Mancip
9!!
10!!
11MODULE solar
12
13  USE defprec
14  USE constantes
15  USE parallel
16
17  IMPLICIT NONE
18
19
20
21CONTAINS
22  !
23  !f90doc CONTAINS
24  !
25  !
26  !
27  !-
28  !=====================================================================
29  SUBROUTINE solarang (julian, julian0, iim, jjm, lon, lat, csang)
30    !---------------------------------------------------------------------
31    !- This subroutine computes the solar angle according to the method
32    !- used by GSWP and developed by J.C. Morill.
33    !- See for further details :
34    !- http://www.atmo.arizona.edu/~morrill/swrad.html
35    !---------------------------------------------------------------------
36    !
37    USE calendar
38    !
39    IMPLICIT NONE
40    !-
41    INTEGER,INTENT(in) :: iim, jjm
42    REAL,INTENT(in)  :: julian
43    REAL,INTENT(in)  :: julian0
44    REAL,DIMENSION(iim,jjm), INTENT(in)  :: lon, lat
45    REAL,DIMENSION(iim,jjm), INTENT(out) :: csang
46    !-
47    REAL :: pi, gamma, dec, decd
48    REAL :: et, gmt, le, ls, lcorr, latime, omega, omegad
49    REAL :: llatd, llat
50    INTEGER :: igmt, ilon, ilat
51    INTEGER,SAVE,ALLOCATABLE :: zone(:), lhour(:)
52    !
53    LOGICAL :: check = .FALSE.
54    !---------------------------------------------------------------------
55    !
56    pi = 4.*ATAN(1.)
57    IF (check) WRITE(numout,*) 'We get the right calendar information'
58    !-
59    !- 1) Day angle gamma
60    !-
61    !   gamma = 2.*pi*MOD(julian,one_year)/one_year
62    gamma = 2.*pi*(julian-julian0)/one_year
63    !-
64    !- 2) Solar declination (assumed constant for a 24 hour period)  page 7
65    !-    in radians
66    !-
67    IF (check) WRITE(numout,*) 'Solar declination'
68    !
69    dec = ( 0.006918-0.399912*COS(gamma)+0.070257*SIN(gamma) &
70         &       -0.006758*COS(2*gamma)+0.000907*SIN(2*gamma)      &
71         &       -0.002697*COS(3*gamma)+0.00148*SIN(3*gamma))
72    decd = dec*(180/pi)
73    !-
74    !- maximum error 0.0006 rad (<3'), leads to error
75    !- of less than 1/2 degree in zenith angle
76    !-
77    IF (check) WRITE(numout,*) 'Equation of time'
78    !- 3)  Equation of time  page 11
79    !-
80    et = ( 0.000075+0.001868*COS(gamma)-0.032077*SIN(gamma)&
81         &      -0.014615*COS(2*gamma)-0.04089*SIN(2*gamma))*229.18
82    !-
83    !- Get the time zones for the current time
84    !-
85    gmt = 24.*(julian-INT(julian))
86    IF (.NOT.ALLOCATED(zone))  ALLOCATE(zone(iim))
87    IF (.NOT.ALLOCATED(lhour)) ALLOCATE(lhour(iim))
88    !-
89    igmt = NINT(gmt)
90    IF (check) WRITE(numout,*) 'Get time zone'
91    CALL time_zone(igmt, iim, jjm, lon, zone, lhour)
92    !-
93    !- Start a loop through the grid
94    !-
95    IF (check) WRITE(numout,*) 'Start a loop through the grid'
96    DO ilon=1,iim
97       !---
98       !--- 4) Local apparent time  page 13
99       !---
100       !--- ls     standard longitude (nearest 15 degree meridian)
101       !--- le     local longtitude
102       !--- lhour  local standard time
103       !--- latime local apparent time
104       !--- lcorr  longitudunal correction (minutes)
105       !---
106       le = lon(ilon,1)
107       ls = ((zone(ilon)-1)*15)-180.
108       lcorr = 4.*(ls-le)*(-1)
109       latime = lhour(ilon)+lcorr/60.+et/60.
110       IF (latime <  0.) latime = latime+24
111       IF (latime > 24.) latime = latime-24
112       !---
113       !--- 5) Hour angle omega  page 15
114       !---
115       !--- hour angle is zero at noon, positive in the morning
116       !--- It ranges from 180 to -180
117       !---
118       !--- omegad is omega in degrees, omega is in radians
119       !---
120       omegad = (latime-12.)*(-15.)
121       omega  = omegad*pi/180.
122       !---
123       DO ilat=1,jjm
124          !-----
125          !----- 6)  Zenith angle  page 15
126          !-----
127          !----- csang cosine of zenith angle (radians)
128          !----- llatd =  local latitude in degrees
129          !----- llat  =  local latitude in radians
130          !-----
131          llatd = lat(1,ilat)
132          llat  = llatd*pi/180.
133          csang(ilon,ilat) = &
134               &  MAX(0.,SIN(dec)*SIN(llat)+COS(dec)*COS(llat)*COS(omega))
135       ENDDO
136    ENDDO
137    !----------------------
138  END SUBROUTINE solarang
139  !-
140  !=====================================================================
141  SUBROUTINE time_zone (gmt, iim, jjm, lon, zone, lhour)
142    !---------------------------------------------------------------------
143    !
144    IMPLICIT NONE
145    !
146    INTEGER :: gmt, iim, jjm, zone(iim), lhour(iim)
147    REAL :: lon(iim,jjm)
148    !-
149    INTEGER :: deg
150    !!??   REAL :: deg
151    INTEGER :: i, ilon, change
152    !---------------------------------------------------------------------
153    DO ilon=1,iim
154       !---
155       !--- Convert longitude index to longtitude in degrees
156       !---
157       deg = lon(ilon,1)
158       !---
159       !--- Determine into which time zone (15 degree interval) the
160       !--- longitude falls.
161       !---
162       DO i=1,25
163          IF (deg < (-187.5+(15*i))) THEN
164             zone(ilon) = i
165             IF (zone(ilon) == 25)   zone(ilon) = 1
166             EXIT
167          ENDIF
168       ENDDO
169       !---
170       !--- Calculate change (in number of hours) from GMT time to
171       !--- local hour.  Change will be negative for zones < 13 and
172       !--- positive for zones > 13.
173       !---
174       !--- There is also a correction for lhour < 0 and lhour > 23
175       !--- to lhour between 0 and 23.
176       !---
177       IF (zone(ilon) < 13) THEN
178          change = zone(ilon)-13
179          lhour(ilon) = gmt+change
180       ENDIF
181       !---
182       IF (zone(ilon) == 13) THEN
183          lhour(ilon) = gmt
184       ENDIF
185       !---
186       IF (zone(ilon) > 13) THEN
187          change = zone(ilon)-13
188          lhour(ilon) = gmt+change
189       ENDIF
190       IF (lhour(ilon) <  0) lhour(ilon) = lhour(ilon)+24
191       IF (lhour(ilon) > 23) lhour(ilon) = lhour(ilon)-24
192       !---
193    ENDDO
194    !-----------------------
195  END SUBROUTINE time_zone
196
197!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
198 
199END MODULE solar
Note: See TracBrowser for help on using the repository browser.