Blame view

src/extras/bio/bio_gsj.F90 24.7 KB
84ed20a5   Dany Dumont   Construction d'un...
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
!$Id: bio_ismer.F90,v 1.11 2007-01-06 11:49:15 kbk Exp $
#include"cppdefs.h"
!-----------------------------------------------------------------------
!BOP
!
! !MODULE: bio_gsj --- Modified from Fasham et al. biological model \label{sec:bio-fasham}
!
! !INTERFACE:
   module bio_gsj
!
! !DESCRIPTION:
!  The model developed by \cite{Fashametal1990} 
!  uses nitrogen as 'currency' according to the evidence that in
!  most cases nitrogen is the limiting macronutrient. It consists of
!  seven state variables: phytoplankton, zooplankton, bacteria,
!  particulate organic matter (detritus), dissolved organic matter
!  and the nutrients nitrate and ammonium.
!  The structure of the \cite{Fashametal1990} biogeochemical model
!  is given in figure \ref{fig_fasham}.
! \begin{figure}
! \begin{center}
! \scalebox{0.5}{\includegraphics{figures/fasham_structure.eps}}
! \caption{Structure of the \cite{Fashametal1990} model with bacteria (bac),
! phytoplankton (phy), detritus (det), zooplankton (zoo), labile dissolved
! organic nitrogen (don), ammonium (amm) and nitrate (nit) as the seven
! state variables.
! The concentrations are in mmol N\,m$^{-3}$,
! all fluxes (green arrows) are conservative.
! }\label{fig_fasham}
! \end{center}
! \end{figure}
!  A detailed mathematical description of all
!  processes is given in section \ref{sec:bio-fasham-rhs}.
!  The version of the \cite{Fashametal1990} model which is implemented includes
!  slight modifications by \cite{KuehnRadach1997} and has been 
!  included into GOTM by \cite{Burchardetal05}. 

! !USES:
!  default: all is private.
   use bio_var
   use output
   use observations, only : aa,g2
   private
!
! !PUBLIC MEMBER FUNCTIONS:
   public init_bio_gsj, init_var_gsj, var_info_gsj, &
          light_gsj, do_bio_gsj, end_bio_gsj
!
! !PRIVATE DATA MEMBERS:
!
! !REVISION HISTORY:!
!  Original author(s): Hans Burchard & Karsten Bolding
!
!
! !LOCAL VARIABLES:
!  from a namelist
   REALTYPE                  ::  p1_init = 0.05
   REALTYPE                  ::  p2_init = 0.05
   REALTYPE                  ::  z1_init = 0.05
   REALTYPE                  ::  z2_init = 0.05
   REALTYPE                  ::  b_init  = 0.001
   REALTYPE                  ::  d_init  = 0.4
   REALTYPE                  ::  l_init  = 0.14
   REALTYPE                  ::  p0      = 0.0
   REALTYPE                  ::  z0      = 0.0
   REALTYPE                  ::  b0      = 0.0
   LOGICAL                   ::  mte     = .true.
   REALTYPE                  ::  ca1     = 3.61
   REALTYPE                  ::  ca2     = 14.58    
   REALTYPE                  ::  ch1     = 3.265
   REALTYPE                  ::  ch2     = 24.923
   REALTYPE                  ::  amratio = 200.0
   REALTYPE                  ::  hmratio = 1000.0
   REALTYPE                  ::  vp1     = 1.5
   REALTYPE                  ::  alpha1  = 0.065
   REALTYPE                  ::  inib1   = 0.05
   REALTYPE                  ::  vp2     = 1.5
   REALTYPE                  ::  alpha2  = 0.065
   REALTYPE                  ::  inib2   = 0.05
   REALTYPE                  ::  theta   = 0.0
   REALTYPE                  ::  w_p1min = -0.06
   REALTYPE                  ::  w_p1max = -0.38
   REALTYPE                  ::  w_p2min = -0.06
   REALTYPE                  ::  w_p2max = -0.38
   REALTYPE                  ::  kn1     = 0.2
   REALTYPE                  ::  ka1     = 0.8
   REALTYPE                  ::  kn2     = 0.2
   REALTYPE                  ::  ka2     = 0.8
   REALTYPE                  ::  mu11    = 0.05
   REALTYPE                  ::  mu12    = 0.05
   REALTYPE                  ::  k5      = 0.2
   REALTYPE                  ::  gamma   = 0.05
   REALTYPE                  ::  w_p1    = -0.5
   REALTYPE                  ::  w_p2    = -0.5
   REALTYPE                  ::  g1max   = 1.0
   REALTYPE                  ::  g2max   = 1.0
   REALTYPE                  ::  k3      = 1.0
   REALTYPE                  ::  beta    = 0.625
   REALTYPE                  ::  mu21    = 0.3
   REALTYPE                  ::  mu22    = 0.3
   REALTYPE                  ::  k6      = 0.2
   REALTYPE                  ::  delta   = 0.1
   REALTYPE                  ::  epsi    = 0.70
   REALTYPE                  ::  r11     = 0.55
   REALTYPE                  ::  r12     = 0.30
   REALTYPE                  ::  r13     = 0.05
   REALTYPE                  ::  r14     = 0.10
   REALTYPE                  ::  r21     = 0.50
   REALTYPE                  ::  r22     = 0.30
   REALTYPE                  ::  r23     = 0.05
   REALTYPE                  ::  r24     = 0.15
