/[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.2 - (show annotations) (download)
Sun Mar 24 02:14:21 2002 UTC (22 years, 2 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint46g_pre, checkpoint46f_post, checkpoint46b_post, checkpoint46d_pre, checkpoint45d_post, checkpoint46a_post, checkpoint46b_pre, checkpoint45a_post, checkpoint46e_pre, checkpoint45b_post, checkpoint46c_pre, checkpoint46, checkpoint46h_pre, checkpoint46a_pre, checkpoint45c_post, checkpoint44h_post, checkpoint46g_post, checkpoint46c_post, checkpoint46e_post, checkpoint45, checkpoint46d_post
Branch point for: release1
Changes since 1.1: +30 -16 lines
Modified test cost functions.

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 kmaxdepth
29 integer i, j, k
30 integer ig, jg
31 integer bi, bj
32 _RL locfc
33 _RL vVel_bar(Nr), theta_bar(Nr)
34 _RL countT(Nr), countV(Nr)
35 _RL petawatt
36 _RL sum
37 parameter( petawatt = 1.e+15 )
38
39 C 80W - 0W at 24N
40 cph parameter( isecbeg = 66, isecend = 86, jsec = 27 )
41 parameter( isecbeg = 70, isecend = 90, jsec = 27 )
42 parameter ( kmaxdepth = 12 )
43 C 80W - 0W at 48N
44 C parameter( isecbeg = 70, isecend = 90, jsec = 33 )
45 C parameter ( kmaxdepth = 14 )
46
47
48
49 C------------------------------------------------------
50 C Accumulate meridionally integrated transports
51 C Note bar(V)*bar(T) not bar(VT)
52 C Attention pYFaceA [m^2*gravity*rhonil]
53 C------------------------------------------------------
54
55 DO bj=myByLo(myThid),myByHi(myThid)
56 DO bi=myBxLo(myThid),myBxHi(myThid)
57
58 locfc = 0.0
59 sum = 0.0
60
61 do j=1,sNy
62 jg = myYGlobalLo-1+(bj-1)*sNy+j
63 if (jg .eq. jsec) then
64
65 #undef ENERGYNORM
66
67 #ifdef ENERGYNORM
68
69 do i=1,sNx
70 ig = myXGlobalLo-1+(bi-1)*sNx+i
71 if ((ig .ge. isecbeg) .and. (ig .le. isecend)) then
72 sum = 0.0
73 do k = 1, kmaxdepth
74 sum = sum
75 & + (vVel(i,j,k,bi,bj) * maskS(i,j,k,bi,bj)
76 & * drF(k))**2
77 end do
78 locfc = locfc + sum*dxG(i,j,bi,bj)
79 end if
80 end do
81
82 #else
83
84 do k = 1, Nr
85 vVel_bar(k) = 0.0
86 theta_bar(k) = 0.0
87 countT(k) = 0.0
88 countV(k) = 0.0
89 do i=1,sNx
90 ig = myXGlobalLo-1+(bi-1)*sNx+i
91 c
92 if ((ig .ge. isecbeg) .and. (ig .le. isecend)) then
93 vVel_bar(k) = vVel_bar(k)
94 & + vVel(i,j,k,bi,bj)
95 & *maskS(i,j,k,bi,bj)
96 theta_bar(k) = theta_bar(k)
97 & + theta(i,j,k,bi,bj)*dxG(i,j,bi,bj)
98 & *maskC(i,j,k,bi,bj)
99 cph & 0.5*(theta(i,j,k,bi,bj)+theta(i,j-1,k,bi,bj) )
100 cph & *maskS(i,j,k,bi,bj)*dxG(i,j,bi,bj)
101 countT(k) = countT(k) + maskC(i,j,k,bi,bj)
102 countV(k) = countV(k) + maskS(i,j,k,bi,bj)
103 end if
104
105 enddo
106 enddo
107 c
108 do k = 1, Nr
109 if ( k .LE. kmaxdepth .AND.
110 & countT(k) .NE. 0 .AND. countV(k) .NE. 0) then
111 sum = sum
112 & + vVel_bar(k) * theta_bar(k) * drF(k)
113 & / ( countT(k) * countV(k) )
114 end if
115 end do
116
117 #endif
118
119 end if
120 end do
121
122 objf_atl(bi,bj) =
123 & sum*HeatCapacity_Cp*rhonil/petawatt
124
125 c-- end of bi,bj loop
126 end do
127 end do
128
129 #endif
130
131 end

  ViewVC Help
Powered by ViewVC 1.1.22