kpp.F90 61.3 KB
Newer Older
dumoda01's avatar
dumoda01 committed
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
112
113
114
115
116
117
118
119
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
164
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
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
!$Id: kpp.F90,v 1.4 2007-01-06 11:49:15 kbk Exp $
#include"cppdefs.h"
!-----------------------------------------------------------------------
!BOP
!
! !MODULE: kpp: the KPP-turbulence model \label{sec:kpp}
!
! !INTERFACE:
   module kpp
!
! !DESCRIPTION:
! This implentation of the KPP turbulence parameterisation is based on the
! publications of \cite{Largeetal94} and \cite{Durskietal2004}.
! The general expression for the turbulent fluxes used in the KPP model is identical to
! that suggested in \eq{fluxes}. It assumes that the turbulent flux is the sum of a
! down-gradient flux and a non-local contribution,
! \begin{equation}
!   \label{kppFluxes}
!   \mean{w'\phi'} = - \nu_t^\phi  \partder{\mean{\phi}}{z}  + \tilde{\Gamma}_\phi
!  \comma
! \end{equation}
! where the super- or subscript $\phi$ is a placeholder for the symbols $m$, $h$, and $s$, indicating
! whether a quantity relates to momentum, heat, or salinity (or any other tracer), respectively.
! Note that turbulence parameters due to salinity stratification are updated only if the
! pre-processor macro {\tt KPP\_SALINITY} has been defined in {\tt cppdefs.h}.
!
! In the notation of the KPP model, the non-local flux is expressed as
! \begin{equation}
!   \label{kppNonlocal}
!   \tilde{\Gamma}_\phi = \nu_t^\phi \gamma_\phi
!  \comma
! \end{equation}
! where independent models are valid for $\nu_t^\phi$ and $\gamma_\phi$.
! The KPP model assumes that the turbulent diffusivity, $\nu_t^\phi$, inside the surface or
! bottom boundary layer is determined by a relation of the form
! \begin{equation}
!   \label{Kx}
!   \nu_t^\phi = h w_\phi(\sigma) G(\sigma)
!  \comma
! \end{equation}
! where $h$ denotes the thickness of the boundary layer, computed according
! to the algorithm discussed below.
! The non-dimensional boundary layer coordinate $\sigma$ is defined according to
! \begin{equation}
!   \label{defSigma}
!   \sigma = \dfrac{d}{h}
!  \comma
! \end{equation}
! where $d$ is the distance from the free surface (or the bottom boundary).
! The velocity scale, $w_\phi$, in \eq{Kx} is computed as described in \sect{sec:wscale}. The
! dimensionless shape function $G$ is a cubic polynomial,
! \begin{equation}
!   \label{GsigmaA}
!   G(\sigma) = a_0 + a_1 \sigma + a_2 \sigma^2 + a_3 \sigma^3
! \point
! \end{equation}
! Physical arguments discussed in \cite{Largeetal94} require  $a_0=0$, $a_1=1$. The remaining
! two parameters $a_2$ and $a_3$ may be re-expressed in terms of the value of $G$ and its
! derivative, $G'$,  at the edge of the boundary layer, $\sigma=1$.
! Then, \eq{GsigmaA} can be re-expressed as
! \begin{equation}
!   \label{GsigmaB}
!   G(\sigma) = \sigma \Big[
!   1 + \sigma \Big(
!     \left( \sigma - 2 \right)
!   + \left( 3 - 2 \sigma  \right) G(1)
!   + \left( \sigma - 1 \right) G'(1)
!   \Big) \Big]
! \end{equation}
! Apart from the boundary layer diffusivities, the KPP model also computes "interior"
! diffusivities, here denoted by $\nu_{ti}^{\phi}$. The function $G$ and its derivative
! can be evaluted from the requirement that, at the edge of the boundary layer, the
! boundary layer diffusivity and its
! derivative correspond exactly to the interior diffusivity and its derivative,
! respectively.
!
! Continuity of the boundary and interior diffusivites is obviously insured, see
! \eq{Kx}, if we require that
! \begin{equation}
!   \label{GOfOne}
!   G(1) = \dfrac{1}{h w_\phi(1)} \nu_{ti}^\phi(z_{bl})
!  \comma
! \end{equation}
! where $z_{bl}$ denotes the vertical coordinate of the surface (or bottom) boundary layer.
!
! A condition for the continuity of the derivatives of $\nu_t^\phi$ and $\nu_{ti}^\phi$
! can be obtained by carrying out the derivative with respect to $z$ of \eq{Kx}, and
! setting it equal to the $z$-derivative of $\nu_{ti}^\phi$. For the surface layer
! this results in
! \begin{equation}
!   \label{GPrimeOfOne}
!   G'(1) = - \dfrac{G(1)}{w(1)} \partder{w}{\sigma} \Big|_{\sigma=1}
!           - \dfrac{1}{w(1)} \partder{\nu_{ti}^\phi}{z} \Big|_{z=z_{bl}}
! \comma
! \end{equation}
! where we used the relation
! \begin{equation}
!   \label{Dsigma}
!   \partder{}{z} = - \frac{1}{h} \partder{}{\sigma}
! \comma
! \end{equation}
! if the motion of the free surface is ignored.
!
! The derivative of $w_\phi$ appearing in \eq{GPrimeOfOne} can be evaluted with
! the help of the formulae given in \sect{sec:wscale}. As discussed in \sect{sec:wscale},
! at $\sigma=1$, the derivative of $w_\phi$ is different from zero only for stably
! stratified flows. Then, the non-dimensional function $\Phi_\phi$ appearing
! in \eq{wscale} is given by \eq{PhiStable}, and it is easy to show that
! \begin{equation}
!   \label{dWdS}
!     \partder{w}{\sigma} \Big|_{\sigma=1} =
!    - 5 h w_\phi(1) \dfrac{B_f}{u_*^4}
! \comma
! \end{equation}
! valid for both, bottom and surface boundary layers.
! Note that in the original publication of \cite{Largeetal94}, erroneously,
! there appears an additional factor $\kappa$ in this relation.
!
! With the help of \eq{dWdS}, one can re-write \eq{GPrimeOfOne} as
! \begin{equation}
!   \label{GPrimeOfOneB}
!   G'(1) = \dfrac{B_f}{u_*^4} \nu_{ti}^\phi \Big|_{z=z_{bl}}
!           - \dfrac{1}{w(1)} \partder{\nu_{ti}^\phi}{z} \Big|_{z=z_{bl}}
! \comma
! \end{equation}
! valid only for the surface boundary layer. For the bottom boundary layer,
! the minus sign in \eq{Dsigma} disappears, with the consequence that
! the minus sign in \eq{GPrimeOfOneB} has to be replaced by a plus. Note that
! if the pre-processor macro {\tt KPP\_CLIP\_GS} is defined in {\tt cppdef.h},
! the slope of  $G$ is set to zero for negative slopes. For stably stratified
! flows with a stabilizing buoyancy flux, this limiter breaks the continuity of
! the first derivatives.
!
! The non-local transport term defined in \eq{kppNonlocal} is computed
! as described in \cite{Largeetal94}, if the pre-processor macro
! {\tt NONLOCAL} is defined. Otherwise, non-local transport is ignored.