63c62f83   Philippe Klotz   Colimitation par ...
112
113
   REALTYPE                  ::  vb1     = 1.2
   REALTYPE                  ::  vb2     = 1.2
84ed20a5   Dany Dumont   Construction d'un...
114
   REALTYPE                  ::  k4      = 0.5
63c62f83   Philippe Klotz   Colimitation par ...
115
   REALTYPE                  ::  k10     = 0.15
84ed20a5   Dany Dumont   Construction d'un...
116
117
   REALTYPE                  ::  w_h     = 0.0
   REALTYPE                  ::  mu3     = 0.15
63c62f83   Philippe Klotz   Colimitation par ...
118
119
   REALTYPE                  ::  etaa    = 0.0
   REALTYPE                  ::  etah    = 0.0
84ed20a5   Dany Dumont   Construction d'un...
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
   REALTYPE                  ::  mu4     = 0.02
   REALTYPE                  ::  mu5     = 0.00
   REALTYPE                  ::  w_d     = -2.0
   REALTYPE, public          ::  kc      = 0.03
   integer                   ::  out_unit
   integer, parameter        ::  n=1,p1=2,p2=3,z1=4,z2=5,d=6,l=7,b1=8,a=9,hc=10,b2=11
!EOP
!-----------------------------------------------------------------------

   contains

!-----------------------------------------------------------------------
!BOP
!
! !IROUTINE: Initialise the bio module
!
! !INTERFACE:
   subroutine init_bio_gsj(namlst,fname,unit)
!
! !DESCRIPTION:
!  Here, the bio namelist {\tt bio\_fasham.nml} is read and
!  various variables (rates and settling velocities)
!  are transformed into SI units.
!
! !USES:
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
   integer,          intent(in)   :: namlst
   character(len=*), intent(in)   :: fname
   character(len=20)              :: pfile
   integer,          intent(in)   :: unit
!
! !REVISION HISTORY:
!  Original author(s): Hans Burchard & Karsten Bolding
!
! !LOCAL VARIABLES:
   namelist /bio_gsj_nml/ numc, &
                        p1_init,p2_init,z1_init,z2_init,                 &
                        b_init,d_init,l_init,                            &
                        p0,z0,b0,vp1,alpha1,inib1,vp2,alpha2,inib2,      &
                        kn1,ka1,kn2,ka2,mu11,mu12,k5,gamma,w_p1,w_p2,    &
                        g1max,g2max,k3,beta,mu21,mu22,k6,delta,epsi,     &
                        r11,r12,r13,r14,r21,r22,r23,r24,                 &
63c62f83   Philippe Klotz   Colimitation par ...
164
                        vb1,vb2,k4,k10,w_h,mu3,etaa,etah,mu4,w_d,kc,mu5, &
84ed20a5   Dany Dumont   Construction d'un...
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
                        theta,w_p1max,w_p1min,w_p2min,w_p2max,           &
                        mte,ca1,ca2,ch1,ch2,amratio,hmratio

!EOP
!-----------------------------------------------------------------------
!BOC
   LEVEL2 'init_bio_gsj'

   open(namlst,file=fname,action='read',status='old',err=98)
   read(namlst,nml=bio_gsj_nml,err=99)
   close(namlst)

   numcc=numc

