1 | MODULE zdfmxl |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE zdfmxl *** |
---|
4 | !! Ocean physics: mixed layer depth |
---|
5 | !!====================================================================== |
---|
6 | !! History : |
---|
7 | !! 9.0 ! 03-08 (G. Madec) OpenMP/autotasking optimization |
---|
8 | !!---------------------------------------------------------------------- |
---|
9 | !! zdf_mxl : Compute the turbocline and mixed layer depths. |
---|
10 | !!---------------------------------------------------------------------- |
---|
11 | !! * Modules used |
---|
12 | USE oce ! ocean dynamics and tracers variables |
---|
13 | USE dom_oce ! ocean space and time domain variables |
---|
14 | USE zdf_oce ! ocean vertical physics |
---|
15 | USE in_out_manager ! I/O manager |
---|
16 | USE prtctl ! Print control |
---|
17 | USE iom |
---|
18 | |
---|
19 | IMPLICIT NONE |
---|
20 | PRIVATE |
---|
21 | |
---|
22 | !! * Routine accessibility |
---|
23 | PUBLIC zdf_mxl ! called by step.F90 |
---|
24 | |
---|
25 | !! * Shared module variables |
---|
26 | INTEGER, PUBLIC, DIMENSION(jpi,jpj) :: & !: |
---|
27 | nmln !: number of level in the mixed layer |
---|
28 | REAL(wp), PUBLIC, DIMENSION(jpi,jpj) :: & !: |
---|
29 | hmld , & !: mixing layer depth (turbocline) (m) |
---|
30 | hmlp , & !: mixed layer depth (rho=rho0+zdcrit) (m) |
---|
31 | hmlpt !: mixed layer depth at t-points (m) |
---|
32 | |
---|
33 | !! * module variables |
---|
34 | REAL(wp) :: & |
---|
35 | avt_c = 5.e-4_wp, & ! Kz criterion for the turbocline depth |
---|
36 | rho_c = 0.01_wp ! density criterion for mixed layer depth |
---|
37 | |
---|
38 | !! * Substitutions |
---|
39 | # include "domzgr_substitute.h90" |
---|
40 | !!---------------------------------------------------------------------- |
---|
41 | !! OPA 9.0 , LOCEAN-IPSL (2005) |
---|
42 | !! $Id$ |
---|
43 | !! This software is governed by the CeCILL licence see modipsl/doc/NEMO_CeCILL.txt |
---|
44 | !!---------------------------------------------------------------------- |
---|
45 | |
---|
46 | CONTAINS |
---|
47 | |
---|
48 | SUBROUTINE zdf_mxl( kt ) |
---|
49 | !!---------------------------------------------------------------------- |
---|
50 | !! *** ROUTINE zdfmxl *** |
---|
51 | !! |
---|
52 | !! ** Purpose : Compute the turbocline depth and the mixed layer depth |
---|
53 | !! with density criteria. |
---|
54 | !! |
---|
55 | !! ** Method : The turbocline depth is the depth at which the vertical |
---|
56 | !! eddy diffusivity coefficient (resulting from the vertical physics |
---|
57 | !! alone, not the isopycnal part, see trazdf.F) fall below a given |
---|
58 | !! value defined locally (avt_c here taken equal to 5 cm/s2) |
---|
59 | !! |
---|
60 | !! ** Action : |
---|
61 | !! |
---|
62 | !! History : |
---|
63 | !! ! 94-11 (M. Imbard) Original code |
---|
64 | !! 8.0 ! 96-01 (E. Guilyardi) sub mixed layer temp. |
---|
65 | !! 8.1 ! 97-07 (G. Madec) optimization |
---|
66 | !! 8.5 ! 02-06 (G. Madec) F90: Free form and module |
---|
67 | !!---------------------------------------------------------------------- |
---|
68 | !! * Arguments |
---|
69 | INTEGER, INTENT( in ) :: kt ! ocean time-step index |
---|
70 | |
---|
71 | !! * Local declarations |
---|
72 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
73 | INTEGER :: ik ! temporary integer |
---|
74 | INTEGER, DIMENSION(jpi,jpj) :: & |
---|
75 | imld ! temporary workspace |
---|
76 | !!---------------------------------------------------------------------- |
---|
77 | |
---|
78 | IF( kt == nit000 ) THEN |
---|
79 | IF(lwp) WRITE(numout,*) |
---|
80 | IF(lwp) WRITE(numout,*) 'zdf_mxl : mixed layer depth' |
---|
81 | IF(lwp) WRITE(numout,*) '~~~~~~~ ' |
---|
82 | ENDIF |
---|
83 | |
---|
84 | |
---|
85 | ! 1. Turbocline depth |
---|
86 | ! ------------------- |
---|
87 | ! last w-level at which avt<avt_c (starting from the bottom jk=jpk) |
---|
88 | ! (since avt(.,.,jpk)=0, we have jpk=< imld =< 2 ) |
---|
89 | DO jk = jpk, 2, -1 |
---|
90 | DO jj = 1, jpj |
---|
91 | DO ji = 1, jpi |
---|
92 | IF( avt(ji,jj,jk) < avt_c ) imld(ji,jj) = jk |
---|
93 | END DO |
---|
94 | END DO |
---|
95 | END DO |
---|
96 | |
---|
97 | ! Turbocline depth and sub-turbocline temperature |
---|
98 | DO jj = 1, jpj |
---|
99 | DO ji = 1, jpi |
---|
100 | ik = imld(ji,jj) |
---|
101 | hmld (ji,jj) = fsdepw(ji,jj,ik) * tmask(ji,jj,1) |
---|
102 | END DO |
---|
103 | END DO |
---|
104 | CALL iom_put( "mldturb", hmld ) ! turbocline depth |
---|
105 | |
---|
106 | !!gm idea |
---|
107 | !! |
---|
108 | !!gm DO jk = jpk, 2, -1 |
---|
109 | !!gm DO jj = 1, jpj |
---|
110 | !!gm DO ji = 1, jpi |
---|
111 | !!gm IF( avt(ji,jj,jk) < avt_c ) hmld(ji,jj) = fsdepw(ji,jj,jk) * tmask(ji,jj,1) |
---|
112 | !!gm END DO |
---|
113 | !!gm END DO |
---|
114 | !!gm END DO |
---|
115 | !!gm |
---|
116 | |
---|
117 | ! 2. Mixed layer depth |
---|
118 | ! -------------------- |
---|
119 | ! Initialization to the number of w ocean point mbathy |
---|
120 | nmln(:,:) = mbathy(:,:) |
---|
121 | |
---|
122 | ! Last w-level at which rhop>=rho surf+rho_c (starting from jpk-1) |
---|
123 | ! (rhop defined at t-point, thus jk-1 for w-level just above) |
---|
124 | DO jk = jpkm1, 2, -1 |
---|
125 | DO jj = 1, jpj |
---|
126 | DO ji = 1, jpi |
---|
127 | IF( rhop(ji,jj,jk) > rhop(ji,jj,1) + rho_c ) nmln(ji,jj) = jk |
---|
128 | END DO |
---|
129 | END DO |
---|
130 | END DO |
---|
131 | |
---|
132 | ! Mixed layer depth |
---|
133 | DO jj = 1, jpj |
---|
134 | DO ji = 1, jpi |
---|
135 | ik = nmln(ji,jj) |
---|
136 | hmlp (ji,jj) = fsdepw(ji,jj,ik) * tmask(ji,jj,1) |
---|
137 | hmlpt(ji,jj) = fsdept(ji,jj,ik-1) |
---|
138 | END DO |
---|
139 | END DO |
---|
140 | CALL iom_put( "mld010", hmlp ) ! mixed layer depth |
---|
141 | |
---|
142 | IF(ln_ctl) CALL prt_ctl( tab2d_1=REAL(nmln,wp), clinfo1=' nmln : ', tab2d_2=hmld, clinfo2=' hmld : ', ovlap=1 ) |
---|
143 | |
---|
144 | END SUBROUTINE zdf_mxl |
---|
145 | |
---|
146 | !!====================================================================== |
---|
147 | END MODULE zdfmxl |
---|