/[MITgcm]/MITgcm/pkg/layers/layers_calc.F
ViewVC logotype

Annotation of /MITgcm/pkg/layers/layers_calc.F

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


Revision 1.16 - (hide annotations) (download)
Thu Oct 18 12:59:17 2012 UTC (11 years, 8 months ago) by jmc
Branch: MAIN
Changes since 1.15: +79 -86 lines
- remove ALLOW_LAYERS_OUTPUT
- finish to rename variables

1 jmc 1.16 C $Header: /u/gcmpack/MITgcm/pkg/layers/layers_calc.F,v 1.15 2012/10/17 18:49:15 rpa Exp $
2 rpa 1.1 C $Name: $
3    
4     #include "LAYERS_OPTIONS.h"
5 jmc 1.6 #ifdef ALLOW_GMREDI
6 dfer 1.5 #include "GMREDI_OPTIONS.h"
7 jmc 1.6 #endif
8 rpa 1.1
9 jmc 1.4 CBOP 0
10     C !ROUTINE: LAYERS_CALC
11 rpa 1.1
12 jmc 1.4 C !INTERFACE:
13 rpa 1.1 SUBROUTINE LAYERS_CALC(
14 jmc 1.4 I myTime, myIter, myThid )
15 rpa 1.1
16 jmc 1.4 C !DESCRIPTION:
17 rpa 1.1 C ===================================================================
18     C Calculate the transport in isopycnal layers.
19 gforget 1.13 C This was the meat of the LAYERS package, which
20     C has been moved to S/R LAYERS_FLUXCALC.F
21 rpa 1.1 C ===================================================================
22    
23 jmc 1.4 C !USES:
24 rpa 1.1 IMPLICIT NONE
25     #include "SIZE.h"
26 jmc 1.4 #include "EEPARAMS.h"
27     #include "PARAMS.h"
28 rpa 1.1 #include "GRID.h"
29     #include "DYNVARS.h"
30     #include "LAYERS_SIZE.h"
31     #include "LAYERS.h"
32 jmc 1.6 #ifdef ALLOW_GMREDI
33     # include "GMREDI.h"
34 dfer 1.5 #endif
35 rpa 1.1
36 jmc 1.4 C !INPUT PARAMETERS:
37     C myTime :: Current time in simulation
38     C myIter :: Current iteration number
39     C myThid :: my Thread Id number
40     _RL myTime
41     INTEGER myIter
42     INTEGER myThid
43     CEOP
44 rpa 1.1
45     #ifdef ALLOW_LAYERS
46    
47 jmc 1.4 C !LOCAL VARIABLES:
48     C bi, bj :: tile indices
49 rpa 1.1 C i,j :: horizontal indices
50 gforget 1.14 C iLa :: layer coordinate index
51 rpa 1.1 C k :: vertical index for model grid
52    
53 gforget 1.14 INTEGER bi, bj, iLa
54 gforget 1.13 #ifdef LAYERS_PRHO_REF
55     INTEGER i, j, k
56 jmc 1.6 #endif
57 rpa 1.1
58 jmc 1.4 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
59    
60 gforget 1.14 DO iLa=1,layers_maxNum
61 jmc 1.16
62     IF ( layers_num(iLa) .EQ. 1 ) THEN
63     CALL LAYERS_FLUXCALC( uVel,vVel,theta,iLa,
64     & layers_UH(1-OLx,1-OLy,1,1,1,iLa),
65     & layers_VH(1-OLx,1-OLy,1,1,1,iLa),
66     & layers_Hw(1-OLx,1-OLy,1,1,1,iLa),
67     & layers_Hs(1-OLx,1-OLy,1,1,1,iLa),
68     & layers_PIw(1-OLx,1-OLy,1,1,1,iLa),
69     & layers_PIs(1-OLx,1-OLy,1,1,1,iLa),
70     & layers_U(1-OLx,1-OLy,1,1,1,iLa),
71     & layers_V(1-OLx,1-OLy,1,1,1,iLa),
72     & myThid )
73     ELSEIF ( layers_num(iLa) .EQ. 2 ) THEN
74     CALL LAYERS_FLUXCALC( uVel,vVel,salt,iLa,
75     & layers_UH(1-OLx,1-OLy,1,1,1,iLa),
76     & layers_VH(1-OLx,1-OLy,1,1,1,iLa),
77     & layers_Hw(1-OLx,1-OLy,1,1,1,iLa),
78     & layers_Hs(1-OLx,1-OLy,1,1,1,iLa),
79     & layers_PIw(1-OLx,1-OLy,1,1,1,iLa),
80     & layers_PIs(1-OLx,1-OLy,1,1,1,iLa),
81     & layers_U(1-OLx,1-OLy,1,1,1,iLa),
82     & layers_V(1-OLx,1-OLy,1,1,1,iLa),
83     & myThid )
84     ELSEIF ( layers_num(iLa) .EQ. 3 ) THEN
85 dfer 1.10 #ifdef LAYERS_PRHO_REF
86 gforget 1.14 C For layers_num(iLa) = 3, calculate the potential density referenced to
87     C the model level given by layers_krho.
88 jmc 1.16 DO bj=myByLo(myThid),myByHi(myThid)
89     DO bi=myBxLo(myThid),myBxHi(myThid)
90     DO k = 1,Nr
91     CALL FIND_RHO_2D( 1-OLx, sNx+OLx, 1-OLy, sNy+OLy,
92     & layers_krho(iLa),
93     & theta(1-OLx,1-OLy,k,bi,bj),
94     & salt(1-OLx,1-OLy,k,bi,bj),
95     & prho(1-OLx,1-OLy,k,bi,bj,iLa),
96     & k, bi, bj, myThid )
97     DO j = 1-OLy,sNy+OLy
98     DO i = 1-OLx,sNx+OLx
99     prho(i,j,k,bi,bj,iLa) = rhoConst + prho(i,j,k,bi,bj,iLa)
100     ENDDO
101     ENDDO
102     ENDDO
103 jmc 1.11 ENDDO
104     ENDDO
105 jmc 1.16 CALL LAYERS_FLUXCALC( uVel,vVel,
106     & prho(1-OLx,1-OLy,1,1,1,iLa),iLa,
107     & layers_UH(1-OLx,1-OLy,1,1,1,iLa),
108     & layers_VH(1-OLx,1-OLy,1,1,1,iLa),
109     & layers_Hw(1-OLx,1-OLy,1,1,1,iLa),
110     & layers_Hs(1-OLx,1-OLy,1,1,1,iLa),
111     & layers_PIw(1-OLx,1-OLy,1,1,1,iLa),
112     & layers_PIs(1-OLx,1-OLy,1,1,1,iLa),
113     & layers_U(1-OLx,1-OLy,1,1,1,iLa),
114     & layers_V(1-OLx,1-OLy,1,1,1,iLa),
115     & myThid )
116 gforget 1.14 #endif
117 jmc 1.16 ENDIF
118 rpa 1.1
119     #ifdef ALLOW_TIMEAVE
120     C-- Time-average
121 gforget 1.14 cgf layers_maxNum loop and dimension would be needed for
122     cgf the following and tave output to work beyond iLa.EQ.1
123 jmc 1.16 IF ( layers_taveFreq.GT.0. .AND. iLa.EQ.1 ) THEN
124 gforget 1.13 C --- The tile loops
125 jmc 1.16 DO bj=myByLo(myThid),myByHi(myThid)
126     DO bi=myBxLo(myThid),myBxHi(myThid)
127 rpa 1.1
128     #ifdef LAYERS_UFLUX
129 jmc 1.16 CALL TIMEAVE_CUMULATE( layers_UH_T, layers_UH, Nlayers,
130     & deltaTclock, bi, bj, myThid )
131 rpa 1.1 #ifdef LAYERS_THICKNESS
132 jmc 1.16 CALL TIMEAVE_CUMULATE( layers_Hw_T, layers_Hw, Nlayers,
133     & deltaTclock, bi, bj, myThid )
134     CALL TIMEAVE_CUMULATE( layers_PIw_T, layers_PIw, Nlayers,
135     & deltaTclock, bi, bj, myThid )
136     CALL TIMEAVE_CUMULATE( layers_U_T, layers_U, Nlayers,
137     & deltaTclock, bi, bj, myThid )
138 rpa 1.1 #endif /* LAYERS_THICKNESS */
139     #endif /* LAYERS_UFLUX */
140     #ifdef LAYERS_VFLUX
141 jmc 1.16 CALL TIMEAVE_CUMULATE( layers_VH_T, layers_VH, Nlayers,
142     & deltaTclock, bi, bj, myThid )
143 rpa 1.1 #ifdef LAYERS_THICKNESS
144 jmc 1.16 CALL TIMEAVE_CUMULATE( layers_Hs_T, layers_Hs, Nlayers,
145     & deltaTclock, bi, bj, myThid )
146     CALL TIMEAVE_CUMULATE( layers_PIs_T, layers_PIs, Nlayers,
147     & deltaTclock, bi, bj, myThid )
148     CALL TIMEAVE_CUMULATE( layers_V_T, layers_V, Nlayers,
149     & deltaTclock, bi, bj, myThid )
150 rpa 1.1 #endif /* LAYERS_THICKNESS */
151     #endif /* LAYERS_VFLUX */
152    
153 dfer 1.10 #ifdef LAYERS_PRHO_REF
154 jmc 1.16 IF ( layers_num(iLa) .EQ. 3 )
155     & CALL TIMEAVE_CUMULATE( prho_tave, prho, Nr,
156     & deltaTclock, bi, bj, myThid )
157 dfer 1.10 #endif /* LAYERS_PRHO_REF */
158    
159 jmc 1.16 layers_TimeAve(bi,bj)=layers_TimeAve(bi,bj)+deltaTclock
160 rpa 1.1
161     C --- End bi,bj loop
162 jmc 1.16 ENDDO
163     ENDDO
164     ENDIF
165 gforget 1.13 #endif /* ALLOW_TIMEAVE */
166 rpa 1.1
167 jmc 1.16 ENDDO !DO iLa=1,layers_maxNum
168 gforget 1.14
169 rpa 1.1 #endif /* ALLOW_LAYERS */
170    
171     RETURN
172     END

  ViewVC Help
Powered by ViewVC 1.1.22