/[MITgcm]/MITgcm/pkg/ecco/cost_salt0.F
ViewVC logotype

Contents of /MITgcm/pkg/ecco/cost_salt0.F

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


Revision 1.1 - (show annotations) (download)
Thu Nov 6 22:10:08 2003 UTC (20 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint52l_pre, checkpoint52e_pre, hrcube4, checkpoint52n_post, checkpoint52j_post, checkpoint53d_post, checkpoint54a_pre, checkpoint55c_post, checkpoint54e_post, checkpoint52e_post, checkpoint53c_post, checkpoint55d_pre, checkpoint52j_pre, checkpoint54a_post, branch-netcdf, checkpoint52d_pre, checkpoint52l_post, checkpoint52k_post, checkpoint52b_pre, checkpoint54b_post, checkpoint53b_pre, checkpoint55b_post, checkpoint54d_post, checkpoint53f_post, hrcube_1, checkpoint52m_post, checkpoint55, checkpoint53a_post, checkpoint54, checkpoint54f_post, checkpoint53b_post, checkpoint53, checkpoint52, checkpoint52d_post, checkpoint52a_post, checkpoint52b_post, checkpoint53g_post, checkpoint52f_post, checkpoint52c_post, ecco_c52_e35, hrcube5, checkpoint52a_pre, checkpoint52i_post, checkpoint55a_post, checkpoint53d_pre, checkpoint54c_post, checkpoint52i_pre, checkpoint51u_post, checkpoint52h_pre, checkpoint52f_pre, hrcube_2, hrcube_3
Branch point for: netcdf-sm0
o merging from ecco-branch
o pkg/ecco now containes ecco-specific part of cost function
o top level routines the_main_loop, forward_step
  supersede those in model/src/
  previous input data.cost now in data.ecco
  (new namelist ecco_cost_nml)

1 C $Header: /u/gcmpack/MITgcm/pkg/cost/Attic/cost_salt0.F,v 1.1.2.2 2003/07/16 16:38:58 heimbach Exp $
2
3 #include "COST_CPPOPTIONS.h"
4
5
6 subroutine cost_salt0(
7 I myiter,
8 I mytime,
9 I mythid
10 & )
11
12 c ==================================================================
13 c SUBROUTINE cost_salt0
14 c ==================================================================
15 c
16 c o Calculate the zonal wind stress contribution to the cost function.
17 c
18 c started: Christian Eckert eckert@mit.edu 30-Jun-1999
19 c
20 c changed: Christian Eckert eckert@mit.edu 25-Feb-2000
21 c
22 c - Restructured the code in order to create a package
23 c for the MITgcmUV.
24 c
25 c ==================================================================
26 c SUBROUTINE cost_salt0
27 c ==================================================================
28
29 implicit none
30
31 c == global variables ==
32
33 #include "EEPARAMS.h"
34 #include "SIZE.h"
35 #include "GRID.h"
36
37 #include "ecco_cost.h"
38 #include "ctrl.h"
39 #include "ctrl_dummy.h"
40 #include "optim.h"
41
42 c == routine arguments ==
43
44 integer myiter
45 _RL mytime
46 integer mythid
47
48 c == local variables ==
49
50 integer bi,bj
51 integer i,j,k
52 integer itlo,ithi
53 integer jtlo,jthi
54 integer jmin,jmax
55 integer imin,imax
56 integer nrec
57 integer irec
58 integer ilfld
59
60 _RL fctile
61 _RL fcthread
62 _RL tmpx
63
64 logical doglobalread
65 logical ladinit
66
67 character*(80) fnamefld
68
69 character*(MAX_LEN_MBUF) msgbuf
70
71 c == external functions ==
72
73 integer ilnblnk
74 external ilnblnk
75
76 c == end of interface ==
77
78 jtlo = mybylo(mythid)
79 jthi = mybyhi(mythid)
80 itlo = mybxlo(mythid)
81 ithi = mybxhi(mythid)
82 jmin = 1
83 jmax = sny
84 imin = 1
85 imax = snx
86
87 c-- Read state record from global file.
88 doglobalread = .false.
89 ladinit = .false.
90
91 irec = 1
92
93 #ifdef ALLOW_SALT0_COST_CONTRIBUTION
94
95 if (optimcycle .ge. 0) then
96 ilfld = ilnblnk( xx_salt_file )
97 write(fnamefld(1:80),'(2a,i10.10)')
98 & xx_salt_file(1:ilfld),'.',optimcycle
99 endif
100
101 fcthread = 0. _d 0
102
103 call active_read_xyz_loc( fnamefld, tmpfld3d, irec, doglobalread,
104 & ladinit, optimcycle, mythid
105 & , xx_salt_dummy )
106
107 c-- Loop over this thread's tiles.
108 do bj = jtlo,jthi
109 do bi = itlo,ithi
110
111 c-- Determine the weights to be used.
112
113 fctile = 0. _d 0
114 do k = 1,nr
115 do j = jmin,jmax
116 do i = imin,imax
117 if (_hFacC(i,j,k,bi,bj) .ne. 0.) then
118 tmpx = tmpfld3d(i,j,k,bi,bj)
119 fctile = fctile
120 & + wsalt(k,bi,bj)*cosphi(i,j,bi,bj)
121 & *tmpx*tmpx
122 endif
123 enddo
124 enddo
125 enddo
126
127 objf_salt0(bi,bj) = objf_salt0(bi,bj) + fctile
128 fcthread = fcthread + fctile
129
130 #ifdef ECCO_VERBOSE
131 c-- Print cost function for each tile in each thread.
132 write(msgbuf,'(a)') ' '
133 call print_message( msgbuf, standardmessageunit,
134 & SQUEEZE_RIGHT , mythid)
135 write(msgbuf,'(a,i8.8,1x,i3.3,1x,i3.3)')
136 & ' cost_salt0: irec,bi,bj = ',irec,bi,bj
137 call print_message( msgbuf, standardmessageunit,
138 & SQUEEZE_RIGHT , mythid)
139 write(msgbuf,'(a,d22.15)')
140 & ' cost function (dS(0)) = ',
141 & fctile
142 call print_message( msgbuf, standardmessageunit,
143 & SQUEEZE_RIGHT , mythid)
144 #endif
145 enddo
146 enddo
147
148 #ifdef ECCO_VERBOSE
149 c-- Print cost function for all tiles.
150 _GLOBAL_SUM_R8( fcthread , myThid )
151 write(msgbuf,'(a)') ' '
152 call print_message( msgbuf, standardmessageunit,
153 & SQUEEZE_RIGHT , mythid)
154 write(msgbuf,'(a,i8.8)')
155 & ' cost_zonstress: irec = ',irec
156 call print_message( msgbuf, standardmessageunit,
157 & SQUEEZE_RIGHT , mythid)
158 write(msgbuf,'(a,d22.15)')
159 & ' global cost function value = ',
160 & fcthread
161 call print_message( msgbuf, standardmessageunit,
162 & SQUEEZE_RIGHT , mythid)
163 write(msgbuf,'(a)') ' '
164 call print_message( msgbuf, standardmessageunit,
165 & SQUEEZE_RIGHT , mythid)
166 #endif
167
168 #else
169 c-- Do not enter the calculation of the salinity increment
170 c-- contribution to the final cost function.
171
172 fctile = 0. _d 0
173 fcthread = 0. _d 0
174
175 #ifdef ECCO_VERBOSE
176 _BEGIN_MASTER( mythid )
177 write(msgbuf,'(a)') ' '
178 call print_message( msgbuf, standardmessageunit,
179 & SQUEEZE_RIGHT , mythid)
180 write(msgbuf,'(a,a)')
181 & ' cost_salt0 : no contribution of the I.C. in temp. ',
182 & ' to cost function.'
183 call print_message( msgbuf, standardmessageunit,
184 & SQUEEZE_RIGHT , mythid)
185 write(msgbuf,'(a)') ' '
186 call print_message( msgbuf, standardmessageunit,
187 & SQUEEZE_RIGHT , mythid)
188 _END_MASTER( mythid )
189 #endif
190
191 #endif
192
193 end

  ViewVC Help
Powered by ViewVC 1.1.22