/[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.9 - (show annotations) (download)
Tue Oct 14 18:04:53 2008 UTC (15 years, 6 months ago) by utke
Branch: MAIN
CVS Tags: checkpoint62c, checkpoint62a, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62w, checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62z, checkpoint62y, checkpoint62x, checkpoint62, checkpoint63, checkpoint63p, checkpoint63q, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63g, checkpoint63a, checkpoint63b, checkpoint63c, checkpoint62b, checkpoint61f, checkpoint61n, checkpoint61q, checkpoint61e, checkpoint61g, checkpoint61o, checkpoint61l, checkpoint61m, checkpoint61j, checkpoint61k, checkpoint61h, checkpoint61i, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p, checkpoint61z, checkpoint61x, checkpoint61y
Changes since 1.8: +2 -2 lines
changed per the following comment from Chris Hill:
... this specific line of code is wrong, as it will produce annoying differences with different compilers. Ideally it should say
          parameter ( petawatt =  1. _d +15 )

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

  ViewVC Help
Powered by ViewVC 1.1.22