! The position of the surface boundary layer depth, $z_{bl}$, corresponds
! to the position where the bulk Richardson number,
! \begin{equation}
!   \label{Rib}
!   Ri_b = \dfrac{(B_r - B(z_{bl})) d}
!                {\magn{\V U_r - \V U(z_{bl})}^2 + V_t^2(z_{bl})}
! \comma
! \end{equation}
! defined by \cite{Largeetal94}, reaches the critical value $Ri_c$.
! The subscript "r" in \eq{Rib} denotes a certain reference value of the buoyancy
! and velocity close to the surface. The choice of this reference value is
! not unique, and several possibilities have been implemented in numerical
! models. Presently, GOTM uses the uppermost grid point as the reference value.
! The bulk Richardson-number is then computed at the grid faces by linear
! interpolation of quantities defined at the centers (if
! ${\tt KPP\_TWOPOINT\_REF}$ is defined) or by simply identifying the 
! neighbouring center-value with the value at the face.
! The "turbulent velocity shear", $V_t$, is computed as described
! by \cite{Largeetal94}. The value of $z_{bl}$ is then found from
! \eq{Rib} by linear interpolation.
!
! To check the boundary layer limit according to the condition
! \begin{equation}
!   \label{RiConditionA}
!   Ri_b(z_{bl})
!   = \dfrac{Ri_\text{top}(z_{bl})}{Ri_\text{bot}(z_{bl})} = Ri_c
! \comma
! \end{equation}
! two methods have been implemented in GOTM. The first method simply evaluates
! \eq{RiConditionA} with a linear interpolation scheme.
! The second method is activated if the pre-processor
! macro ${\tt KPP\_IP\_FC}$ is defined. Then, the condition \eq{RiConditionA}
! is reformulated as
! \begin{equation}
!   \label{RiConditionB}
!   F_c(z_{bl}) = Ri_\text{top}(z_{bl})  - Ri_c Ri_\text{bot}(z_{bl}) = 0
! \point
! \end{equation}
! The position where the function $F_c$ changes sign is computed
! from linear interpolation.  This method has been suggested in the ROMS
! code as the numerically more stable alternative.
! Clearly, all approaches are grid-depending, a difficulty that cannot
! be overcome with the KPP model.
!
! Finally, provided {\tt clip\_mld=.true.} in {\tt kpp.nml}, the boundary layer is cut
! if it exceeds the Ekman or the Monin-Obukhov length scale, see \cite{Largeetal94}.
!
! !USES:

  use turbulence,   only: num,nuh,nus
  use turbulence,   only: gamu,gamv,gamh,gams
  use turbulence,   only: Rig
  use turbulence,   only: kappa

