/[MITgcm]/MITgcm/pkg/cost/cost_atlantic_heat.F
ViewVC logotype

Annotation of /MITgcm/pkg/cost/cost_atlantic_heat.F

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


Revision 1.1 - (hide annotations) (download)
Thu Jan 17 17:03:34 2002 UTC (22 years, 4 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint44f_post, chkpt44d_post, checkpoint44e_pre, chkpt44a_post, checkpoint44h_pre, chkpt44c_pre, checkpoint44e_post, checkpoint44g_post, checkpoint44b_post, chkpt44a_pre, checkpoint44b_pre, checkpoint44, chkpt44c_post, checkpoint44f_pre
Branch point for: release1_final
* added Atlantic heat transport cost function
* added vector-valued cost function
...and corresponding options.

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

  ViewVC Help
Powered by ViewVC 1.1.22