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

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

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


Revision 1.10 - (show annotations) (download)
Tue Sep 4 15:04:51 2012 UTC (11 years, 9 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f, checkpoint63s, checkpoint64, checkpoint65, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65f, checkpoint65d, checkpoint65e
Changes since 1.9: +4 -4 lines
- remove #ifdef ALLOW_SMOOTH_CORREL3D brackets.
- add more relevant #ifdef ALLOW_SMOOTH ones.
- sort out useAtmWind, useSMOOTH, ctrlSmoothCorrel2D.

1 C $Header: /u/gcmpack/MITgcm/pkg/ecco/cost_kapgm.F,v 1.9 2012/08/10 19:45:26 jmc Exp $
2 C $Name: $
3
4 #include "ECCO_OPTIONS.h"
5
6
7 subroutine cost_kapgm(
8 I myiter,
9 I mytime,
10 I mythid
11 & )
12
13 C o==========================================================o
14 C | subroutine cost_kapgm |
15 C | o GM 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_GMREDI
27 # include "GMREDI.h"
28 #endif
29
30 #include "ecco_cost.h"
31 #include "CTRL_SIZE.h"
32 #include "ctrl.h"
33 #include "ctrl_dummy.h"
34 #include "optim.h"
35 c == routine arguments ==
36
37 integer myiter
38 _RL mytime
39 integer mythid
40
41 c == local variables ==
42
43 integer bi,bj
44 integer i,j,k
45 integer itlo,ithi
46 integer jtlo,jthi
47 integer jmin,jmax
48 integer imin,imax
49 integer nrec
50 integer irec
51 integer ilfld
52
53 _RL fctile
54 _RL fcthread
55 _RL tmpx
56
57 logical doglobalread
58 logical ladinit
59
60 character*(80) fnamefld
61
62 character*(MAX_LEN_MBUF) msgbuf
63
64 c == external functions ==
65
66 integer ilnblnk
67 external ilnblnk
68
69 c == end of interface ==
70
71 jtlo = mybylo(mythid)
72 jthi = mybyhi(mythid)
73 itlo = mybxlo(mythid)
74 ithi = mybxhi(mythid)
75 jmin = 1
76 jmax = sny
77 imin = 1
78 imax = snx
79
80 c-- Read state record from global file.
81 doglobalread = .false.
82 ladinit = .false.
83
84 irec = 1
85
86 #ifdef ALLOW_KAPGM_COST_CONTRIBUTION
87
88 if (optimcycle .ge. 0) then
89 ilfld = ilnblnk( xx_kapgm_file )
90 write(fnamefld(1:80),'(2a,i10.10)')
91 & xx_kapgm_file(1:ilfld),'.',optimcycle
92 endif
93
94 fcthread = 0. _d 0
95
96 call active_read_xyz( fnamefld, tmpfld3d, irec, doglobalread,
97 & ladinit, optimcycle, mythid
98 & , xx_kapgm_dummy )
99
100 c-- Loop over this thread tiles.
101 do bj = jtlo,jthi
102 do bi = itlo,ithi
103
104 c-- Determine the weights to be used.
105
106 fctile = 0. _d 0
107 do k = 1,nr
108 do j = jmin,jmax
109 do i = imin,imax
110 if (_hFacC(i,j,k,bi,bj) .ne. 0.) then
111 c tmpx = (tmpfld3d(i,j,k,bi,bj)-GM_background_K)
112 tmpx = tmpfld3d(i,j,k,bi,bj)
113 IF ( .NOT.ctrlSmoothCorrel3D ) THEN
114 fctile = fctile
115 & + wkapgmFld(i,j,k,bi,bj)*cosphi(i,j,bi,bj)
116 & *tmpx*tmpx
117 ELSE !IF ( .NOT.ctrlSmoothCorrel3D ) THEN
118 fctile = fctile + tmpx*tmpx
119 ENDIF !IF ( .NOT.ctrlSmoothCorrel3D ) THEN
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_RL( fcthread , myThid )
149 write(msgbuf,'(a)') ' '
150 call print_message( msgbuf, standardmessageunit,
151 & SQUEEZE_RIGHT , mythid)
152 write(msgbuf,'(a,i8.8)')
153 & ' cost_kapgm: 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
168 fctile = 0. _d 0
169 fcthread = 0. _d 0
170
171 #ifdef ECCO_VERBOSE
172 _BEGIN_MASTER( mythid )
173 write(msgbuf,'(a)') ' '
174 call print_message( msgbuf, standardmessageunit,
175 & SQUEEZE_RIGHT , mythid)
176 write(msgbuf,'(a)')
177 & ' cost_kapgm : no contribution to cost function'
178 call print_message( msgbuf, standardmessageunit,
179 & SQUEEZE_RIGHT , mythid)
180 write(msgbuf,'(a)') ' '
181 call print_message( msgbuf, standardmessageunit,
182 & SQUEEZE_RIGHT , mythid)
183 _END_MASTER( mythid )
184 #endif
185
186 #endif
187
188 return
189 end
190
191

  ViewVC Help
Powered by ViewVC 1.1.22