#ifdef EXTRA_OUTPUT
  use turbulence,   only: turb1,turb2,turb3,turb4,turb5
#endif

  use eqstate

  IMPLICIT NONE

  private

! !PUBLIC MEMBER FUNCTIONS:
!
  public init_kpp, do_kpp

! !PUBLIC DATA MEMBERS:
!
! z-position of surface boundary layer depth
  REALTYPE, public                 :: zsbl

! z-position of bottom boundary layer depth
  REALTYPE, public                 :: zbbl


! !DEFINED PARAMETERS:
!
!  non-dimensional extent of the surface layer (epsilon=0.1)
   REALTYPE, parameter :: epsilon = 0.1

!  critical gradient Richardson number below which turbulent
!  mixing occurs (Ri0=0.7)
   REALTYPE, parameter :: Ri0     = 0.7

!  value of double-diffusive density ratio where mixing goes
!  to zero in salt fingering (Rrho0=1.9)
   REALTYPE, parameter :: Rrho0    = 1.9

!  buoancy frequency (1/s2) limit for convection (bvfcon=-2.0E-5)
   REALTYPE, parameter :: bvfcon  = -2.0E-5

!  scaling factor for double diffusion of temperature in salt
!  fingering case (fdd=0.7)
   REALTYPE, parameter :: fdd     = 0.7

!  maximum interior convective viscosity and diffusivity
!  due to shear instability (nu0c=0.01)
   REALTYPE, parameter :: nu0c    = 0.01

!  maximum interior viscosity (m2/s) due to shear
!  instability (nu0m=10.0E-4)
   REALTYPE, parameter :: nu0m    = 10.0E-4

!  maximum interior diffusivity (m2/s) due to shear
!  instability (nu0s=10.0E-4)
   REALTYPE, parameter :: nu0s    = 10.0E-4

!  scaling factor for double diffusion in salt fingering (nu=1.5E-6)
   REALTYPE, parameter :: nu      = 1.5E-6

!  scaling factor for double diffusion in salt fingering (nuf=10.0E-4)
   REALTYPE, parameter :: nuf     = 10.0E-4

!  interior viscosity (m2/s) due to wave breaking (nuwm=1.0E-5)
   REALTYPE, parameter :: nuwm    = 1.0E-5

!  interior diffusivity (m2/s) due to wave breaking (nuwm=1.0E-6)
257
   REALTYPE, parameter :: nuws    = 1.0E-6
258
259
260
!   REALTYPE, parameter :: nuws    = 3.0E-6  ! Bourgault et al. (2011) - Amundsen Gulf
!   REALTYPE, parameter :: nuws    = 1.4E-5  ! Cyr et al. (submitted)  - LSLE 80m-bottom
!   REALTYPE, parameter :: nuws    = 4.3E-5  ! Cyr et al. (submitted)  - LSLE 30m-80m
261
!   REALTYPE, parameter :: nuws    = 5.5E-5  ! Cyr et al. (submitted)  - LSLE 30m-80m
dumoda01's avatar
dumoda01 committed
262
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
339
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
589
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
678
679
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
730
731
732
733
734
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

!  double diffusion constant for salinity in diffusive
!  convection case (sdd1=0.15)
   REALTYPE, parameter :: sdd1    = 0.15

!  double diffusion constant for salinity in diffusive convection
!  (sdd2=1.85)
   REALTYPE, parameter :: sdd2    = 1.85

!  double diffusion constant for salinity in diffusive convection
!  (sdd3=0.85)
   REALTYPE, parameter :: sdd3    = 0.85

!  double diffusion constant for temperature in diffusive convection
!  (tdd1=0.909)
   REALTYPE, parameter :: tdd1    = 0.909

!  double diffusion constant for temperature in diffusive convection
!  (tdd2=4.6)
   REALTYPE, parameter :: tdd2    = 4.6

!  double diffusion constant for temperature in diffusive convection case
!  (tdd3=0.54).
   REALTYPE, parameter :: tdd3    = 0.54

!  proportionality coefficient parameterizing nonlocal  transport
!  (Cstar=10.0)
   REALTYPE, parameter :: Cstar   = 10.0

!  ratio of interior Brunt-Vaisala frequency to that
!  at entrainment depth (Cv=1.5-1.6)
   REALTYPE, parameter :: Cv      = 1.6

!  ratio of entrainment flux to surface buoyancy flux (betaT=-0.2)
   REALTYPE, parameter :: betaT   = -0.2

!  constant for computation of Ekman scale (cekman=0.7)
   REALTYPE, parameter :: cekman  = 0.7

!  constant for computation of Monin-Obukhov scale (cmonob  = 1.0)
   REALTYPE, parameter :: cmonob  = 1.0

!  coefficient of flux profile for momentum in their
!  1/3 power law regimes (am=1.26)
   REALTYPE, parameter :: am      = 1.257

