/[MITgcm]/MITgcm_contrib/bling/pkg/bling_light.F
ViewVC logotype

Contents of /MITgcm_contrib/bling/pkg/bling_light.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.2 - (show annotations) (download)
Sun Feb 28 21:49:24 2016 UTC (9 years, 4 months ago) by mmazloff
Branch: MAIN
Changes since 1.1: +41 -63 lines
Update to BLING version 2

1 C $Header: $
2 C $Name: $
3
4 #include "BLING_OPTIONS.h"
5
6 CBOP
7 subroutine BLING_LIGHT(
8 I mld,
9 U irr_inst, irr_eff,
10 I bi, bj, imin, imax, jmin, jmax,
11 I myIter, myTime, myThid )
12
13 C =================================================================
14 C | subroutine bling_light
15 C | o calculate effective light for phytoplankton growth
16 C | There are multiple types of light.
17 C | - irr_inst is the instantaneous irradiance field.
18 C | - irr_mix is the same, but with the irr_inst averaged throughout
19 C | the mixed layer. This quantity is intended to represent the
20 C | light to which phytoplankton subject to turbulent transport in
21 C | the mixed-layer would be exposed.
22 C | - irr_mem is a temporally smoothed field carried between
23 C | timesteps, to represent photoadaptation.
24 C | - irr_eff is the effective irradiance for photosynthesis,
25 C | given either by irr_inst or irr_mix, depending on model
26 C | options and location.
27 C =================================================================
28
29 implicit none
30
31 C === Global variables ===
32 C irr_inst :: Instantaneous irradiance
33 C irr_mem :: Phyto irradiance memory
34
35 #include "SIZE.h"
36 #include "EEPARAMS.h"
37 #include "PARAMS.h"
38 #include "FFIELDS.h"
39 #include "GRID.h"
40 #include "DYNVARS.h"
41 #include "BLING_VARS.h"
42 #include "PTRACERS_SIZE.h"
43 #include "PTRACERS_PARAMS.h"
44 #ifdef ALLOW_AUTODIFF
45 # include "tamc.h"
46 #endif
47
48 C === Routine arguments ===
49 C bi,bj :: tile indices
50 C iMin,iMax :: computation domain: 1rst index range
51 C jMin,jMax :: computation domain: 2nd index range
52 C myTime :: current time
53 C myIter :: current timestep
54 C myThid :: thread Id. number
55 INTEGER bi, bj, imin, imax, jmin, jmax
56 INTEGER myThid
57 INTEGER myIter
58 _RL myTime
59 C === Input ===
60 _RL mld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
61 C === Output ===
62 C irr_inst :: instantaneous light
63 C irr_eff :: effective light for photosynthesis
64 _RL irr_inst (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
65 _RL irr_eff (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
66
67 C === Local variables ===
68 _RL solar, albedo
69 _RL dayfrac, yday, delta
70 _RL lat, sun1, dayhrs
71 _RL cosz, frac, fluxi
72 _RL atten
73 _RL irr_surf (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
74 #ifdef ML_MEAN_LIGHT
75 _RL irr_mix (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
76 _RL SumMLIrr
77 _RL tmp_ML
78 #endif
79 #ifndef READ_PAR
80 #ifndef USE_QSW
81 _RL sfac (1-OLy:sNy+OLy)
82 #endif
83 #endif
84 integer i,j,k
85 CEOP
86
87 DO k=1,Nr
88 DO j=jmin,jmax
89 DO i=imin,imax
90 irr_eff(i,j,k) = 0. _d 0
91 ENDDO
92 ENDDO
93 ENDDO
94
95 c ---------------------------------------------------------------------
96 c Surface insolation
97
98 #ifndef USE_EXFQSW
99 c From pkg/dic/dic_insol
100 c find light as function of date and latitude
101 c based on paltridge and parson
102
103 solar = 1360. _d 0 !solar constant
104 albedo = 0.6 _d 0 !planetary albedo
105
106 C Case where a 2-d output array is needed: for now, stop here.
107 IF ( usingCurvilinearGrid .OR. rotateGrid ) THEN
108 STOP 'ABNORMAL END: S/R INSOL: 2-D output not implemented'
109 ENDIF
110
111 C find day (****NOTE for year starting in winter*****)
112 dayfrac=mod(myTime,360. _d 0*86400. _d 0)
113 & /(360. _d 0*86400. _d 0) !fraction of year
114 yday = 2. _d 0*PI*dayfrac !convert to radians
115 delta = (0.006918 _d 0
116 & -(0.399912 _d 0*cos(yday)) !cosine zenith angle
117 & +(0.070257 _d 0*sin(yday)) !(paltridge+platt)
118 & -(0.006758 _d 0*cos(2. _d 0*yday))
119 & +(0.000907 _d 0*sin(2. _d 0*yday))
120 & -(0.002697 _d 0*cos(3. _d 0*yday))
121 & +(0.001480 _d 0*sin(3. _d 0*yday)) )
122 DO j=1-OLy,sNy+OLy
123 C latitude in radians
124 lat=YC(1,j,1,bj)*deg2rad
125 C latitute in radians, backed out from coriolis parameter
126 C (makes latitude independent of grid)
127 IF ( usingCartesianGrid .OR. usingCylindricalGrid )
128 & lat = asin( fCori(1,j,1,bj)/(2. _d 0*omega) )
129 sun1 = -sin(delta)/cos(delta) * sin(lat)/cos(lat)
130 IF (sun1.LE.-0.999 _d 0) sun1=-0.999 _d 0
131 IF (sun1.GE. 0.999 _d 0) sun1= 0.999 _d 0
132 dayhrs = abs(acos(sun1))
133 cosz = ( sin(delta)*sin(lat)+ !average zenith angle
134 & (cos(delta)*cos(lat)*sin(dayhrs)/dayhrs) )
135 IF (cosz.LE.5. _d -3) cosz= 5. _d -3
136 frac = dayhrs/PI !fraction of daylight in day
137 C daily average photosynthetically active solar radiation just below surface
138 fluxi = solar*(1. _d 0-albedo)*cosz*frac*parfrac
139
140 C convert to sfac
141 sfac(j) = MAX(1. _d -5,fluxi)
142 ENDDO !j
143
144 #endif
145
146 c ---------------------------------------------------------------------
147 c instantaneous light, mixed layer averaged light
148
149 C$TAF LOOP = parallel
150 DO j=jmin,jmax
151 C$TAF LOOP = parallel
152 DO i=imin,imax
153
154 c Photosynthetically-available radiations (PAR)
155 #ifdef USE_EXFQSW
156 irr_surf(i,j) = max(epsln,
157 & -parfrac*Qsw(i,j,bi,bj)*maskC(i,j,1,bi,bj))
158 #else
159 irr_surf(i,j) = sfac(j)
160 #endif
161 cav IF ( .NOT. QSW_underice ) THEN
162 c if using Qsw but not seaice/thsice or coupled, then
163 c ice fraction needs to be taken into account
164 cav irr_surf(i,j) = irr_surf(i,j)*(1. _d 0 - FIce(i,j,bi,bj))
165 cav ENDIF
166
167 #ifdef ML_MEAN_LIGHT
168 SumMLIrr = 0. _d 0
169 tmp_ML = 0. _d 0
170 #endif
171
172 c C$TAF init ml_stuff = static, Nr
173 DO k=1,Nr
174 c C$TAF STORE SumMLDepth = ml_stuff
175
176 IF (hFacC(i,j,k,bi,bj).gt.0) THEN
177
178 IF (k.eq.1) THEN
179 c Light attenuation in middle of top layer
180 atten = k0*drF(1)/2. _d 0*hFacC(i,j,1,bi,bj)
181 irr_inst(i,j,1) = irr_surf(i,j)*exp(-atten)
182 ELSE
183 c Attenuation from one more layer
184 atten = k0*drF(k)/2. _d 0*hFacC(i,j,k,bi,bj)
185 & + k0*drF(k-1)/2. _d 0*hFacC(i,j,k-1,bi,bj)
186 irr_inst(i,j,k) =
187 & irr_inst(i,j,k-1)*exp(-atten)
188 ENDIF
189
190 #ifdef ML_MEAN_LIGHT
191 c Mean irradiance in the mixed layer
192 IF ((-rf(k+1) .le. mld(i,j)).and.
193 & (-rf(k+1).lt.200. _d 0)) THEN
194 SumMLIrr = SumMLIrr+drF(k)*irr_inst(i,j,k)
195 tmp_ML = tmp_ML + drF(k)
196 irr_mix(i,j) = SumMLIrr/tmp_ML
197 ENDIF
198 #endif
199
200 ENDIF
201
202 ENDDO
203 ENDDO
204 ENDDO
205
206
207 DO k=1,Nr
208 DO j=jmin,jmax
209 DO i=imin,imax
210
211 IF (hFacC(i,j,k,bi,bj) .gt. 0. _d 0) THEN
212
213 irr_eff(i,j,k) = irr_inst(i,j,k)
214 #ifdef ML_MEAN_LIGHT
215 c Inside mixed layer, effective light is set to mean mixed layer light
216 IF ((-rf(k+1) .le. mld(i,j)).and.
217 & (-rf(k+1).lt.200. _d 0)) THEN
218 irr_eff(i,j,k) = irr_mix(i,j)
219 ENDIF
220 #endif
221
222 ENDIF
223
224 ENDDO
225 ENDDO
226 ENDDO
227
228 #ifdef ALLOW_DIAGNOSTICS
229 IF ( useDiagnostics ) THEN
230 CALL DIAGNOSTICS_FILL(Qsw,'BLGQSW ',0,1,1,bi,bj,myThid)
231 CALL DIAGNOSTICS_FILL(irr_inst,'BLGIRRIS',0,Nr,2,bi,bj,myThid)
232 ENDIF
233 #endif
234
235 RETURN
236 END
237

  ViewVC Help
Powered by ViewVC 1.1.22