! Print some parameter values in standard output
! and save them in a separate file [out_fn]_ismer.par
   pfile = trim(out_fn) // '_gsj.par'
   open(10,status='unknown',action='write',file=pfile)
   LEVEL3 'GSJ parameters saved in ', pfile
   write(10,901) vp1
   write(10,901) vp2
   write(10,901) alpha1
   write(10,901) alpha2
   write(10,901) inib1
   write(10,901) inib2
   write(10,901) kn1
   write(10,901) ka1
   write(10,901) kn2
   write(10,901) ka2
   write(10,901) mu11
   write(10,901) mu12
   write(10,901) k5
   write(10,901) gamma
   write(10,901) w_p1
   write(10,901) w_p2
   write(10,901) g1max
   write(10,901) g2max
   write(10,901) k3 
   write(10,901) beta
   write(10,901) mu21
   write(10,901) mu22
   write(10,901) k6
   write(10,901) delta
   write(10,901) epsi
   write(10,901) r11
   write(10,901) r12
   write(10,901) r13
   write(10,901) r14
   write(10,901) r21
   write(10,901) r22
   write(10,901) r23
   write(10,901) r24
63c62f83   Philippe Klotz   Colimitation par ...
217
218
   write(10,901) vb1
   write(10,901) vb2
84ed20a5   Dany Dumont   Construction d'un...
219
220
   write(10,901) k4
   write(10,901) mu3
63c62f83   Philippe Klotz   Colimitation par ...
221
222
   write(10,901) etaa
   write(10,901) etah
84ed20a5   Dany Dumont   Construction d'un...
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
   write(10,901) mu4
   write(10,901) mu5
   write(10,901) w_d
   write(10,901) kc
   if (mte) then
      write(10,901) ca1
      write(10,901) ca2
      write(10,901) ch1
      write(10,901) ch2
      write(10,901) amratio
      write(10,901) hmratio
   endif
   

900 format (a,f8.5)
901 format (f8.5)

!  Conversion from day to second
   vp1     = vp1     /secs_pr_day
   vp2     = vp2     /secs_pr_day
63c62f83   Philippe Klotz   Colimitation par ...
243
244
   vb1     = vb1     /secs_pr_day
   vb2     = vb2     /secs_pr_day
84ed20a5   Dany Dumont   Construction d'un...
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
   mu11    = mu11    /secs_pr_day
   mu12    = mu12    /secs_pr_day
   mu21    = mu21    /secs_pr_day
   mu22    = mu22    /secs_pr_day
   mu3     = mu3     /secs_pr_day
   mu4     = mu4     /secs_pr_day
   mu5     = mu5     /secs_pr_day
   g1max   = g1max   /secs_pr_day
   g2max   = g2max   /secs_pr_day
   w_p1    = w_p1    /secs_pr_day
   w_p2    = w_p2    /secs_pr_day
   w_p1min = w_p1min /secs_pr_day
   w_p1max = w_p1max /secs_pr_day
   w_p2min = w_p2min /secs_pr_day
   w_p2max = w_p2max /secs_pr_day
   theta   = theta   /secs_pr_day
   w_d     = w_d     /secs_pr_day
c58569d0   Dany Dumont   Modification du t...
262
   w_h     = w_h     /secs_pr_day
84ed20a5   Dany Dumont   Construction d'un...
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
   alpha1  = alpha1  /secs_pr_day
   inib1   = inib1   /secs_pr_day
   alpha2  = alpha2  /secs_pr_day
   inib2   = inib2   /secs_pr_day

   out_unit=unit

   LEVEL3 'GSJ bio module initialised ...'

   return

98 LEVEL2 'I could not open bio_gsj.nml'
   LEVEL2 'If thats not what you want you have to supply bio_gsj.nml'
   LEVEL2 'See the bio example on www.gotm.net for a working bio_gsj.nml'
   return
99 FATAL 'I could not read bio_gsj.nml'
   stop 'init_bio_gsj'
   end subroutine init_bio_gsj
!EOC

!-----------------------------------------------------------------------
!BOP
!
! !IROUTINE: Initialise the concentration variables
!
! !INTERFACE:
   subroutine init_var_gsj(nlev)
!
! !DESCRIPTION:
!  Here, the the initial conditions are set and the settling velocities are
!  transferred to all vertical levels. All concentrations are declared
!  as non-negative variables, and it is defined which variables would be
!  taken up by benthic filter feeders.
!
! !USES:
   use observations,    only: nprof,aprof,hcprof
   use meanflow,        only: nit,amm,hcb,T,S

   IMPLICIT NONE

!
! !INPUT PARAMETERS:
   integer, intent(in)                 :: nlev
!
! !REVISION HISTORY:
!  Original author(s): Hans Burchard & Karsten Bolding

! !LOCAL VARIABLES:
  integer                    :: i