!  coefficient of flux profile for momentum in their
!  1/3 power law regimes (as=-28.86)
   REALTYPE, parameter :: as      = -28.86

!  coefficient of flux profile for momentum in their
!  1/3 power law regimes (cm=8.38)
   REALTYPE, parameter :: cm      = 8.38

!  coefficient of flux profile for momentum in their
!  1/3 power law regimes (cs=98.96)
   REALTYPE, parameter :: cs      = 98.96

!  maximum stability parameter "zeta" value of the 1/3
!  power law regime of flux profile for momentum (zetam=-0.2)
   REALTYPE, parameter :: zetam   = -0.2

!  maximum stability parameter "zeta" value of the 1/3
!  power law regime of flux profile for tracers (zetas=-1.0)
   REALTYPE, parameter :: zetas   = -1.0

! !REVISION HISTORY:
!  Original author(s): Lars Umlauf
!
!  $Log: kpp.F90,v $
!  Revision 1.4  2007-01-06 11:49:15  kbk
!  namelist file extension changed .inp --> .nml
!
!  Revision 1.3  2005/11/15 11:35:02  lars
!  documentation finish for print
!
!  Revision 1.2  2005/07/21 10:20:00  lars
!  polished documentation
!
!  Revision 1.1  2005/06/27 10:54:33  kbk
!  new files needed
!

! !LOCAL VARIABLES:
!
!  proportionality coefficient for
!  parameterizing non-local transport
   REALTYPE                              ::    Cg

!  coefficient from computation of
!  turbulent shear velocity
   REALTYPE                              ::    Vtc

!  acceleration of gravity
   REALTYPE                              ::    g

!  reference density
   REALTYPE                              ::    rho_0

!  g/rho_0
   REALTYPE                              ::    gorho0

!  critical bulk Richardson number
   REALTYPE                              ::    Ric

!  compute surface and bottom BBL
   logical                               ::    kpp_sbl,kpp_bbl

!  compute internal mixing
   logical                               ::    kpp_interior

!  use clipping of MLD at Ekman and Monin-Oboukhov scale
   logical                               ::    clip_mld

!  positions of grid faces and centers
   REALTYPE, dimension(:), allocatable   ::    z_w,z_r

!  distance between centers
   REALTYPE, dimension(:), allocatable   ::    h_r

   integer                               ::    ksblOld
   REALTYPE                              ::    zsblOld

!
!
!EOP
!-----------------------------------------------------------------------

   contains

!-----------------------------------------------------------------------
!BOP
!
! !IROUTINE: Initialise the KPP module
!
! !INTERFACE:
   subroutine init_kpp(namlst,fn,nlev,h0,h,kpp_g,kpp_rho_0)
!
! !DESCRIPTION:
! This routine first reads the namelist {\tt kpp}, which has to be contained
! in a file with filename specified by the string {\tt fn} (typically called
! {\tt kpp.nml}). Since the {\tt kpp} module uses fields defined in the
! {\tt turbulence} module, it has to allocate dynamic memory for them.
! Apart from this, this routine reports the model settings and initialises a
! number of parameters needed later in the time loop.
!
! If you call the GOTM KPP routines from a three-dimensional model, make sure
! that this function is called \emph{after} the call to {\tt init\_turbulence()}.
! Also make sure that you pass the correct arguments.
!
! !USES:
   IMPLICIT NONE
!
! !INPUT PARAMETERS:

!  namelist reference
   integer,          intent(in)        :: namlst

!  filename containing namelist
   character(len=*), intent(in)        :: fn

!  number of grid cells
   integer,          intent(in)        :: nlev

!  bathymetry (m)
   REALTYPE,         intent(in)        :: h0

!  size of grid cells (m)
   REALTYPE,         intent(in)        :: h(0:nlev)

!  acceleration of gravity (m/s^2)
   REALTYPE,         intent(in)        :: kpp_g

!  reference density (kg/m^3)
   REALTYPE,         intent(in)        :: kpp_rho_0

! !REVISION HISTORY:
!  Original author(s): Lars Umlauf
!
!EOP
!
! !LOCAL VARIABLES:
   integer                             :: k
   integer                             :: rc

   namelist /kpp/                      kpp_sbl,kpp_bbl,kpp_interior,    &
                                       clip_mld,Ric
!
!-----------------------------------------------------------------------
!BOC

   LEVEL1 'init_kpp...'

   ! read the variables from the namelist file
   open(namlst,file=fn,status='old',action='read',err=80)

   LEVEL2 'reading kpp namelist...'

   read(namlst,nml=kpp,err=81)
   close (namlst)

   LEVEL2 'done.'

