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

Annotation of /MITgcm/pkg/ecco/cost_ssh_mean.F

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


Revision 1.2 - (hide annotations) (download)
Mon Oct 11 16:38:53 2004 UTC (19 years, 8 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57d_post, checkpoint57b_post, checkpoint57c_pre, checkpoint55j_post, checkpoint56b_post, checkpoint55h_post, checkpoint57e_post, checkpoint56c_post, checkpoint57a_post, checkpoint55g_post, checkpoint55f_post, checkpoint57a_pre, checkpoint55i_post, checkpoint57, checkpoint56, eckpoint57e_pre, checkpoint57c_post, checkpoint55e_post, checkpoint56a_post, checkpoint55d_post
Changes since 1.1: +0 -0 lines
o ECCO specific cost function terms (up-to-date with 1x1 runs)
o ecco_cost_weights is modified to 1x1 runs
o modifs to allow observations to be read in as
  single file or yearly files

1 heimbach 1.1 C $Header: /u/gcmpack/MITgcm/pkg/cost/Attic/cost_ssh_mean.F,v 1.1.2.3 2002/04/04 10:58:59 heimbach Exp $
2    
3     #include "COST_CPPOPTIONS.h"
4    
5    
6     subroutine cost_ssh_mean(
7     I psmean, offset
8     O , costmean
9     I , mythid
10     & )
11    
12     c ==================================================================
13     c SUBROUTINE cost_ssh_mean
14     c ==================================================================
15     c
16     c o Evaluate cost function contribution of sea surface height.
17     c using of geoid error covariances requires regular model grid
18     c o TODO: interpolate model grid to regular grid
19     c check mask
20     c
21     c started: Detlef Stammer, Ralf Giering Jul-1996
22     c changed: Ralf Giering Ralf.Giering@FastOpt.de 12-Jun-2001
23     c - totally rewrite for parallel processing
24     c heimbach@mit.edu 13-Mar-2002
25     c - several wrong i-loop boundaries (spotted by G. Gebbie)
26     c - nprocs need to be replaced by nPx*nPy
27     c - geoid error does not work as of now
28     c
29     c ==================================================================
30     c SUBROUTINE cost_ssh_mean
31     c ==================================================================
32    
33     implicit none
34    
35     c == global variables ==
36    
37     #include "EEPARAMS.h"
38     #include "SIZE.h"
39     #include "PARAMS.h"
40    
41     #include "ecco_cost.h"
42    
43     #ifdef ALLOW_EGM96_ERROR_COV
44     # include "EESUPPORT.h"
45     # include "cost_sph.h"
46     #endif
47    
48     c == routine arguments ==
49    
50     _RL psmean ( 1-olx:snx+olx, 1-oly:sny+oly, nsx, nsy )
51     _RL offset
52     _RL costmean
53     integer mythid
54    
55     c == local variables ==
56    
57     integer bi,bj
58     integer i,j
59     integer itlo,ithi
60     integer jtlo,jthi
61     integer jmin,jmax
62     integer imin,imax
63    
64     #ifdef ALLOW_EGM96_ERROR_COV
65     _RL cphi
66     _RL xmean
67     _RL shc( ncshc )
68     _RL misfit ( 1-olx:snx+olx, 1-oly:sny+oly, nsx,nsy )
69     _RL misfitgl( 1-olx:snx+olx, 1-oly:sny+oly, nsx,nsy, npx,npy )
70    
71     _RL pwork ( (lshcmax+1)*(lshcmax+2)/2 )
72     _RL diffearth ( lonshc, latshc )
73     integer joffset
74     integer ipx, ipy
75     integer iglobe, jglobe
76     integer iproc
77     integer mpirc
78     integer mpicrd(2)
79     #endif
80    
81     _RL diff
82     Real*8 sumc
83    
84     c == end of interface ==
85    
86     c-- Initialise local variables.
87     jtlo = mybylo(mythid)
88     jthi = mybyhi(mythid)
89     itlo = mybxlo(mythid)
90     ithi = mybxhi(mythid)
91    
92     jmin = 1
93     jmax = sny
94     imin = 1
95     imax = snx
96    
97     costmean = 0.
98    
99     #ifdef ALLOW_EGM96_ERROR_COV
100    
101     c-- compute misfit on local domain
102     do bj = jtlo,jthi
103     do bi = itlo,ithi
104     do j = jmin,jmax
105     do i = imin,imax
106     if ( tpmeanmask(i,j,bi,bj) .gt. 0. ) then
107     misfit(i,j,bi,bj) =
108     & tpmean(i,j,bi,bj) - offset - psmean(i,j,bi,bj)
109     else
110     misfit(i,j,bi,bj) = 0.
111     endif
112     enddo
113     enddo
114     enddo
115     enddo
116    
117     c-- communicate to get misfit on full 2d model surface grid
118     call exch_all_2d_rl( misfit, misfitgl, mythid )
119    
120     c-- set meridional index offset,
121     c-- it is non zero if the grid does no reach the poles
122     joffset = (latshc - ny)/2
123    
124     _BEGIN_MASTER( mythid )
125     print *
126     print *, ' cost_SSH: grid starts at point ', joffset
127     & , ' of full earth.'
128     _END_MASTER( mythid )
129    
130     c-- preset field
131     c-- necessary if the grid does no reach the poles
132     do j = 1, latshc
133     do i = 1, lonshc
134     diffearth(i,j) = 0.
135     enddo
136     enddo
137    
138     c-- -interpolate- from model grid to regular grid
139     c-- so far we assume that the model grid is already regular
140     _BEGIN_MASTER( mythid )
141    
142     do iproc = 1, nprocs
143    
144     c-- get coordinates of processor (iporc-1)
145     call MPI_Cart_coords(
146     I MPI_COMM_MODEL, iproc-1, 2, mpicrd
147     O , mpirc
148     & )
149    
150     ipx = 1 + mpicrd(1)
151     ipy = 1 + mpicrd(2)
152    
153     do bj = jtlo,jthi
154     do bi = itlo,ithi
155    
156     do j = jmin,jmax
157     do i = imin,imax
158     jglobe = joffset+ j + sNy*(bj-1) + sNy*nSy*(ipy-1)
159     iglobe = i + sNx*(bi-1) + sNx*nSx*(ipx-1)
160    
161     diffearth(iglobe,jglobe) = misfitgl(i,j,bi,bj,ipx,ipy)
162     enddo
163     enddo
164    
165     enddo
166     enddo
167    
168     enddo
169    
170     c-- Project regular grid misfit onto sperical harmonics
171     call shc4grid(
172     I lshcmax
173     O , shc
174     I , latshc, lonshc
175     I , diffearth
176     W , pwork
177     & )
178    
179     c-- Remove the C(0,0) component, i.e. the global mean.
180     shc(1) = 0.
181    
182     C-- Compute the cost function for the mean SSH.
183     call cost_geoid(
184     O costmean
185     I , shc
186     I , mythid
187     & )
188    
189     _END_MASTER( mythid )
190    
191     _BARRIER
192    
193     #else
194    
195     c-- Compute cost function for SSH by using the diagonal
196     c-- of the error covariance only.
197     c--
198     c-- Note: wp is assumed to include latitude dependence
199     c-- of error due to meridian convergence;
200     c-- --> no weighting by cosphi.
201    
202     sumc = 0.
203     do bj = jtlo,jthi
204     do bi = itlo,ithi
205     do j = jmin,jmax
206     do i = imin,imax
207     diff = psmean(i,j,bi,bj) - tpmean(i,j,bi,bj) + offset
208     sumc = sumc + diff*diff
209     & * wp(i,j,bi,bj)*tpmeanmask(i,j,bi,bj)
210     ce --> check masks!
211     enddo
212     enddo
213     enddo
214     enddo
215    
216     c-- Do the global summation.
217     _GLOBAL_SUM_R8( sumc, mythid )
218     costmean = sumc
219     #endif
220    
221     end

  ViewVC Help
Powered by ViewVC 1.1.22