!EOP
!-----------------------------------------------------------------------
!BOC
   do i=1,nlev
      cc(n,i) = nprof(i)                                  !CHG3
      cc(p1,i)= p1_init
      cc(p2,i)= p2_init
      cc(z1,i)= z1_init
      cc(z2,i)= z2_init
      cc(d,i) = d_init
      cc(l,i) = l_init
      cc(b1,i)= b_init
      cc(a,i) = aprof(i)                                  !CHG5
      cc(hc,i)= hcprof(i)
      cc(b2,i)= b_init
   end do

   do i=0,nlev
      ws(n,i)  = _ZERO_
      ws(p1,i) = w_p1
      ws(p2,i) = w_p2
      ws(z1,i) = _ZERO_
      ws(z2,i) = _ZERO_
      ws(d,i)  = w_d
      ws(l,i)  = _ZERO_
      ws(b1,i) = _ZERO_
      ws(a,i)  = _ZERO_
c58569d0   Dany Dumont   Modification du t...
339
      ws(hc,i) = w_h
84ed20a5   Dany Dumont   Construction d'un...
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
      ws(b2,i) = _ZERO_
   end do

   posconc(n)  = 1
   posconc(p1) = 1
   posconc(p2) = 1
   posconc(z1) = 1
   posconc(z2) = 1
   posconc(d)  = 1
   posconc(l)  = 1
   posconc(b1) = 1
   posconc(a)  = 1
   posconc(hc) = 1
   posconc(b2) = 1

   LEVEL3 'GSJ variables initialised ...'

   return

   end subroutine init_var_gsj
!EOC

!-----------------------------------------------------------------------
!BOP
!
! !IROUTINE: Providing info on variables
!
! !INTERFACE:
   subroutine var_info_gsj()
!
! !DESCRIPTION:
!  This subroutine provides information about the variable names as they
!  will be used when storing data in NetCDF files.
!
! !USES:
   IMPLICIT NONE
!
! !REVISION HISTORY:
!  Original author(s): Hans Burchard & Karsten Bolding
!
! !LOCAL VARIABLES:
!EOP
!-----------------------------------------------------------------------
!BOC
   var_names(1) = 'nit'
   var_units(1) = 'mmol/m**3'
   var_long(1)  = 'nitrate'

   var_names(2) = 'fla'
   var_units(2) = 'mmol/m**3'
   var_long(2)  = 'flagellates'

   var_names(3) = 'dia'
   var_units(3) = 'mmol/m**3'
   var_long(3)  = 'diatoms'

   var_names(4) = 'mcz'
   var_units(4) = 'mmol/m**3'
   var_long(4)  = 'micro-zooplankton'

   var_names(5) = 'msz'
   var_units(5) = 'mmol/m**3'
   var_long(5)  = 'meso-zooplankton'

   var_names(6) = 'det'
   var_units(6) = 'mmol/m**3'
   var_long(6)  = 'detritus'

   var_names(7) = 'ldn'
   var_units(7) = 'mmol/m**3'
   var_long(7)  = 'labile_dissolved_organic_nitrogen'

   var_names(8) = 'bac1'
   var_units(8) = 'mmol/m**3'
   var_long(8)  = 'bacteria'

   var_names(9) = 'amm'
   var_units(9) = 'mmol/m**3'
   var_long(9)  = 'ammonium'

   var_names(10)= 'hcb'
   var_units(10)= 'mmol/m**3'
   var_long(10) = 'hydrocarbon'

   var_names(11)= 'bac2'
   var_units(11)= 'mmol/m**3'
   var_long(11) = 'bacteria'

   return
   end subroutine var_info_gsj
!EOC

!-----------------------------------------------------------------------
!BOP
!
! !IROUTINE: Light properties for the ISMER model
!
! !INTERFACE:
   subroutine light_gsj(nlev,bioshade_feedback)
!
! !DESCRIPTION:
! Here, the photosynthetically available radiation is calculated
! by simply assuming that the short wave part of the total
! radiation is available for photosynthesis. 
! The photosynthetically
! available radiation, $I_{PAR}$, follows from (\ref{light}).
! The user should make
! sure that this is consistent with the light class given in the
! {\tt extinct} namelist of the {\tt obs.nml} file.
! The self-shading effect is also calculated in this subroutine,
! which may be used to consider the effect of bio-turbidity also
! in the temperature equation (if {\tt bioshade\_feedback} is set
! to true in {\tt bio.nml}).
! For details, see section \ref{sec:do-bio}.
!
! !USES:
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
  integer, intent(in)                  :: nlev
  logical, intent(in)                  :: bioshade_feedback