!  allocate memory for variables defined in other modules
!
   allocate(num(0:nlev),stat=rc)
   if (rc /= 0) stop 'init_turbulence: Error allocating (num)'
   num = _ZERO_

   allocate(nuh(0:nlev),stat=rc)
   if (rc /= 0) stop 'init_turbulence: Error allocating (nuh)'
   nuh = _ZERO_

   allocate(nus(0:nlev),stat=rc)
   if (rc /= 0) stop 'init_turbulence: Error allocating (nus)'
   nus = _ZERO_

   allocate(gamu(0:nlev),stat=rc)
   if (rc /= 0) stop 'init_turbulence: Error allocating (gamu)'
   gamu = _ZERO_

   allocate(gamv(0:nlev),stat=rc)
   if (rc /= 0) stop 'init_turbulence: Error allocating (gamv)'
   gamv = _ZERO_

   allocate(gamh(0:nlev),stat=rc)
   if (rc /= 0) stop 'init_turbulence: Error allocating (gamh)'
   gamh = _ZERO_

   allocate(gams(0:nlev),stat=rc)
   if (rc /= 0) stop 'init_turbulence: Error allocating (gams)'
   gams = _ZERO_

   allocate(Rig(0:nlev),stat=rc)
   if (rc /= 0) stop 'init_turbulence: Error allocating (Rig)'
   Rig = _ZERO_

   allocate(z_w(0:nlev),stat=rc)
   if (rc /= 0) stop 'init_turbulence: Error allocating (z_w)'
   z_w = _ZERO_

   allocate(z_r(0:nlev),stat=rc)
   if (rc /= 0) stop 'init_turbulence: Error allocating (z_r)'
   z_r = _ZERO_

   allocate(h_r(0:nlev),stat=rc)
   if (rc /= 0) stop 'init_turbulence: Error allocating (h_r)'
   h_r = _ZERO_

# ifdef EXTRA_OUTPUT

   allocate(turb1(0:nlev),stat=rc)
   if (rc /= 0) stop 'init_turbulence: Error allocating (turb1)'
   turb1 = _ZERO_

   allocate(turb2(0:nlev),stat=rc)
   if (rc /= 0) stop 'init_turbulence: Error allocating (turb2)'
   turb2 = _ZERO_

   allocate(turb3(0:nlev),stat=rc)
   if (rc /= 0) stop 'init_turbulence: Error allocating (turb3)'
   turb3 = _ZERO_

   allocate(turb4(0:nlev),stat=rc)
   if (rc /= 0) stop 'init_turbulence: Error allocating (turb4)'
   turb4 = _ZERO_

   allocate(turb5(0:nlev),stat=rc)
   if (rc /= 0) stop 'init_turbulence: Error allocating (turb5)'
   turb5 = _ZERO_

# endif


!  report model parameters

   LEVEL2 '--------------------------------------------------------'
   LEVEL3 'You are using the KPP turbulence model          '
   LEVEL3 'with the following specifications:                      '
   LEVEL3 '                                                        '
   if (kpp_interior) then
      LEVEL3 'Interior mixing algorithm                  - active -   '
# ifdef KPP_SHEAR
      LEVEL4 'KPP shear instability mixing           - active -   '
# else
      LEVEL4 'KPP shear instability mixing       - not active -   '
# endif
# ifdef KPP_INTERNAL_WAVE
      LEVEL4 'KPP internal wave mixing               - active -   '
# else
      LEVEL4 'KPP internal wave mixing           - not active -   '
# endif
# ifdef KPP_CONVEC
      LEVEL4 'KPP convective instability mixing      - active -   '
# else
      LEVEL4 'KPP convective instability mixing  - not active -   '
# endif
# ifdef KPP_DDMIX
      LEVEL4 'KPP double diffusive mixing            - active -   '
# else
      LEVEL4 'KPP double diffusive mixing        - not active -   '
# endif
   else
      LEVEL3 'Interior mixing algorithm              - not active -   '
   endif

   if (kpp_sbl) then
      LEVEL3 'Surface layer mixing algorithm             - active -   '
      if (clip_mld) then
         LEVEL4 'Clipping at Ekman/Oboukhov scale       - active -   '
      else
         LEVEL4 'Clipping at Ekman/Oboukhov scale   - not active -   '
      endif
# ifdef KPP_SALINITY
      LEVEL4 'Compute salinity fluxes                - active -   '
# else
      LEVEL4 'Compute salinity fluxes            - not active -   '
# endif
# ifdef NONLOCAL
      LEVEL4 'Nonlocal fluxes                        - active -   '
# else
      LEVEL4 'Nonlocal fluxes                    - not active -   '
# endif
# ifdef KPP_TWOPOINT_REF
      LEVEL4 'Ri_b from 2-point interpolation        - active -   '
# else
      LEVEL4 'Ri_b from 2-point interpolation    - not active -   '
# endif
# ifdef KPP_IP_FC
      LEVEL4 'F_c =0 criterion for SL-depth          - active -   '
# else
      LEVEL4 'Ri_b - Ri_c =0 criterion for SL-depth  - active -   '
# endif
# ifdef KPP_CLIP_GS
      LEVEL4 'Clipping G''(sigma) for matching        - active -   '
