| 1 |
C $Header: $ |
| 2 |
C $Name: $ |
| 3 |
|
| 4 |
#include "BLING_OPTIONS.h" |
| 5 |
|
| 6 |
CBOP |
| 7 |
subroutine BLING_MIXEDLAYER( |
| 8 |
U sumMLDepth, |
| 9 |
I bi, bj, imin, imax, jmin, jmax, |
| 10 |
I myIter, myTime, myThid ) |
| 11 |
|
| 12 |
C ================================================================= |
| 13 |
C | subroutine bling_mixedlayer |
| 14 |
C | o |
| 15 |
C | - irr_mix is the same as irr_inst, but averaged throughout |
| 16 |
C | the mixed layer. This quantity is intended to represent the |
| 17 |
C | light to which phytoplankton subject to turbulent transport in |
| 18 |
C | the mixed-layer would be exposed. |
| 19 |
C ================================================================= |
| 20 |
|
| 21 |
implicit none |
| 22 |
|
| 23 |
C === Global variables === |
| 24 |
|
| 25 |
#include "SIZE.h" |
| 26 |
#include "EEPARAMS.h" |
| 27 |
#include "PARAMS.h" |
| 28 |
#include "FFIELDS.h" |
| 29 |
#include "GRID.h" |
| 30 |
#include "DYNVARS.h" |
| 31 |
#include "BLING_VARS.h" |
| 32 |
#include "PTRACERS_SIZE.h" |
| 33 |
#include "PTRACERS_PARAMS.h" |
| 34 |
#ifdef ALLOW_AUTODIFF |
| 35 |
# include "tamc.h" |
| 36 |
#endif |
| 37 |
|
| 38 |
C === Routine arguments === |
| 39 |
C bi,bj :: tile indices |
| 40 |
C iMin,iMax :: computation domain: 1rst index range |
| 41 |
C jMin,jMax :: computation domain: 2nd index range |
| 42 |
C myTime :: current time |
| 43 |
C myIter :: current timestep |
| 44 |
C myThid :: thread Id. number |
| 45 |
INTEGER bi, bj, imin, imax, jmin, jmax |
| 46 |
INTEGER myThid |
| 47 |
INTEGER myIter |
| 48 |
_RL myTime |
| 49 |
C === Output === |
| 50 |
C sumMLDepth :: mixed layer depth |
| 51 |
_RL sumMLDepth(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
| 52 |
|
| 53 |
C === Local variables === |
| 54 |
_RL dens_surf (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
| 55 |
_RL dens_z (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
| 56 |
_RL delta_dens(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
| 57 |
integer i,j,k |
| 58 |
CEOP |
| 59 |
|
| 60 |
c --------------------------------------------------------------------- |
| 61 |
c Mixed layer depth |
| 62 |
|
| 63 |
DO j=jmin,jmax |
| 64 |
DO i=imin,imax |
| 65 |
SumMLDepth(i,j) = drf(1) |
| 66 |
ENDDO |
| 67 |
ENDDO |
| 68 |
|
| 69 |
c Surface density |
| 70 |
CALL FIND_RHO_2D( |
| 71 |
I 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1, |
| 72 |
I theta(1-OLx,1-OLy,1,bi,bj), salt(1-OLx,1-OLy,1,bi,bj), |
| 73 |
O dens_surf, |
| 74 |
I 1, bi, bj, myThid ) |
| 75 |
|
| 76 |
DO k=1,Nr |
| 77 |
DO j=jmin,jmax |
| 78 |
DO i=imin,imax |
| 79 |
if (k.eq.1) then |
| 80 |
delta_dens(i,j,1) = 0. _d 0 |
| 81 |
else |
| 82 |
delta_dens(i,j,k) = 9999. _d 0 |
| 83 |
endif |
| 84 |
ENDDO |
| 85 |
ENDDO |
| 86 |
ENDDO |
| 87 |
|
| 88 |
DO k = 2,Nr |
| 89 |
|
| 90 |
c Potential density |
| 91 |
CALL FIND_RHO_2D( |
| 92 |
I 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1, |
| 93 |
I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj), |
| 94 |
O dens_z, |
| 95 |
I k, bi, bj, myThid ) |
| 96 |
|
| 97 |
DO j=jmin,jmax |
| 98 |
DO i=imin,imax |
| 99 |
|
| 100 |
c SumMLDepth(i,j) = 0. _d 0 |
| 101 |
|
| 102 |
IF (hFacC(i,j,k,bi,bj) .gt. 0. _d 0) THEN |
| 103 |
delta_dens(i,j,k) = dens_z(i,j)-dens_surf(i,j) |
| 104 |
IF (delta_dens(i,j,k) .LT. 0.03 _d 0) THEN |
| 105 |
SumMLDepth(i,j) = SumMLDepth(i,j)+drF(k) |
| 106 |
ENDIF |
| 107 |
ENDIF |
| 108 |
|
| 109 |
ENDDO |
| 110 |
ENDDO |
| 111 |
ENDDO |
| 112 |
|
| 113 |
RETURN |
| 114 |
END |
| 115 |
|