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