# else
      LEVEL4 'Clipping G''(sigma) for matching    - not active -   '
# endif

      LEVEL4 'Ri_c  = ', Ric

   else
      LEVEL3 'Surface layer mixing algorithm         - not active -   '
   endif

   if (kpp_bbl) then
      LEVEL3 'Bottom layer mixing algorithm              - active -   '
      LEVEL4 '(Same parameters as surface layer mixing)'

   else
      LEVEL3 'Bottom layer mixing algorithm          - not active -   '
   endif

   LEVEL2 '--------------------------------------------------------'



!  pre-compute coefficient for turbulent shear velocity
   Vtc=Cv*sqrt(-betaT)/(sqrt(cs*epsilon)*Ric*kappa*kappa)


!  pre-compute proportionality coefficient for
!  boundary layer non-local transport
   Cg=Cstar*kappa*(cs*kappa*epsilon)**(1.0/3.0)


!  set acceleration of gravity and reference density
   g        = kpp_g
   rho_0    = kpp_rho_0
   gorho0   = g/rho_0

   LEVEL1 'done.'

   return

80 FATAL 'I could not open "kpp.nml"'
   stop 'init_kpp'
81 FATAL 'I could not read "kpp" namelist'
   stop 'init_kpp'

 end subroutine init_kpp
!EOC

!-----------------------------------------------------------------------
!BOP
!
! !IROUTINE: Loop over the KPP-algorithm
!
! !INTERFACE:
   subroutine do_kpp(nlev,h0,h,rho,u,v,NN,NNT,NNS,SS,u_taus,u_taub,  &
                     tFlux,btFlux,sFlux,bsFlux,tRad,bRad,f)
!
! !DESCRIPTION:
! Here, the time step for the KPP model is managed. If {\tt kpp\_interior=.true.}
! in {\tt kpp.nml}, the mixing algorithm for the computation of the interior
! diffusivities is called first. This algorithm is described in \sect{sec:kppInterior}.
! Then, if {\tt kpp\_sbl=.true.} and {\tt kpp\_bbl=.true.}, the algorithms
! for the surface and bottom boundary layer are called. They are described in
! \sect{sec:kppSurface} and \sect{sec:kppBottom}, respectively.
!
! If this routine is called from a three-dimensional code, it is essential to
! pass the correct arguments. The first 3 parameters relate to the numerical
! grid, discussed in \sect{SectionNumericsMean}. Note that {\tt h0} denotes
! the local bathymetry, i.e.\ the positive distance between the reference level
! $z=0$ and the bottom.

! The next three parameters denote the potential density, $\rho$,
! and the two mean velocity components, $U$ and $V$. The buoyancy frequency, $N^2$,
! and the different contributions to it, $N_\Theta^2$ and $N_S^2$, have to be computed
! from the potential density as discussed in \sect{sec:stratification}. The shear frequency,
! $M^2$, is defined in \eq{MSquared}. The vertical discretisation does not necessarly
! have to follow \eq{shearsquared}, since in the KPP model no TKE equation is solved
! and thus energy conservation is not an issue. All three-dimensional fields have to
! be interpolated "in a smart way" to the water column defined by GOTM. The corresponding
! interpolation schemes may be quite different for the different staggered grids,
! finite volume, and finite element approaches used in the horizontal. Therefore,
! we cannot offer a general recipe here.

! The bottom friction velocity is computed as described in \sect{sec:friction}. If this
! parameter is passed from a three-dimensional code, it has to be insured that the parameter
! $r$ in \eq{uStar} is computed consistently, see \eq{rParam}.
!
! All fluxes without exception are counted positive, if they enter the water body.
! Note that for consistency, the equations of state in GOTM cannot
! be used if the KPP routines are called from a 3-D model. Therefore,
! it is necessary to pass the temperature and salinity fluxes, as well as the
! corresponding buoyancy fluxes. The same applies
! to the radiative fluxes. The user is responsible for
! performing the flux conversions in the correct way. To get an idea
! have a look at \sect{sec:convertFluxes}.
!
! The last argument is the Coriolis parameter, $f$. It is only used for clippling the mixing
! depth at the Ekman depth.
!
!
!
! !USES:
   IMPLICIT NONE
!
! !INPUT PARAMETERS:
!  number of grid cells
   integer                                       :: nlev

!  bathymetry (m)
   REALTYPE                                      :: h0

!  thickness of grid cells (m)
   REALTYPE                                      :: h(0:nlev)

!  potential density at grid centers (kg/m^3)
   REALTYPE                                      :: rho(0:nlev)

!  velocity components at grid centers (m/s)
   REALTYPE                                      :: u(0:nlev),v(0:nlev)

!  square of buoyancy frequency (1/s^2)
   REALTYPE                                      :: NN(0:nlev)

!  square of buoyancy frequency caused by
!  temperature and salinity stratification
   REALTYPE                                      :: NNT(0:nlev),NNS(0:nlev)

