/[MITgcm]/MITgcm_contrib/heimbach/OpenAD/code_ad_moc/cost_atlantic_heat.F
ViewVC logotype

Contents of /MITgcm_contrib/heimbach/OpenAD/code_ad_moc/cost_atlantic_heat.F

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


Revision 1.1 - (show annotations) (download)
Mon Sep 15 22:16:00 2008 UTC (16 years, 10 months ago) by utke
Branch: MAIN
CVS Tags: HEAD
keep this version around

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, jdel
28 integer jsecbeg , jsecend , isec, idel
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 Meridional:
43 C 80W - 0W at 24N
44 parameter( isecbeg = 278, isecend = 350, jsec = 106 )
45 #if defined (ALLOW_COST_ATLANTIC_HEAT_DOMASS)
46 parameter ( kmaxdepth = 14 )
47 #else
48 parameter ( kmaxdepth = 23 )
49 #endif
50
51 C Zonal:
52 parameter( jsecbeg = 10, jsecend = 27, isec = 59 )
53 C 80W - 0W at 48N
54 C parameter( isecbeg = 70, isecend = 90, jsec = 33 )
55 C parameter ( kmaxdepth = 14 )
56
57 C------------------------------------------------------
58 C Accumulate meridionally integrated transports
59 C Note bar(V)*bar(T) not bar(VT)
60 C Attention pYFaceA [m^2*gravity*rhoConst]
61 C------------------------------------------------------
62
63 #ifdef ALLOW_COST_ATLANTIC_HEAT_DOTHETA
64 jdel = 5
65 #else
66 jdel = 0
67 #endif
68
69 DO bj=myByLo(myThid),myByHi(myThid)
70 DO bi=myBxLo(myThid),myBxHi(myThid)
71
72 locfc = 0.0
73 sum = 0.0
74
75 #define MERID_TRANSPORT
76
77 #ifdef MERID_TRANSPORT
78
79 do j=1,sNy
80 jg = myYGlobalLo-1+(bj-1)*sNy+j
81 if ( jg.GE.(jsec-jdel) .AND.
82 & jg.LE.(jsec+jdel) ) then
83
84 do k = 1, Nr
85 vVel_bar(k) = 0.0
86 theta_bar(k) = 0.0
87 thetaVvel_bar(k) = 0.0
88 countV(k) = 0.0
89 countT(k) = 0.0
90 countTV(k) = 0.0
91 do i=1,sNx
92 ig = myXGlobalLo-1+(bi-1)*sNx+i
93 if ((ig .ge. isecbeg) .and. (ig .le. isecend)) then
94 c--
95 c-- accumulate volume transport
96 vVel_bar(k) = vVel_bar(k)
97 & + cMeanVVel(i,j,k,bi,bj)
98 & *maskS(i,j,k,bi,bj)*dxG(i,j,bi,bj)
99 c-- accumulate heat content
100 theta_bar(k) = theta_bar(k)
101 & + cMeanTheta(i,j,k,bi,bj)
102 & *maskC(i,j,k,bi,bj)*rA(i,j,bi,bj)
103 c-- accumulate heat transport
104 thetaVvel_bar(k) = thetaVvel_bar(k)
105 & + cMeanThetaVVel(i,j,k,bi,bj)*dxG(i,j,bi,bj)
106 & *maskS(i,j,k,bi,bj)*maskC(i,j,k,bi,bj)
107 c-- fill counters
108 countV(k) = countV(k) + maskS(i,j,k,bi,bj)
109 countT(k) = countT(k) + maskC(i,j,k,bi,bj)
110 countTV(k) = countTV(k) +
111 & maskS(i,j,k,bi,bj)*maskC(i,j,k,bi,bj)
112 c--
113 end if
114 enddo
115 enddo
116 c
117 do k = 1, Nr
118 #if defined (ALLOW_COST_ATLANTIC_HEAT_DOMASS)
119 if ( k .LE. kmaxdepth .AND. countV(k) .NE. 0) then
120 sum = sum
121 & + vVel_bar(k)*drF(k)/countV(k)
122 end if
123 #elif defined (ALLOW_COST_ATLANTIC_HEAT_DOTHETA)
124 if ( k .LE. kmaxdepth .AND. countT(k) .NE. 0) then
125 sum = sum
126 & + theta_bar(k)*drF(k)/countT(k)
127 end if
128 #else
129 if ( k .LE. kmaxdepth .AND. countTV(k) .NE. 0) then
130 sum = sum + thetaVVel_bar(k)*drF(k)
131 cph & /countTV(k)
132 end if
133 #endif
134 end do
135
136 #else
137
138 cph need to change this part to go from
139 cph \bar{u}*\bar{T} to \bar{u*T}
140 cph (required store dir. are now in place)
141
142 do i=1,sNx
143 ig = myXGlobalLo-1+(bi-1)*sNx+i
144 if (ig .eq. isec) then
145
146 do k = 1, Nr
147 uVel_bar(k) = 0.0
148 theta_bar(k) = 0.0
149 countT(k) = 0.0
150 countU(k) = 0.0
151 do j=1,sNy
152 jg = myYGlobalLo-1+(bj-1)*sNy+j
153 c
154 if ((jg .ge. jsecbeg) .and. (jg .le. jsecend)) then
155 uVel_bar(k) = uVel_bar(k)
156 & + cMeanUVel(i,j,k,bi,bj)
157 & *maskW(i,j,k,bi,bj)
158 & *maskC(i,j,k,bi,bj)*maskC(i-1,j,k,bi,bj)
159 theta_bar(k) = theta_bar(k) +
160 & 0.5*( cMeanTheta(i,j,k,bi,bj)
161 & +cMeanTheta(i-,j,k,bi,bj) )
162 & *maskW(i,j,k,bi,bj)*dyG(i,j,bi,bj)
163 & *maskC(i,j,k,bi,bj)*maskC(i-1,j,k,bi,bj)
164 countT(k) = countT(k) + maskW(i,j,k,bi,bj)
165 & *maskC(i,j,k,bi,bj)*maskC(i-1,j,k,bi,bj)
166 countU(k) = countU(k) + maskW(i,j,k,bi,bj)
167 & *maskC(i,j,k,bi,bj)*maskC(i-1,j,k,bi,bj)
168 end if
169
170 enddo
171 enddo
172 c
173 do k = 1, Nr
174 if ( k .LE. kmaxdepth .AND.
175 & countT(k) .NE. 0 .AND. countU(k) .NE. 0) then
176 sum = sum
177 & + uVel_bar(k) * theta_bar(k) * drF(k)
178 & / ( countT(k) * countU(k) )
179 end if
180 end do
181
182 #endif
183
184 end if
185 end do
186
187 #if defined (ALLOW_COST_ATLANTIC_HEAT_DOMASS)
188 objf_atl(bi,bj) =
189 & sum*1.E-6
190 #elif defined (ALLOW_COST_ATLANTIC_HEAT_DOTHETA)
191 objf_atl(bi,bj) =
192 & sum*HeatCapacity_Cp*rhoConst
193 #else
194 objf_atl(bi,bj) =
195 & sum*HeatCapacity_Cp*rhoConst/petawatt
196 #endif
197
198 c-- end of bi,bj loop
199 end do
200 end do
201
202 #endif
203
204 end

  ViewVC Help
Powered by ViewVC 1.1.22