!
! !REVISION HISTORY:
!  Original author(s): Hans Burchard, Karsten Bolding
!
! !LOCAL VARIABLES:
   integer                   :: i
   REALTYPE                  :: zz,add
!EOP
!-----------------------------------------------------------------------
!BOC
   zz = _ZERO_
   add = _ZERO_
   do i=nlev,1,-1
      add=add+0.5*h(i)*(cc(p1,i)+cc(p2,i)+2*p0)
      zz=zz+0.5*h(i)
      par(i)=rad(nlev)*(1.-aa)*exp(-zz/g2)*exp(-kc*add)
      add=add+0.5*h(i)*(cc(p1,i)+cc(p2,i)+2*p0)
      zz=zz+0.5*h(i)
      if (bioshade_feedback) bioshade_(i)=exp(-kc*add)
   end do


   return
   end subroutine light_gsj
!EOC

!-----------------------------------------------------------------------
!BOP
!
! !IROUTINE: Right hand sides of geobiochemical model \label{sec:bio-fasham-rhs}
!
! !INTERFACE:
   subroutine do_bio_gsj(first,numc,nlev,cc,pp,dd)
!
! !DESCRIPTION:
! 
! The \cite{Fashametal1990} model consisting of the $I=7$
! state variables phytoplankton, bacteria, detritus, zooplankton, 
! nitrate, ammonium and dissolved organic nitrogen is described here
! in detail.
! 
! Phytoplankton mortality and zooplankton grazing loss of phytoplankton:
! \begin{equation}\label{d13}
! d_{1,3} = \mu_1 \frac{c_1+c_{1}^{\min}}{K_5+c_1+c_{1}^{\min}}c_1+
! (1-\beta)\frac{g\rho_1 c_1^2}{K_3 \sum_{j=1}^3 \rho_jc_j
! + \sum_{j=1}^3 \rho_jc_j^2} (c_4+c_{4}^{\min}).
! \end{equation}
! Phytoplankton loss to LDON (labile dissolved organic nitrogen):
! \begin{equation}\label{d17}
! d_{1,7} = \gamma
! F(I_{PAR})\frac{\frac{c_5}{K_1}
! +\frac{c_6}{K_2}}{1+\frac{c_5}{K_1}+\frac{c_6}{K_2}}c_1,
! \end{equation}
! with
! \begin{equation}\label{FI}
!  F(I_{PAR}) = \frac{V_p\alpha I_{PAR}(z)}{\left(V_p^2+\alpha^2(I_{PAR}(z))^2 
! \right)^{1/2}}.
! \end{equation}
! With $I_{PAR}$ from (\ref{light}). 
! 
! Zooplankton grazing loss:
! \begin{equation}\label{di3}
! d_{2,3} = (1-\beta)\frac{g\rho_2 c_2^2}{K_3 \sum_{j=1}^3 \rho_jc_j 
! + \sum_{j=1}^3 \rho_jc_j^2} (c_4+c_{4}^{\min}).
! \end{equation}
! Zooplankton grazing:
! \begin{equation}\label{di4}
! d_{i,4} = \beta\frac{g\rho_i c_i^2}{K_3 \sum_{j=1}^3 \rho_jc_j 
! + \sum_{j=1}^3 \rho_jc_j^2} (c_4+c_{4}^{\min}), \quad i=1,\dots,3.
! \end{equation}
! Bacteria excretion rate:
! \begin{equation}\label{d26}
! d_{2,6} = \mu_3 c_2.
! \end{equation}
! Detritus breakdown rate:
! \begin{equation}\label{d37}
! d_{3,7} = \mu_4 c_3.
! \end{equation}
! Zooplankton losses to detritus, ammonium and LDON:
! \begin{equation}\label{d43}
! d_{4,3} = (1-\epsilon-\delta)\mu_2 
! \frac{c_4+c_{4}^{\min}}{K_6+c_4+c_{4}^{\min}}c_4.
! \end{equation}
! \begin{equation}\label{d46}
! d_{4,6} = \epsilon\mu_2 \frac{c_4+c_{4}^{\min}}{K_6+c_4+c_{4}^{\min}}c_4.
! \end{equation}
! \begin{equation}\label{d47}
! d_{4,7} = \delta\mu_2 \frac{c_4+c_{4}^{\min}}{K_6+c_4+c_{4}^{\min}}c_4.
! \end{equation}
! Nitrate uptake by phytoplankton:
! \begin{equation}\label{d51}
! d_{5,1} = F(I_{PAR})\frac{\frac{c_5}{K_1}}{1+\frac{c_5}{K_1}
! +\frac{c_6}{K_2}}(c_1+c_{1}^{\min}).
! \end{equation}
! Ammonium uptake by phytoplankton:
! \begin{equation}\label{d61}
! d_{6,1} = F(I_{PAR})\frac{\frac{c_6}{K_2}}{1+\frac{c_5}{K_1}
! +\frac{c_6}{K_2}}(c_1+c_{1}^{\min}).
! \end{equation}
! Ammonium uptake by bacteria:
! \begin{equation}\label{d62}
! d_{6,2} = V_b \frac{\min(c_6,\eta c_7)}{K_4+\min(c_6,\eta c_7)+c_7} 
! (c_2+c_{2}^{\min}).
! \end{equation}
! LDON uptake by bacteria:
! \begin{equation}\label{d72}
! d_{7,2} = V_b \frac{c_7}{K_4+\min(c_6,\eta c_7)+c_7} (c_2+c_{2}^{\min}).
! \end{equation}
!
! !USES:
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
   logical, intent(in)                 :: first
   integer, intent(in)                 :: numc,nlev
   REALTYPE, intent(in)                :: cc(1:numc,0:nlev)