!  square of shear frequency (1/s^2)
   REALTYPE                                      :: SS(0:nlev)

!  surface and bottom friction velocities (m/s)
   REALTYPE                                      :: u_taus,u_taub

!  surface temperature flux (K m/s) and
!  salinity flux (psu m/s) (negative for loss)
   REALTYPE                                      :: tFlux,sFlux

!  surface buoyancy fluxes (m^2/s^3) due to
!  heat and salinity fluxes
   REALTYPE                                      :: btFlux,bsFlux

!  radiative flux [ I(z)/(rho Cp) ] (K m/s)
!  and associated buoyancy flux (m^2/s^3)
   REALTYPE                                      :: tRad(0:nlev),bRad(0:nlev)

!  Coriolis parameter (rad/s)
   REALTYPE                                      :: f


! !REVISION HISTORY:
!  Original author(s): Lars Umlauf
!
!EOP
!
! !LOCAL VARIABLES:
   REALTYPE                            :: cff
   integer                             :: k

!
!-----------------------------------------------------------------------
!BOC
!
!-----------------------------------------------------------------------
! Update model grid
!-----------------------------------------------------------------------

!  Compute distance between centers (between rho-points)
!  Note that h is the distance between faces (between w-points)
765
766
   do k=1,nlev-1
      h_r(k) = 0.5*(h(k) + h(k+1))
dumoda01's avatar
dumoda01 committed
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
   enddo

!  Compute position of interfaces (w-points)
   z_w(0) = - h0
   do k=1,nlev
      z_w(k) = z_w(k-1) + h(k)
   enddo

!  Compute position of centers (rho-points)
   z_r(1) = - h0 + 0.5*h(1)
   do k=2,nlev
      z_r(k) = z_r(k-1) + h_r(k-1)
   enddo



!-----------------------------------------------------------------------
! compute interior mixing
!-----------------------------------------------------------------------

   if (kpp_interior) then
      call interior(nlev,NN,NNT,NNS,SS)
   else
      num = _ZERO_
      nuh = _ZERO_
      nus = _ZERO_
   endif


!-----------------------------------------------------------------------
! compute surface boundary layer mixing
!-----------------------------------------------------------------------

   if (kpp_sbl) then
      call surface_layer(nlev,h0,h,rho,u,v,NN,u_taus,u_taub,            &
                         tFlux,btFlux,sFlux,bsFlux,tRad,bRad,f)
   endif

!-----------------------------------------------------------------------
! compute bottom boundary layer mixing
!-----------------------------------------------------------------------

   if (kpp_bbl) then
      call bottom_layer(nlev,h0,h,rho,u,v,NN,u_taus,u_taub,             &
                        _ZERO_,_ZERO_,_ZERO_,_ZERO_,tRad,bRad,f)
   endif



 end subroutine do_kpp
!EOC



!-----------------------------------------------------------------------
!BOP
!
! !IROUTINE: Compute interior fluxes \label{sec:kppInterior}
!
! !INTERFACE:
   subroutine interior(nlev,NN,NNT,NNS,SS)
!
! !DESCRIPTION:
! Here, the interior diffusivities (defined as the diffusivities outside the
! surface and bottom boundary layers) are computed. The algorithms are identical
! to those suggested by \cite{Largeetal94}. For numerical efficiency, the
! algorithms for different physical processes are active only if certain
! pre-processor macros are defined in {\tt cppdefs.h}.
! \begin{itemize}
!  \item The shear instability algorithm is active if the macro
!        {\tt KPP\_SHEAR} is defined.
!  \item The internal wave algorithm is active if the macro
!        {\tt KPP\_INTERNAL\_WAVE} is defined.
!  \item The convective instability algorithm is active if the macro
!        {\tt KPP\_CONVEC} is defined.
!  \item The double-diffusion algorithm is active if the macro
!        {\tt KPP\_DDMIX} is defined. Note that in this case, the
!        macro {\tt SALINITY} has to be defined as well.
! \end{itemize}

!
! !USES:
   IMPLICIT NONE
!
! !INPUT PARAMETERS:

!  number of grid cells
   integer                                       :: nlev

!  square of buoyancy frequency (1/s^2)
   REALTYPE                                      :: NN(0:nlev)

!  square of buoyancy frequencies caused by
!  temperature and salinity stratification
   REALTYPE                                      :: NNT(0:nlev),NNS(0:nlev)

!  square of shear frequency (1/s^2)
   REALTYPE                                      :: SS(0:nlev)


! !REVISION HISTORY:
!  Original author(s): Lars Umlauf
!
!EOP
!
! !LOCAL VARIABLES:
   REALTYPE , parameter       :: eps=1.0E-14

   integer                    :: i
   REALTYPE                   :: cff,shear2
   REALTYPE                   :: nu_sx,nu_sxc
   REALTYPE                   :: iwm,iws
   REALTYPE                   :: drhoT,drhoS,Rrho,nu_dds,nu_ddt

