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