!
! !OUTPUT PARAMETERS:
   REALTYPE, intent(out)               :: pp(1:numc,1:numc,0:nlev)
   REALTYPE, intent(out)               :: dd(1:numc,1:numc,0:nlev)
!
! !REVISION HISTORY:
!  Original author(s): Hans Burchard, Karsten Bolding
!
! !LOCAL VARIABLES:
   REALTYPE                   :: amr1,amr2
   REALTYPE                   :: hmr1,hmr2
   REALTYPE                   :: fac1,fac2,fac3,fac4
63c62f83   Philippe Klotz   Colimitation par ...
589
   REALTYPE                   :: minal,minhl,qn1,qa1,qn2,qa2
84ed20a5   Dany Dumont   Construction d'un...
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
   REALTYPE                   :: ps1,ps2,ff1,ff2
   REALTYPE                   :: Ea,Eh,kBeV,T0
   integer                    :: i,j,ci
!EOP
!-----------------------------------------------------------------------
!BOC
!KBK - is it necessary to initialise every time - expensive in a 3D model
   pp = _ZERO_
   dd = _ZERO_

   do ci=1,nlev

      if (mte) then
         ! Temperature-dependent metabolic rates
         ! Gillooly et al. (2002)
         Eh = 0.65        ! Activation energy for heterotrophs (eV)
         Ea = 0.32        ! Activation energy for autotrophs (eV)
         kBeV = 8.62e-5   ! Boltzmann constant (eV K-1)
         T0 = 273.15-1.9  ! Temperature at which metabolism stops
      
         ! metabolic rate coefficient (a(T0) in Gillooly et al. 2002 in kg^1/4 s^-1)
         ! they are tuned to 
         ca1 = 3.61       ! pico-phytoplankton
         ca2 = 14.58      ! micro-phytoplankton	
         ch1 = 3.265      ! micro-zooplankton
         ch2 = 24.923     ! meso-zooplankton

         ! The mass ratio between pico-phytoplankton (p1) and
         ! nano- and micro-phytoplankton (p2) is amratio,
         ! which translates into a ratio of fac3=amratio**0.25 between
         ! their respective metabolic rates, according to the MTE.
         ! The same applies to micro-zooplankton (z1) versus
         ! meso-zooplankton (z2) through fac4=hmratio**0.25.

         fac3 = amratio**0.25
         fac4 = hmratio**0.25

         ! Autotroph metabolic rate
         amr1 = max(0.0,ca1*0.25*exp(Ea/(kBeV*T0**2)*(T(ci)/(1+T(ci)/T0)))     )
         amr2 = max(0.0,ca2*0.25*exp(Ea/(kBeV*T0**2)*(T(ci)/(1+T(ci)/T0)))/fac3)

         ! Heterotroph metabolic rate
         hmr1 = max(0.0,ch1*0.25*exp(Eh/(kBeV*T0**2)*(T(ci)/(1+T(ci)/T0)))     )
         hmr2 = max(0.0,ch2*0.25*exp(Eh/(kBeV*T0**2)*(T(ci)/(1+T(ci)/T0)))/fac4)
      else
         amr1 = 1.0
         amr2 = 1.0
         hmr1 = 1.0
         hmr2 = 1.0
      endif

      ! Smith (1936) - saturation (default)
      ! ff= vp*alpha*par(ci)/sqrt(vp**2+alpha**2*par(ci)**2)
      ! Blackman (1919)
      ! if (par(ci) .lt. vp/alpha) then
      !    ff=alpha*par(ci)
      ! else
      !    ff=vp
      ! endif       
      ! Steele (1962) - inhibition
      !    ff= vp*((par(ci)/I_opt)*exp(1-(par(ci)/I_opt)))
      ! Parker (1974) - inhibition
      !    ff= vp*((par(ci)/I_opt)*exp(1-(par(ci)/I_opt)))**2
   
      ! Platt et al. (1980) - inhibition
      ! Light limitation factor (PI curve)
      ps1= vp1/((alpha1/(alpha1+inib1))*(alpha1/(alpha1+inib1))**(inib1/alpha1))
      ff1= ps1*(1.-exp(-1.*alpha1*par(ci)/ps1))*exp(-1.*inib1*par(ci)/ps1)

      ps2= vp2/((alpha2/(alpha2+inib2))*(alpha2/(alpha2+inib2))**(inib2/alpha2))
      ff2= ps2*(1.-exp(-1.*alpha2*par(ci)/ps2))*exp(-1.*inib2*par(ci)/ps2)

      ! Nutrient limitation factors
      qn1=(cc(n,ci)/kn1)/(1.+cc(n,ci)/kn1+cc(a,ci)/ka1)
      qa1=(cc(a,ci)/ka1)/(1.+cc(n,ci)/kn1+cc(a,ci)/ka1)

      qn2=(cc(n,ci)/kn2)/(1.+cc(n,ci)/kn2+cc(a,ci)/ka2)
      qa2=(cc(a,ci)/ka2)/(1.+cc(n,ci)/kn2+cc(a,ci)/ka2)

      ! Grazing preference normalization factors
      fac1=(cc(z1,ci)+z0)/(k3*(r11*cc(p1,ci)+r12*cc(p2,ci)+0.5*r13*(cc(b1,ci)+cc(b2,ci))+r14*cc(d,ci))  &
                     +r11*cc(p1,ci)**2+r12*cc(p2,ci)**2+0.5*r13*(cc(b1,ci)+cc(b2,ci))**2+r14*cc(d,ci)**2)

      fac2=(cc(z2,ci)+z0)/(k3*(r21*cc(p1,ci)+r22*cc(p2,ci)                &
                     +r23*cc(d,ci)+r24*cc(z1,ci))                         &
                     +r21*cc(p1,ci)**2+r22*cc(p2,ci)**2                   &
                     +r23*cc(d,ci)**2+r24*cc(z1,ci)**2)