!
!-----------------------------------------------------------------------
!BOC
!
!-----------------------------------------------------------------------
! Compute gradient Richardson number
!-----------------------------------------------------------------------
!
   do i=1,nlev-1
      Rig(i) = NN(i)/(SS(i) + eps)
   enddo

   Rig(0)    = _ZERO_
   Rig(nlev) = _ZERO_
!
!-----------------------------------------------------------------------
!  Compute "interior" viscosities and diffusivities everywhere as
!  the superposition of three processes: local Richardson number
!  instability due to resolved vertical shear, internal wave
!  breaking, and double diffusion.
!-----------------------------------------------------------------------
!
   do i=1,nlev-1
!
!     Smooth gradient Richardson number
      Rig(i)=0.25*Rig(i-1) + 0.50*Rig(i) + 0.25*Rig(i+1)
!
!     Compute interior diffusivity due to shear instability mixing.
# ifdef KPP_SHEAR
      cff=min(_ONE_,max(_ZERO_,Rig(i))/Ri0)
      nu_sx  = _ONE_-cff*cff
      nu_sx  = nu_sx*nu_sx*nu_sx
!
!     The shear mixing should be also a function of the actual magnitude
!     of the shear, see Polzin (1996, JPO, 1409-1425).
      shear2 = SS(i)
      cff    = shear2*shear2/(shear2*shear2+16.0E-10)
      nu_sx  = cff*nu_sx
# else
      nu_sx=_ZERO_
# endif
!
# ifdef KPP_INTERNAL_WAVE
!
!      Compute interior diffusivity due to wave breaking
!
!      Version A, see Gargett and Holloway (1984)
!      cff  =  _ONE_/sqrt(max(NN(i),1.0d-7))
!      iwm  =  1.0E-6*cff
!      iws  =  1.0E-7*cff

!     Version B, see Large et al. (1994)
      iwm  =  nuwm
      iws  =  nuws
# else
      iwm  =  _ZERO_
      iws  =  _ZERO_
# endif
!
# ifdef KPP_CONVEC
!     Compute interior convective diffusivity due to static instability
!     mixing
      cff    =  max(NN(i),bvfcon)
      cff    =  min(_ONE_,(bvfcon-cff)/bvfcon)
      nu_sxc =  _ONE_-cff*cff
      nu_sxc =  nu_sxc*nu_sxc*nu_sxc
# else
      nu_sxc =  _ZERO_
# endif
!
!     Sum contributions due to internal wave breaking, shear instability
!     and convective diffusivity due to shear instability.
      num(i)=iwm+nu0m*nu_sx+nu0c*nu_sxc
      nuh(i)=iws+nu0s*nu_sx+nu0c*nu_sxc
      nus(i)=nuh(i)
!
# ifdef KPP_DDMIX
!
!-----------------------------------------------------------------------
!  Compute double-diffusive mixing.  It can occur when vertical
!  gradient of density is stable but the vertical gradient of
!  salinity (salt figering) or temperature (diffusive convection)
!  is unstable.
!-----------------------------------------------------------------------
!
!     Compute double-diffusive density ratio, Rrho.
      drhoT = NNT(i)
      drhoS = NNS(i)
      Rrho= - drhoT/drhoS
!
!
!     Salt fingering case.
      if ((Rrho.gt._ONE_).and.(drhoS.gt._ZERO_)) then
!
!        Compute interior diffusivity for double diffusive mixing of
!        salinity.  Upper bound "Rrho" by "Rrho0"; (Rrho0=1.9, nuf=0.001).
         Rrho=min(Rrho,Rrho0)
         nu_dds=_ONE_-((Rrho-_ONE_)/(Rrho0-_ONE_))**2.0
         nu_dds=nuf*nu_dds*nu_dds*nu_dds
!
!        Compute interior diffusivity for double diffusive mixing
!        of temperature (fdd=0.7).
         nu_ddt=fdd*nu_dds
!
!
!     Diffusive convection case.
      elseif ((_ZERO_.lt.Rrho).and.(Rrho.lt._ONE_).and.(drhoS.lt._ZERO_)) then
!
!        Compute interior diffusivity for double diffusive mixing of
!        temperature (Marmorino and Caldwell, 1976); (nu=1.5e-6,
!        tdd1=0.909, tdd2=4.6, tdd3=0.54).
         nu_ddt=nu*tdd1*exp(tdd2*exp(-tdd3*((_ONE_/Rrho)-_ONE_)))

!        Compute interior diffusivity for double diffusive mixing
!        of salinity (sdd1=0.15, sdd2=1.85, sdd3=0.85).
!
         if (Rrho.lt.0.5) then
            nu_dds=nu_ddt*sdd1*Rrho
         else
            nu_dds=nu_ddt*(sdd2*Rrho-sdd3)
For faster browsing, not all history is shown. View entire blame