/[MITgcm]/MITgcm/pkg/cost/cost_atlantic_heat.F
ViewVC logotype

Contents of /MITgcm/pkg/cost/cost_atlantic_heat.F

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


Revision 1.6 - (show annotations) (download)
Sat Mar 24 02:23:09 2007 UTC (17 years, 1 month ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint59, checkpoint58y_post, checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59a, checkpoint59c, checkpoint59b, checkpoint59h, checkpoint58x_post
Changes since 1.5: +34 -20 lines
o Changes in cost to compute \bar{T*V} instead of \bar{T}*{V}
o one fix in pkg/flt

1
2 #include "COST_CPPOPTIONS.h"
3
4 subroutine cost_atlantic_heat( myThid )
5 C /==========================================================\
6 C | subroutine cost_atlantic_heat |
7 C | o This routine computes the meridional heat transport. |
8 C | The current indices are for North Atlantic 29N |
9 C | 2x2 global setup. |
10 C \==========================================================/
11 implicit none
12
13 C == Global variables ===
14 #include "SIZE.h"
15 #include "EEPARAMS.h"
16 #include "PARAMS.h"
17 #include "GRID.h"
18 #include "DYNVARS.h"
19 #include "cost.h"
20
21 C ======== Routine arguments ======================
22 C myThid - Thread number for this instance of the routine.
23 integer myThid
24
25 #ifdef ALLOW_COST_ATLANTIC_HEAT
26 C ========= Local variables =========================
27 integer isecbeg , isecend , jsec
28 integer jsecbeg , jsecend , isec
29 integer kmaxdepth
30 integer i, j, k
31 integer ig, jg
32 integer bi, bj
33 _RL locfc
34 _RL uVel_bar(Nr), vVel_bar(Nr), theta_bar(Nr)
35 _RL thetaUvel_bar(Nr), thetaVvel_bar(Nr)
36 _RL countU(Nr), countV(Nr), countT(Nr)
37 _RL countTU(Nr), countTV(Nr)
38 _RL petawatt
39 _RL sum
40 parameter( petawatt = 1.e+15 )
41
42 C 80W - 0W at 24N
43 parameter( isecbeg = 69, isecend = 87, jsec = 28 )
44 cph parameter( isecbeg = 70, isecend = 90, jsec = 30 )
45 parameter( jsecbeg = 10, jsecend = 27, isec = 59 )
46 parameter ( kmaxdepth = 5 )
47 C 80W - 0W at 48N
48 C parameter( isecbeg = 70, isecend = 90, jsec = 33 )
49 C parameter ( kmaxdepth = 14 )
50
51
52
53 C------------------------------------------------------
54 C Accumulate meridionally integrated transports
55 C Note bar(V)*bar(T) not bar(VT)
56 C Attention pYFaceA [m^2*gravity*rhoConst]
57 C------------------------------------------------------
58
59 DO bj=myByLo(myThid),myByHi(myThid)
60 DO bi=myBxLo(myThid),myBxHi(myThid)
61
62 locfc = 0.0
63 sum = 0.0
64
65 #define MERID_TRANSPORT
66
67 #ifdef MERID_TRANSPORT
68
69 #undef ENERGYNORM
70 #ifdef ENERGYNORM
71
72 do i=1,sNx
73 ig = myXGlobalLo-1+(bi-1)*sNx+i
74 if ((ig .ge. isecbeg) .and. (ig .le. isecend)) then
75 sum = 0.0
76 do k = 1, kmaxdepth
77 sum = sum
78 & + (vVel(i,j,k,bi,bj) * maskS(i,j,k,bi,bj)
79 & * drF(k))**2
80 end do
81 locfc = locfc + sum*dxG(i,j,bi,bj)
82 end if
83 end do
84
85 #else
86
87 do j=1,sNy
88 jg = myYGlobalLo-1+(bj-1)*sNy+j
89 if (jg .eq. jsec) then
90
91 do k = 1, Nr
92 vVel_bar(k) = 0.0
93 thetaVvel_bar(k) = 0.0
94 countV(k) = 0.0
95 countTV(k) = 0.0
96 do i=1,sNx
97 ig = myXGlobalLo-1+(bi-1)*sNx+i
98 c
99 if ((ig .ge. isecbeg) .and. (ig .le. isecend)) then
100 vVel_bar(k) = vVel_bar(k)
101 & + cMeanVVel(i,j,k,bi,bj)*maskS(i,j,k,bi,bj)
102
103 thetaVvel_bar(k) = thetaVvel_bar(k)
104 & + cMeanThetaVVel(i,j,k,bi,bj)*dxG(i,j,bi,bj)
105 & *maskS(i,j,k,bi,bj)*maskC(i,j,k,bi,bj)
106
107 countTV(k) = countTV(k) +
108 & maskS(i,j,k,bi,bj)*maskC(i,j,k,bi,bj)
109 countV(k) = countV(k) +
110 & maskS(i,j,k,bi,bj)
111 end if
112
113 enddo
114 enddo
115 c
116 do k = 1, Nr
117 #ifdef ALLOW_COST_ATLANTIC_HEAT_DOMASS
118 if ( k .LE. kmaxdepth .AND. countV(k) .NE. 0) then
119 sum = sum
120 & + vVel_bar(k)*drF(k)/countV(k)
121 end if
122 #else
123 if ( k .LE. kmaxdepth .AND. countTV(k) .NE. 0) then
124 sum = sum
125 & + thetaVVel_bar(k)*drF(k)/countTV(k)
126 end if
127 #endif
128 end do
129
130 #endif /* ENERGYNORM */
131
132 #else
133
134 cph need to change this part to go from
135 cph \bar{u}*\bar{T} to \bar{u*T}
136 cph (required store dir. are now in place)
137
138 do i=1,sNx
139 ig = myXGlobalLo-1+(bi-1)*sNx+i
140 if (ig .eq. isec) then
141
142 do k = 1, Nr
143 uVel_bar(k) = 0.0
144 theta_bar(k) = 0.0
145 countT(k) = 0.0
146 countU(k) = 0.0
147 do j=1,sNy
148 jg = myYGlobalLo-1+(bj-1)*sNy+j
149 c
150 if ((jg .ge. jsecbeg) .and. (jg .le. jsecend)) then
151 uVel_bar(k) = uVel_bar(k)
152 & + cMeanUVel(i,j,k,bi,bj)
153 & *maskW(i,j,k,bi,bj)
154 & *maskC(i,j,k,bi,bj)*maskC(i-1,j,k,bi,bj)
155 theta_bar(k) = theta_bar(k) +
156 & 0.5*( cMeanTheta(i,j,k,bi,bj)
157 & +cMeanTheta(i-,j,k,bi,bj) )
158 & *maskW(i,j,k,bi,bj)*dyG(i,j,bi,bj)
159 & *maskC(i,j,k,bi,bj)*maskC(i-1,j,k,bi,bj)
160 countT(k) = countT(k) + maskW(i,j,k,bi,bj)
161 & *maskC(i,j,k,bi,bj)*maskC(i-1,j,k,bi,bj)
162 countU(k) = countU(k) + maskW(i,j,k,bi,bj)
163 & *maskC(i,j,k,bi,bj)*maskC(i-1,j,k,bi,bj)
164 end if
165
166 enddo
167 enddo
168 c
169 do k = 1, Nr
170 if ( k .LE. kmaxdepth .AND.
171 & countT(k) .NE. 0 .AND. countU(k) .NE. 0) then
172 sum = sum
173 & + uVel_bar(k) * theta_bar(k) * drF(k)
174 & / ( countT(k) * countU(k) )
175 end if
176 end do
177
178 #endif
179
180 end if
181 end do
182
183 #ifdef ALLOW_COST_ATLANTIC_HEAT_DOMASS
184 objf_atl(bi,bj) =
185 & sum*1.E-6
186 #else
187 objf_atl(bi,bj) =
188 & sum*HeatCapacity_Cp*rhoConst/petawatt
189 #endif
190
191 c-- end of bi,bj loop
192 end do
193 end do
194
195 #endif
196
197 end

  ViewVC Help
Powered by ViewVC 1.1.22