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

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

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


Revision 1.2 - (show annotations) (download)
Mon Oct 11 16:38:53 2004 UTC (19 years, 8 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57d_post, checkpoint57b_post, checkpoint57c_pre, checkpoint55j_post, checkpoint56b_post, checkpoint55h_post, checkpoint57e_post, checkpoint56c_post, checkpoint57a_post, checkpoint55g_post, checkpoint55f_post, checkpoint57a_pre, checkpoint55i_post, checkpoint57, checkpoint56, eckpoint57e_pre, checkpoint57c_post, checkpoint55e_post, checkpoint56a_post, checkpoint55d_post
Changes since 1.1: +3 -1 lines
o ECCO specific cost function terms (up-to-date with 1x1 runs)
o ecco_cost_weights is modified to 1x1 runs
o modifs to allow observations to be read in as
  single file or yearly files

1 C $Header: /u/gcmpack/MITgcm/pkg/cost/Attic/cost_theta0.F,v 1.1.2.1 2002/04/04 10:58:59 heimbach Exp $
2
3 #include "COST_CPPOPTIONS.h"
4
5
6 subroutine cost_theta0(
7 I myiter,
8 I mytime,
9 I mythid
10 & )
11
12 c ==================================================================
13 c SUBROUTINE cost_zonstress
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_zonstress
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_THETA0_COST_CONTRIBUTION
94
95 if (optimcycle .ge. 0) then
96 ilfld = ilnblnk( xx_theta_file )
97 write(fnamefld(1:80),'(2a,i10.10)')
98 & xx_theta_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_theta_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 if ( ABS(R_low(i,j,bi,bj)) .LT. 1000. )
120 & tmpx = tmpx*ABS(R_low(i,j,bi,bj))/1000.
121 fctile = fctile
122 & + wtheta(k,bi,bj)*cosphi(i,j,bi,bj)
123 & *tmpx*tmpx
124 endif
125 enddo
126 enddo
127 enddo
128
129 objf_temp0(bi,bj) = objf_temp0(bi,bj) + fctile
130 fcthread = fcthread + fctile
131
132 #ifdef ECCO_VERBOSE
133 c-- Print cost function for each tile in each thread.
134 write(msgbuf,'(a)') ' '
135 call print_message( msgbuf, standardmessageunit,
136 & SQUEEZE_RIGHT , mythid)
137 write(msgbuf,'(a,i8.8,1x,i3.3,1x,i3.3)')
138 & ' cost_theta0: irec,bi,bj = ',irec,bi,bj
139 call print_message( msgbuf, standardmessageunit,
140 & SQUEEZE_RIGHT , mythid)
141 write(msgbuf,'(a,d22.15)')
142 & ' cost function (dT(0)) = ',
143 & fctile
144 call print_message( msgbuf, standardmessageunit,
145 & SQUEEZE_RIGHT , mythid)
146 #endif
147 enddo
148 enddo
149
150 #ifdef ECCO_VERBOSE
151 c-- Print cost function for all tiles.
152 _GLOBAL_SUM_R8( fcthread , myThid )
153 write(msgbuf,'(a)') ' '
154 call print_message( msgbuf, standardmessageunit,
155 & SQUEEZE_RIGHT , mythid)
156 write(msgbuf,'(a,i8.8)')
157 & ' cost_: irec = ',irec
158 call print_message( msgbuf, standardmessageunit,
159 & SQUEEZE_RIGHT , mythid)
160 write(msgbuf,'(a,d22.15)')
161 & ' global cost function value = ',
162 & fcthread
163 call print_message( msgbuf, standardmessageunit,
164 & SQUEEZE_RIGHT , mythid)
165 write(msgbuf,'(a)') ' '
166 call print_message( msgbuf, standardmessageunit,
167 & SQUEEZE_RIGHT , mythid)
168 #endif
169
170 #else
171 c-- Do not enter the calculation of the salinity increment
172 c-- contribution to the final cost function.
173
174 fctile = 0. _d 0
175 fcthread = 0. _d 0
176
177 #ifdef ECCO_VERBOSE
178 _BEGIN_MASTER( mythid )
179 write(msgbuf,'(a)') ' '
180 call print_message( msgbuf, standardmessageunit,
181 & SQUEEZE_RIGHT , mythid)
182 write(msgbuf,'(a,a)')
183 & ' cost_theta0 : no contribution of the I.C. in salin. ',
184 & ' to cost function.'
185 call print_message( msgbuf, standardmessageunit,
186 & SQUEEZE_RIGHT , mythid)
187 write(msgbuf,'(a)') ' '
188 call print_message( msgbuf, standardmessageunit,
189 & SQUEEZE_RIGHT , mythid)
190 _END_MASTER( mythid )
191 #endif
192
193 #endif
194
195 return
196 end
197
198

  ViewVC Help
Powered by ViewVC 1.1.22