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