/[MITgcm]/MITgcm/pkg/cost/cost_initvaria.F
ViewVC logotype

Contents of /MITgcm/pkg/cost/cost_initvaria.F

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


Revision 1.1.2.8 - (show annotations) (download)
Sat Jun 21 23:46:36 2003 UTC (21 years ago) by heimbach
Branch: ecco-branch
CVS Tags: ecco_c50_e33a, ecco_c51_e34, ecco_c51_e34e, ecco_c51_e34d, ecco_c51_e34f, ecco_c51_e34g, ecco_c51_e34a, ecco_c51_e34b, ecco_c51_e34c, ecco_c50_e33
Changes since 1.1.2.7: +3 -1 lines
Modifications, corrections and additions to obcs code.

1 C $Header: /u/gcmpack/MITgcm/pkg/cost/Attic/cost_initvaria.F,v 1.1.2.7 2003/06/20 05:30:14 heimbach Exp $
2
3 #include "COST_CPPOPTIONS.h"
4
5
6 subroutine cost_initvaria( mythid )
7
8 c ==================================================================
9 c SUBROUTINE cost_initvaria
10 c ==================================================================
11 c
12 c o Initialise the variable cost function part.
13 c
14 c started: Christian Eckert eckert@mit.edu 30-Jun-1999
15 c changed: Christian Eckert eckert@mit.edu 18-Apr-2000
16 c - Restructured the code in order to create a package
17 c for the MITgcmUV.
18 c added sea-ice term: menemenlis@jpl.nasa.gov 26-Feb-2003
19 c
20 c ==================================================================
21 c SUBROUTINE cost_initvaria
22 c ==================================================================
23
24 implicit none
25
26 c == global variables ==
27
28 #include "EEPARAMS.h"
29 #include "SIZE.h"
30 #include "GRID.h"
31
32 #include "cost.h"
33
34 c == routine arguments ==
35
36 integer mythid
37
38 c == local variables ==
39
40 integer bi,bj
41 integer itlo,ithi
42 integer jtlo,jthi
43 integer imin, imax
44 integer jmin, jmax
45 integer i,j,k
46
47 logical exst
48
49 c == external functions ==
50
51 c == end of interface ==
52 jtlo = mybylo(mythid)
53 jthi = mybyhi(mythid)
54 itlo = mybxlo(mythid)
55 ithi = mybxhi(mythid)
56 jmin = 1-OLy
57 jmax = sny+OLy
58 imin = 1-OLx
59 imax = snx+OLy
60
61 c-- Initialise adjoint of monthly mean files calculated
62 c-- in cost_averagesfields (and their ad...).
63 call cost_averagesinit( mythid )
64 _BARRIER
65 #ifdef ALLOW_TANGENTLINEAR_RUN
66 cph(
67 cph The following init. shoud not be applied if in the middle
68 cph of a divided adjoint run
69 cph)
70 inquire( file='costfinal', exist=exst )
71 if ( .NOT. exst) then
72 call cost_init_barfiles( mythid )
73 endif
74 #endif
75
76 c-- Initialize the tiled cost function contributions.
77 do bj = jtlo,jthi
78 do bi = itlo,ithi
79 objf_hflux(bi,bj) = 0. _d 0
80 objf_hfluxm(bi,bj) = 0. _d 0
81 objf_hfluxmm(bi,bj) = 0. _d 0
82 objf_sflux(bi,bj) = 0. _d 0
83 objf_sfluxm(bi,bj) = 0. _d 0
84 objf_sfluxmm(bi,bj) = 0. _d 0
85 objf_tauu(bi,bj) = 0. _d 0
86 objf_tauum(bi,bj) = 0. _d 0
87 objf_tauv(bi,bj) = 0. _d 0
88 objf_tauvm(bi,bj) = 0. _d 0
89 objf_temp(bi,bj) = 0. _d 0
90 objf_salt(bi,bj) = 0. _d 0
91 objf_temp0(bi,bj) = 0. _d 0
92 objf_salt0(bi,bj) = 0. _d 0
93 objf_tmi(bi,bj) = 0. _d 0
94 objf_sst(bi,bj) = 0. _d 0
95 objf_sss(bi,bj) = 0. _d 0
96 objf_h(bi,bj) = 0. _d 0
97 objf_atl(bi,bj) = 0. _d 0
98 objf_ctdt(bi,bj) = 0. _d 0
99 objf_ctds(bi,bj) = 0. _d 0
100 objf_ctdtclim(bi,bj) = 0. _d 0
101 objf_ctdsclim(bi,bj) = 0. _d 0
102 objf_xbt(bi,bj) = 0. _d 0
103 objf_argot(bi,bj) = 0. _d 0
104 objf_argos(bi,bj) = 0. _d 0
105 objf_drift(bi,bj) = 0. _d 0
106 objf_wdrift(bi,bj) = 0. _d 0
107 objf_sdrift(bi,bj) = 0. _d 0
108 objf_tdrift(bi,bj) = 0. _d 0
109 objf_scatx(bi,bj) = 0. _d 0
110 objf_scaty(bi,bj) = 0. _d 0
111 objf_scatxm(bi,bj) = 0. _d 0
112 objf_scatym(bi,bj) = 0. _d 0
113 objf_atemp(bi,bj) = 0. _d 0
114 objf_aqh(bi,bj) = 0. _d 0
115 objf_uwind(bi,bj) = 0. _d 0
116 objf_vwind(bi,bj) = 0. _d 0
117 objf_obcsn(bi,bj) = 0. _d 0
118 objf_obcss(bi,bj) = 0. _d 0
119 objf_obcsw(bi,bj) = 0. _d 0
120 objf_obcse(bi,bj) = 0. _d 0
121 objf_curmtr(bi,bj) = 0. _d 0
122 objf_ageos(bi,bj) = 0. _d 0
123 objf_ice(bi,bj) = 0. _d 0
124 objf_diffkr(bi,bj) = 0. _d 0
125 objf_kapgm(bi,bj) = 0. _d 0
126 objf_theta_ini_fin(bi,bj) = 0. _d 0
127 objf_salt_ini_fin(bi,bj) = 0. _d 0
128 enddo
129 enddo
130
131 k = 1
132 do bj = jtlo,jthi
133 do bi = itlo,ithi
134 do j = jmin,jmax
135 do i = imin,imax
136 #ifdef ALLOW_SSH_COST_CONTRIBUTION
137 if (_hFacC(i,j,k,bi,bj) .eq. 0.) then
138 tpmeanmask(i,j,bi,bj) = 0. _d 0
139 else
140 tpmeanmask(i,j,bi,bj) = 1. _d 0
141 endif
142 tpmean(i,j,bi,bj) = 0. _d 0
143 #endif
144 #ifdef ALLOW_SSH_TPANOM_COST_CONTRIBUTION
145 if (_hFacC(i,j,k,bi,bj) .eq. 0.) then
146 tpmask(i,j,bi,bj) = 0. _d 0
147 else
148 tpmask(i,j,bi,bj) = 1. _d 0
149 endif
150 tpobs(i,j,bi,bj) = 0. _d 0
151 #endif
152 #ifdef ALLOW_SSH_ERSANOM_COST_CONTRIBUTION
153 if (_hFacC(i,j,k,bi,bj) .eq. 0.) then
154 ersmask(i,j,bi,bj) = 0. _d 0
155 else
156 ersmask(i,j,bi,bj) = 1. _d 0
157 endif
158 ersobs(i,j,bi,bj) = 0. _d 0
159 #endif
160 #ifdef ALLOW_TMI_SST_COST_CONTRIBUTION
161 if (_hFacC(i,j,k,bi,bj) .eq. 0.) then
162 tmimask(i,j,bi,bj) = 0. _d 0
163 else
164 tmimask(i,j,bi,bj) = 1. _d 0
165 endif
166 tmidat(i,j,bi,bj) = 0. _d 0
167 #endif
168 #ifdef ALLOW_SST_COST_CONTRIBUTION
169 if (_hFacC(i,j,k,bi,bj) .eq. 0.) then
170 sstmask(i,j,bi,bj) = 0. _d 0
171 else
172 sstmask(i,j,bi,bj) = 1. _d 0
173 endif
174 sstdat(i,j,bi,bj) = 0. _d 0
175 #endif
176 #ifdef ALLOW_SSS_COST_CONTRIBUTION
177 if (_hFacC(i,j,k,bi,bj) .eq. 0.) then
178 sssmask(i,j,bi,bj) = 0. _d 0
179 else
180 sssmask(i,j,bi,bj) = 1. _d 0
181 endif
182 sssdat(i,j,bi,bj) = 0. _d 0
183 #endif
184 enddo
185 enddo
186 enddo
187 enddo
188
189 c-- Initialise the "global" parts of the cost function.
190 _BEGIN_MASTER( mythid )
191 objf_obcsvol = 0. _d 0
192 objf_hmean = 0. _d 0
193 fc = 0. _d 0
194 _END_MASTER( mythid )
195
196 _BARRIER
197
198 return
199 end
200

  ViewVC Help
Powered by ViewVC 1.1.22