/[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.5 - (show annotations) (download)
Fri Dec 10 20:15:10 2004 UTC (19 years, 6 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57m_post, checkpoint57g_pre, checkpoint57s_post, checkpoint58b_post, checkpoint57b_post, checkpoint57g_post, checkpoint57y_post, checkpoint57r_post, checkpoint57d_post, checkpoint57i_post, checkpoint58, checkpoint58f_post, checkpoint57n_post, checkpoint58d_post, checkpoint58a_post, checkpoint57z_post, checkpoint58t_post, checkpoint58m_post, checkpoint57l_post, checkpoint57t_post, checkpoint57v_post, checkpoint57f_post, checkpoint57a_post, checkpoint57h_pre, checkpoint58w_post, checkpoint57h_post, checkpoint57y_pre, checkpoint58o_post, checkpoint57c_post, checkpoint58p_post, checkpoint58q_post, checkpoint58e_post, checkpoint57c_pre, checkpoint58r_post, checkpoint58n_post, checkpoint57e_post, checkpoint57p_post, checkpint57u_post, checkpoint57q_post, eckpoint57e_pre, checkpoint58k_post, checkpoint58v_post, checkpoint58l_post, checkpoint57h_done, checkpoint57j_post, checkpoint57f_pre, checkpoint58g_post, checkpoint58h_post, checkpoint58j_post, checkpoint57o_post, checkpoint57k_post, checkpoint57w_post, checkpoint58i_post, checkpoint57x_post, checkpoint58c_post, checkpoint58u_post, checkpoint58s_post
Changes since 1.4: +0 -4 lines
This set of changes restores TAMC(!) compatibility.

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 countU(Nr), countV(Nr), countT(Nr)
36 _RL petawatt
37 _RL sum
38 parameter( petawatt = 1.e+15 )
39
40 C 80W - 0W at 24N
41 parameter( isecbeg = 69, isecend = 87, jsec = 28 )
42 cph parameter( isecbeg = 70, isecend = 90, jsec = 30 )
43 parameter( jsecbeg = 10, jsecend = 27, isec = 59 )
44 parameter ( kmaxdepth = 5 )
45 C 80W - 0W at 48N
46 C parameter( isecbeg = 70, isecend = 90, jsec = 33 )
47 C parameter ( kmaxdepth = 14 )
48
49
50
51 C------------------------------------------------------
52 C Accumulate meridionally integrated transports
53 C Note bar(V)*bar(T) not bar(VT)
54 C Attention pYFaceA [m^2*gravity*rhoConst]
55 C------------------------------------------------------
56
57 DO bj=myByLo(myThid),myByHi(myThid)
58 DO bi=myBxLo(myThid),myBxHi(myThid)
59
60 locfc = 0.0
61 sum = 0.0
62
63 #define MERID_TRANSPORT
64
65 #ifdef MERID_TRANSPORT
66
67 #undef ENERGYNORM
68 #ifdef ENERGYNORM
69
70 do i=1,sNx
71 ig = myXGlobalLo-1+(bi-1)*sNx+i
72 if ((ig .ge. isecbeg) .and. (ig .le. isecend)) then
73 sum = 0.0
74 do k = 1, kmaxdepth
75 sum = sum
76 & + (vVel(i,j,k,bi,bj) * maskS(i,j,k,bi,bj)
77 & * drF(k))**2
78 end do
79 locfc = locfc + sum*dxG(i,j,bi,bj)
80 end if
81 end do
82
83 #else
84
85 do j=1,sNy
86 jg = myYGlobalLo-1+(bj-1)*sNy+j
87 if (jg .eq. jsec) then
88
89 do k = 1, Nr
90 vVel_bar(k) = 0.0
91 theta_bar(k) = 0.0
92 countT(k) = 0.0
93 countV(k) = 0.0
94 do i=1,sNx
95 ig = myXGlobalLo-1+(bi-1)*sNx+i
96 c
97 if ((ig .ge. isecbeg) .and. (ig .le. isecend)) then
98 vVel_bar(k) = vVel_bar(k)
99 & + cMeanVVel(i,j,k,bi,bj)
100 & *maskS(i,j,k,bi,bj)
101 & *maskC(i,j,k,bi,bj)*maskC(i,j-1,k,bi,bj)
102 theta_bar(k) = theta_bar(k) +
103 & 0.5*( cMeanTheta(i,j,k,bi,bj)
104 & +cMeanTheta(i,j-1,k,bi,bj) )
105 & *maskS(i,j,k,bi,bj)*dxG(i,j,bi,bj)
106 & *maskC(i,j,k,bi,bj)*maskC(i,j-1,k,bi,bj)
107 countT(k) = countT(k) + maskS(i,j,k,bi,bj)
108 & *maskC(i,j,k,bi,bj)*maskC(i,j-1,k,bi,bj)
109 countV(k) = countV(k) + maskS(i,j,k,bi,bj)
110 & *maskC(i,j,k,bi,bj)*maskC(i,j-1,k,bi,bj)
111 end if
112
113 enddo
114 enddo
115 c
116 do k = 1, Nr
117 if ( k .LE. kmaxdepth .AND.
118 & countT(k) .NE. 0 .AND. countV(k) .NE. 0) then
119 sum = sum
120 & + vVel_bar(k) * theta_bar(k) * drF(k)
121 & / ( countT(k) * countV(k) )
122 end if
123 end do
124
125 #endif /* ENERGYNORM */
126
127 #else
128
129 do i=1,sNx
130 ig = myXGlobalLo-1+(bi-1)*sNx+i
131 if (ig .eq. isec) then
132
133 do k = 1, Nr
134 uVel_bar(k) = 0.0
135 theta_bar(k) = 0.0
136 countT(k) = 0.0
137 countU(k) = 0.0
138 do j=1,sNy
139 jg = myYGlobalLo-1+(bj-1)*sNy+j
140 c
141 if ((jg .ge. jsecbeg) .and. (jg .le. jsecend)) then
142 uVel_bar(k) = uVel_bar(k)
143 & + cMeanUVel(i,j,k,bi,bj)
144 & *maskW(i,j,k,bi,bj)
145 & *maskC(i,j,k,bi,bj)*maskC(i-1,j,k,bi,bj)
146 theta_bar(k) = theta_bar(k) +
147 & 0.5*( cMeanTheta(i,j,k,bi,bj)
148 & +cMeanTheta(i-,j,k,bi,bj) )
149 & *maskW(i,j,k,bi,bj)*dyG(i,j,bi,bj)
150 & *maskC(i,j,k,bi,bj)*maskC(i-1,j,k,bi,bj)
151 countT(k) = countT(k) + maskW(i,j,k,bi,bj)
152 & *maskC(i,j,k,bi,bj)*maskC(i-1,j,k,bi,bj)
153 countU(k) = countU(k) + maskW(i,j,k,bi,bj)
154 & *maskC(i,j,k,bi,bj)*maskC(i-1,j,k,bi,bj)
155 end if
156
157 enddo
158 enddo
159 c
160 do k = 1, Nr
161 if ( k .LE. kmaxdepth .AND.
162 & countT(k) .NE. 0 .AND. countU(k) .NE. 0) then
163 sum = sum
164 & + uVel_bar(k) * theta_bar(k) * drF(k)
165 & / ( countT(k) * countU(k) )
166 end if
167 end do
168
169 #endif
170
171 end if
172 end do
173
174 objf_atl(bi,bj) =
175 & sum*HeatCapacity_Cp*rhoConst/petawatt
176
177 c-- end of bi,bj loop
178 end do
179 end do
180
181 #endif
182
183 end

  ViewVC Help
Powered by ViewVC 1.1.22