/[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.3 - (show annotations) (download)
Wed Sep 25 19:36:50 2002 UTC (21 years, 7 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint46n_post, checkpoint51k_post, checkpoint47e_post, checkpoint52l_pre, hrcube4, hrcube5, checkpoint46l_post, checkpoint47c_post, checkpoint50c_post, checkpoint52d_pre, checkpoint48e_post, checkpoint50c_pre, checkpoint52j_pre, checkpoint51o_pre, checkpoint54d_post, checkpoint54e_post, checkpoint51l_post, checkpoint48i_post, checkpoint46l_pre, checkpoint50d_pre, checkpoint52k_post, checkpoint55, checkpoint54, checkpoint51, checkpoint50, checkpoint53, checkpoint52, checkpoint50d_post, checkpoint52f_post, checkpoint50b_pre, checkpoint54f_post, checkpoint51f_post, checkpoint48b_post, checkpoint51d_post, checkpoint48c_pre, checkpoint47d_pre, checkpoint51t_post, checkpoint51n_post, checkpoint52i_pre, hrcube_1, hrcube_2, hrcube_3, checkpoint51s_post, checkpoint47a_post, checkpoint55c_post, checkpoint48d_pre, checkpoint51j_post, checkpoint47i_post, checkpoint52e_pre, checkpoint52e_post, checkpoint51n_pre, checkpoint47d_post, checkpoint53d_post, checkpoint48d_post, checkpoint48f_post, checkpoint52b_pre, checkpoint54b_post, checkpoint46j_pre, checkpoint51l_pre, checkpoint52m_post, checkpoint48h_post, checkpoint51q_post, checkpoint51b_pre, checkpoint47g_post, checkpoint52b_post, checkpoint52c_post, checkpoint46j_post, checkpoint51h_pre, checkpoint46k_post, checkpoint48a_post, checkpoint50f_post, checkpoint50a_post, checkpoint50f_pre, checkpoint52f_pre, checkpoint47j_post, checkpoint54a_pre, checkpoint53c_post, branch-exfmods-tag, branchpoint-genmake2, checkpoint54a_post, checkpoint51r_post, checkpoint48c_post, checkpoint51i_post, checkpoint55b_post, checkpoint51b_post, checkpoint51c_post, checkpoint53a_post, checkpoint47b_post, checkpoint52d_post, checkpoint53g_post, checkpoint46m_post, checkpoint50g_post, checkpoint52a_pre, checkpoint50h_post, checkpoint52i_post, checkpoint50e_pre, checkpoint50i_post, checkpoint51i_pre, checkpoint52h_pre, checkpoint53f_post, checkpoint52j_post, checkpoint47f_post, checkpoint50e_post, checkpoint46i_post, branch-netcdf, checkpoint52l_post, checkpoint52n_post, checkpoint53b_pre, checkpoint51e_post, checkpoint55a_post, checkpoint47, checkpoint48, checkpoint49, checkpoint46h_post, checkpoint51o_post, checkpoint51f_pre, checkpoint48g_post, checkpoint53b_post, checkpoint47h_post, checkpoint52a_post, checkpoint51g_post, ecco_c52_e35, checkpoint50b_post, checkpoint51m_post, checkpoint53d_pre, checkpoint54c_post, checkpoint51a_post, checkpoint51p_post, checkpoint51u_post
Branch point for: branch-exfmods-curt, branch-genmake2, branch-nonh, tg2-branch, netcdf-sm0, checkpoint51n_branch
Changes since 1.2: +2 -2 lines
o cleaned up the use of rhoNil and rhoConst.
  - rhoNil should only appear in the LINEAR equation of state, everywhere
    else rhoNil is replaced by rhoConst, e.g. find_rho computes rho-rhoConst
    and the dynamical equations are all divided by rhoConst
o introduced new parameter rhoConstFresh, a reference density of fresh
  water, to remove the fresh water flux's dependence on rhoNil. The default
  value is 999.8 kg/m^3
o cleanup up external_forcing.F and external_forcing_surf.F
  - can now be used by both OCEANIC and OCEANICP

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*rhoConst]
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*rhoConst/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