63c62f83   Philippe Klotz   Colimitation par ...
678
679
      minal=min(cc(a,ci),etaa*cc(l,ci))
      minhl=min(cc(hc,ci),etah*cc(l,ci))
84ed20a5   Dany Dumont   Construction d'un...
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729

      ! Light and nutrient limitation factors
      lumlim1(ci) =amr1*ff1
      nitlim1(ci) =qn1
      ammlim1(ci) =qa1
      lumlim2(ci) =amr2*ff2
      nitlim2(ci) =qn2
      ammlim2(ci) =qa2
      
      ! Nutrient uptake by pico- and nano-phytoplankton 
      dd(n,p1,ci) =amr1*ff1*qn1*(cc(p1,ci)+p0)
      dd(a,p1,ci) =amr1*ff1*qa1*(cc(p1,ci)+p0)
      dd(n,p2,ci) =amr2*ff2*qn2*(cc(p2,ci)+p0)
      dd(a,p2,ci) =amr2*ff2*qa2*(cc(p2,ci)+p0)

      dd(p1,l,ci) =amr1*gamma*ff1*(qn1+qa1)*cc(p1,ci)
      dd(p2,l,ci) =amr2*gamma*ff2*(qn2+qa2)*cc(p2,ci)

      dd(p1,d,ci) =amr1*mu11*(cc(p1,ci)+p0)/(k5+cc(p1,ci)+p0)*cc(p1,ci)  &
                   +(1.-beta)*cc(p1,ci)**2*(hmr1*g1max*r11*fac1+hmr2*g2max*r21*fac2)
      dd(p2,d,ci) =amr2*mu12*(cc(p2,ci)+p0)/(k5+cc(p2,ci)+p0)*cc(p2,ci)  &
                   +(1.-beta)*cc(p2,ci)**2*(hmr1*g1max*r12*fac1+hmr2*g2max*r22*fac2)

      dd(b1,d,ci) =(1.-beta)*g1max*r13*cc(b1,ci)**2*fac1

      dd(p1,z1,ci)=hmr1*beta*g1max*r11*cc(p1,ci)**2*fac1
      dd(p2,z1,ci)=hmr1*beta*g1max*r12*cc(p2,ci)**2*fac1
      dd(b1,z1,ci)=hmr1*beta*g1max*0.5*r13*cc(b1,ci)**2*fac1
      dd(b2,z1,ci)=hmr1*beta*g1max*0.5*r13*cc(b2,ci)**2*fac1
      dd(d,z1,ci) =hmr1*beta*g1max*r14*cc(d,ci)**2*fac1

      dd(p1,z2,ci)=hmr2*beta*g2max*r21*cc(p1,ci)**2*fac2
      dd(p2,z2,ci)=hmr2*beta*g2max*r22*cc(p2,ci)**2*fac2
      dd(d,z2,ci) =hmr2*beta*g2max*r23*cc(d,ci)**2*fac2
      dd(z1,z2,ci)=hmr2*beta*g2max*r24*cc(z1,ci)**2*fac2

      dd(b1,a,ci) =mu3*cc(b1,ci)
      dd(b2,a,ci) =mu3*cc(b2,ci)
      dd(d,l,ci)  =mu4*cc(d,ci)
      dd(a,n,ci)  =mu5*cc(a,ci)

      dd(z1,d,ci) =(1.-beta)*hmr2*g2max*r24*cc(z1,ci)**2*fac2            &
                   +hmr1*(1.-epsi-delta)*mu21*(cc(z1,ci)+z0)/(k6+cc(z1,ci)+z0)*cc(z1,ci)
      dd(z1,a,ci) =hmr1*epsi*mu21*(cc(z1,ci)+z0)/(k6+cc(z1,ci)+z0)*cc(z1,ci)
      dd(z1,l,ci) =hmr1*delta*mu21*(cc(z1,ci)+z0)/(k6+cc(z1,ci)+z0)*cc(z1,ci)

      dd(z2,d,ci) =hmr2*(1.-epsi-delta)*mu22*(cc(z2,ci)+z0)/(k6+cc(z2,ci)+z0)*cc(z2,ci)
      dd(z2,a,ci) =hmr2*epsi*mu22*(cc(z2,ci)+z0)/(k6+cc(z2,ci)+z0)*cc(z2,ci)
      dd(z2,l,ci) =hmr2*delta*mu22*(cc(z2,ci)+z0)/(k6+cc(z2,ci)+z0)*cc(z2,ci)

