[2] | 1 | SUBROUTINE alboc(rjour,rlat,albedo) |
---|
| 2 | IMPLICIT none |
---|
| 3 | c====================================================================== |
---|
| 4 | c Auteur(s): Z.X. Li (LMD/CNRS) (adaptation du GCM du LMD) |
---|
| 5 | c Date: le 16 mars 1995 |
---|
| 6 | c Objet: Calculer l'albedo sur l'ocean |
---|
| 7 | c Methode: Integrer numeriquement l'albedo pendant une journee |
---|
| 8 | c |
---|
| 9 | c Arguments; |
---|
| 10 | c rjour (in,R) : jour dans l'annee (a compter du 1 janvier) |
---|
| 11 | c rlat (in,R) : latitude en degre |
---|
| 12 | c albedo (out,R): albedo obtenu (de 0 a 1) |
---|
| 13 | c====================================================================== |
---|
| 14 | #include "dimensions.h" |
---|
| 15 | #include "dimphy.h" |
---|
| 16 | #include "YOMCST.h" |
---|
| 17 | c |
---|
| 18 | REAL fmagic ! un facteur magique pour regler l'albedo |
---|
| 19 | PARAMETER (fmagic=0.7) |
---|
| 20 | INTEGER npts ! il controle la precision de l'integration |
---|
| 21 | PARAMETER (npts=120) ! 120 correspond a l'interval 6 minutes |
---|
| 22 | c |
---|
| 23 | REAL rlat(klon), rjour, albedo(klon) |
---|
| 24 | REAL zdist, zlonsun, zpi, zdeclin |
---|
| 25 | REAL rmu,alb, srmu, salb, fauxo, aa, bb |
---|
| 26 | INTEGER i, k |
---|
| 27 | c |
---|
| 28 | zpi = 4. * ATAN(1.) |
---|
| 29 | c |
---|
| 30 | c Calculer la longitude vraie de l'orbite terrestre: |
---|
| 31 | CALL orbite(rjour,zlonsun,zdist) |
---|
| 32 | c |
---|
| 33 | c Calculer la declinaison du soleil (qui varie entre + et - R_incl): |
---|
| 34 | zdeclin = ASIN(SIN(zlonsun*zpi/180.0)*SIN(R_incl*zpi/180.0)) |
---|
| 35 | c |
---|
| 36 | DO 999 i=1,klon |
---|
| 37 | aa = SIN(rlat(i)*zpi/180.0) * SIN(zdeclin) |
---|
| 38 | bb = COS(rlat(i)*zpi/180.0) * COS(zdeclin) |
---|
| 39 | c |
---|
| 40 | c Midi local (angle du temps = 0.0): |
---|
| 41 | rmu = aa + bb * COS(0.0) |
---|
| 42 | rmu = MAX(0.0, rmu) |
---|
| 43 | fauxo = (1.47-ACOS(rmu))/.15 |
---|
| 44 | alb = 0.03+0.630/(1.+fauxo*fauxo) |
---|
| 45 | srmu = rmu |
---|
| 46 | salb = alb * rmu |
---|
| 47 | c |
---|
| 48 | c Faire l'integration numerique de midi a minuit (le facteur 2 |
---|
| 49 | c prend en compte l'autre moitie de la journee): |
---|
| 50 | DO k = 1, npts |
---|
| 51 | rmu = aa + bb * COS(FLOAT(k)/FLOAT(npts)*zpi) |
---|
| 52 | rmu = MAX(0.0, rmu) |
---|
| 53 | fauxo = (1.47-ACOS(rmu))/.15 |
---|
| 54 | alb = 0.03+0.630/(1.+fauxo*fauxo) |
---|
| 55 | srmu = srmu + rmu * 2.0 |
---|
| 56 | salb = salb + alb*rmu * 2.0 |
---|
| 57 | ENDDO |
---|
| 58 | IF (srmu .NE. 0.0) THEN |
---|
| 59 | albedo(i) = salb / srmu * fmagic |
---|
| 60 | ELSE ! nuit polaire (on peut prendre une valeur quelconque) |
---|
| 61 | albedo(i) = fmagic |
---|
| 62 | ENDIF |
---|
| 63 | 999 CONTINUE |
---|
| 64 | RETURN |
---|
| 65 | END |
---|
| 66 | c===================================================================== |
---|
| 67 | SUBROUTINE alboc_cd(rmu0,albedo) |
---|
| 68 | IMPLICIT none |
---|
| 69 | c====================================================================== |
---|
| 70 | c Auteur(s): Z.X. Li (LMD/CNRS) |
---|
| 71 | c date: 19940624 |
---|
| 72 | c Calculer l'albedo sur l'ocean en fonction de l'angle zenithal moyen |
---|
| 73 | c Formule due a Larson and Barkstrom (1977) Proc. of the symposium |
---|
| 74 | C on radiation in the atmosphere, 19-28 August 1976, science Press, |
---|
| 75 | C 1977 pp 451-453, ou These de 3eme cycle de Sylvie Joussaume. |
---|
| 76 | c |
---|
| 77 | c Arguments |
---|
| 78 | c rmu0 (in): cosinus de l'angle solaire zenithal |
---|
| 79 | c albedo (out): albedo de surface de l'ocean |
---|
| 80 | c====================================================================== |
---|
| 81 | #include "dimensions.h" |
---|
| 82 | #include "dimphy.h" |
---|
| 83 | REAL rmu0(klon), albedo(klon) |
---|
| 84 | c |
---|
| 85 | REAL fauxo |
---|
| 86 | INTEGER i |
---|
| 87 | c |
---|
| 88 | DO i = 1, klon |
---|
| 89 | fauxo = ( 1.47 - ACOS( rmu0(i) ) )/0.15 |
---|
| 90 | albedo(i) = 1.1*( .03 + .630/( 1. + fauxo*fauxo)) |
---|
| 91 | albedo(i) = MAX(MIN(albedo(i),0.60),0.04) |
---|
| 92 | ENDDO |
---|
| 93 | c |
---|
| 94 | RETURN |
---|
| 95 | END |
---|
| 96 | c======================================================================== |
---|
| 97 | SUBROUTINE albsno(veget, agesno, alb_neig) |
---|
| 98 | IMPLICIT none |
---|
| 99 | c |
---|
| 100 | #include "dimensions.h" |
---|
| 101 | #include "dimphy.h" |
---|
| 102 | INTEGER nvm |
---|
| 103 | PARAMETER (nvm=8) |
---|
| 104 | REAL veget(klon,nvm) |
---|
| 105 | REAL alb_neig(klon) |
---|
| 106 | REAL agesno(klon) |
---|
| 107 | c |
---|
| 108 | INTEGER i, nv |
---|
| 109 | c |
---|
| 110 | REAL init(nvm), decay(nvm), as |
---|
| 111 | SAVE init, decay |
---|
| 112 | DATA init /0.55, 0.14, 0.18, 0.29, 0.15, 0.15, 0.14, 0./ |
---|
| 113 | DATA decay/0.30, 0.67, 0.63, 0.45, 0.40, 0.14, 0.06, 1./ |
---|
| 114 | c |
---|
| 115 | DO i = 1, klon |
---|
| 116 | alb_neig(i) = 0.0 |
---|
| 117 | ENDDO |
---|
| 118 | DO nv = 1, nvm |
---|
| 119 | DO i = 1, klon |
---|
| 120 | as = init(nv)+decay(nv)*EXP(-agesno(i)/5.) |
---|
| 121 | alb_neig(i) = alb_neig(i) + veget(i,nv)*as |
---|
| 122 | ENDDO |
---|
| 123 | ENDDO |
---|
| 124 | c |
---|
| 125 | RETURN |
---|
| 126 | END |
---|