63c62f83   Philippe Klotz   Colimitation par ...
730
731
732
733
734
      dd(a,b1,ci) =vb1*minal/(k4+minal+cc(l,ci))*(cc(b1,ci)+b0)
      dd(l,b1,ci) =vb1*cc(l,ci)/(k4+minal+cc(l,ci))*(cc(b1,ci)+b0)
      dd(hc,b2,ci)=vb2*minhl/(k10+minhl+cc(l,ci))*(cc(b2,ci)+b0) 
      dd(l,b2,ci) =vb2*cc(l,ci)/(k10+minhl+cc(l,ci))*(cc(b2,ci)+b0)
      !dd(hc,b2,ci)=vb*cc(hc,ci)/(k10+cc(hc,ci))*(cc(b2,ci)+b0) 
84ed20a5   Dany Dumont   Construction d'un...
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779


      do i=1,numc
         do j=1,numc
            pp(i,j,ci)=dd(j,i,ci)
         end do
      end do
   end do

   return
   end subroutine do_bio_gsj
!EOC

!-----------------------------------------------------------------------
!BOP
!
! !IROUTINE: Finish the bio calculations
!
! !INTERFACE:
   subroutine end_bio_gsj
!
! !DESCRIPTION:
!  Nothing done yet --- supplied for completeness.
!
! !USES:
   IMPLICIT NONE
!
! !REVISION HISTORY:
!  Original author(s): Hans Burchard & Karsten Bolding
!
!EOP
!-----------------------------------------------------------------------
!BOC

   return
   end subroutine end_bio_gsj
!EOC

!-----------------------------------------------------------------------

   end module bio_gsj

!-----------------------------------------------------------------------
! Copyright by the GOTM-team under the GNU Public License - www.gnu.org
!-----------------------------------------------------------------------