/[MITgcm]/MITgcm_contrib/gael/matlab_class/ecco_v4/cost_altimeter.m
ViewVC logotype

Annotation of /MITgcm_contrib/gael/matlab_class/ecco_v4/cost_altimeter.m

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


Revision 1.15 - (hide annotations) (download)
Sat Nov 7 16:00:59 2015 UTC (9 years, 8 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint65r, checkpoint65q
Changes since 1.14: +32 -0 lines
- add another global mean plot.

1 gforget 1.7 function []=cost_altimeter(dirModel,dirMat);
2     %object: compute or plot the various sea level statistics
3     % (std model-obs, model, obs, leading to cost function terms)
4     %inputs: dirModel is the model run directory
5     % dirMat is the directory where diagnozed .mat files will be saved
6     % -> set it to '' to use the default [dirModel 'mat/']
7     % doComp is the switch from computation to plot
8 gforget 1.13
9 gforget 1.4 doComp=1;
10     if doComp==1;
11 gforget 1.11 doSave=1;
12     doPlot=0;
13 gforget 1.4 else;
14 gforget 1.11 doSave=0;
15     doPlot=1;
16 gforget 1.4 end;
17     doPrint=0;
18    
19 gforget 1.14 %%%%%%%%%%%
20     %load grid:
21     %%%%%%%%%%%
22    
23     gcmfaces_global;
24     if ~isfield(mygrid,'XC'); grid_load('./GRID/',5,'compact'); end;
25     if ~isfield(mygrid,'LATS_MASKS'); gcmfaces_lines_zonal; end;
26    
27     %%%%%%%%%%%%%%%
28     %define pathes:
29     %%%%%%%%%%%%%%%
30    
31     if ~isempty(dir([dirModel 'sigma_MDT_glob_eccollc.bin']));
32     dirSigma=dirModel;
33     else;
34     dirSigma='/net/nares/raid11/ecco-shared/ecco-version-4/input/input_err/';
35     end;
36    
37     if isempty(dir([dirModel 'barfiles']));
38     dirEcco=dirModel;
39     else;
40     dirEcco=[dirModel 'barfiles' filesep];
41     end;
42    
43     %nameSigma={'sigma_MDT_glob_eccollc.bin','sigma_SLA_smooth_eccollc.bin','sigma_SLA_PWS07r2_glob_eccollc.bin'}; maxLatObs=66
44     nameSigma={'sigma_MDT_glob_eccollc.bin','slaerr_largescale_r4.err','slaerr_gridscale_r4.err'}; maxLatObs=90
45    
46     if isempty(dirMat); dirMat=[dirModel 'mat/']; else; dirMat=[dirMat '/']; end;
47     if isempty(dir(dirMat)); mkdir([dirMat]); end;
48     if ~isdir([dirMat 'cost/']); mkdir([dirMat 'cost/']); end;
49     runName=pwd; tmp1=strfind(runName,filesep); runName=runName(tmp1(end)+1:end);
50    
51 gforget 1.15
52     %%%%%%%%%%%%%%%%%
53     %global means:
54     %%%%%%%%%%%%%%%%%
55    
56     tmpName=dir([dirEcco 'm_eta_day*data']);
57     avisoName=which('MSL_aviso.bin');
58     if length(tmpName)==1&~isempty(avisoName);
59     m_eta_day=rdmds2gcmfaces([dirEcco tmpName.name(1:end-5)]);
60     avisoName=which('MSL_aviso.bin');
61     etaglo_aviso=read2memory(avisoName)';
62     %global mean
63     nt=size(m_eta_day{1},3); etaglo=NaN*zeros(1,nt);
64     area=mygrid.mskC(:,:,1).*mygrid.RAC; areatot=nansum(area);
65     for tt=1:nt;
66     etaglo(tt)=nansum(area.*m_eta_day(:,:,tt))/areatot;
67     end;
68     %monthly mean
69     etaglo_mon=NaN*zeros(1,240);
70     etaglo_aviso_mon=NaN*zeros(1,240);
71     etadate=datenum([1992 1 1 12 0 0])+[1:nt]-1;
72     for mm=1:240;
73     t0=datenum([1992 1+(mm-1) 1 0 0 0]);
74     t1=datenum([1992 1+mm 1 0 0 0]);
75     tt=find(etadate>t0&etadate<t1);
76     etaglo_mon(mm)=mean(etaglo(tt));
77     etaglo_aviso_mon(mm)=mean(etaglo_aviso(tt));
78     end;
79     %
80     if doSave; eval(['save ' dirMat 'cost/cost_altimeter_etaglo.mat etaglo*;']); end;
81     end;
82    
83 gforget 1.14 %%%%%%%%%%%%%%%%%
84     %do computations:
85     %%%%%%%%%%%%%%%%%
86    
87 gforget 1.4 for doDifObsOrMod=1:3;
88 gforget 1.14
89 gforget 1.11 if doDifObsOrMod==1; suf='modMobs'; elseif doDifObsOrMod==2; suf='obs'; else; suf='mod'; end;
90    
91     if doComp==0;
92     eval(['load ' dirMat 'cost/cost_altimeter_' suf '.mat myflds;']);
93     else;
94    
95     tic;
96     %mdt cost function term (misfit plot)
97 gforget 1.14 dif_mdt=rdmds2gcmfaces([dirEcco 'mdtdiff_smooth']);
98 gforget 1.12 sig_mdt=read_bin([dirSigma nameSigma{1}],1,0);
99 gforget 1.11 sig_mdt(find(sig_mdt==0))=NaN;
100     %store:
101     myflds.dif_mdt=dif_mdt;
102     myflds.sig_mdt=sig_mdt;
103    
104     %skip blanks:
105 gforget 1.14 tmp1=dir([dirEcco 'sladiff_smooth*data']);
106 gforget 1.11 nrec=tmp1.bytes/90/1170/4;
107     ttShift=17;
108     %skip 1992
109     listRecs=[1+365-ttShift:17+365*18-ttShift];
110     nRecs=length(listRecs); TT=1993+[0:nRecs-1]/365.25;
111    
112     toc; tic;
113     %pointwise/1day terms:
114     sladiff_point=repmat(0*mygrid.RAC,[1 1 nRecs]); count_point=sladiff_point;
115     for ii=1:3;
116     if ii==1; myset='tp'; elseif ii==2; myset='gfo'; else; myset='ers'; end;
117     %topex pointwise misfits:
118     if doDifObsOrMod==1;
119 gforget 1.14 sladiff_tmp=cost_altimeter_read([dirEcco 'sladiff_' myset '_raw'],ttShift+listRecs);
120 gforget 1.11 elseif doDifObsOrMod==2;
121 gforget 1.14 sladiff_tmp=cost_altimeter_read([dirEcco 'slaobs_' myset '_raw'],ttShift+listRecs);
122 gforget 1.11 else;
123 gforget 1.14 sladiff_tmp=cost_altimeter_read([dirEcco 'sladiff_' myset '_raw'],ttShift+listRecs);
124     sladiff_tmp=sladiff_tmp+cost_altimeter_read([dirEcco 'slaobs_' myset '_raw'],ttShift+listRecs);
125 gforget 1.11 end;
126     %add to overall data set:
127     sladiff_point=sladiff_point+sladiff_tmp; count_point=count_point+(sladiff_tmp~=0);
128     %finalize data set:
129     sladiff_tmp(sladiff_tmp==0)=NaN;
130     %compute rms:
131     count_tmp=sum(~isnan(sladiff_tmp),3); count_tmp(find(count_tmp<10))=NaN; msk_tmp=1+0*count_tmp;
132     eval(['myflds.rms_' myset '=sqrt(nansum(sladiff_tmp.^2,3)./count_tmp);']);
133     eval(['myflds.std_' myset '=nanstd(sladiff_tmp,0,3).*msk_tmp;']);
134     end;
135     %finalize overall data set:
136     count_point(count_point==0)=NaN;
137     sladiff_point=sladiff_point./count_point;
138     %compute overall rms,std:
139     count_tmp=sum(~isnan(sladiff_point),3); count_tmp(find(count_tmp<10))=NaN; msk_tmp=1+0*count_tmp;
140     myflds.rms_sladiff_point=sqrt(nansum(sladiff_point.^2,3)./count_tmp);
141     myflds.std_sladiff_point=nanstd(sladiff_point,0,3).*msk_tmp;
142     %
143     msk_point=1*(count_point>0);
144     %fill blanks:
145     warning('off','MATLAB:divideByZero');
146     msk=mygrid.mskC(:,:,1); msk(find(abs(mygrid.YC)>maxLatObs))=NaN;
147     myflds.xtrp_rms_sladiff_point=diffsmooth2D_extrap_inv(myflds.rms_sladiff_point,msk);
148     myflds.xtrp_std_sladiff_point=diffsmooth2D_extrap_inv(myflds.std_sladiff_point,msk);
149     warning('on','MATLAB:divideByZero');
150     %get weight:
151 gforget 1.12 sig_sladiff_point=read_bin([dirSigma nameSigma{3}],1,0);
152 gforget 1.11 sig_sladiff_point(find(sig_sladiff_point==0))=NaN;
153     myflds.sig_sladiff_point=sig_sladiff_point;
154    
155     %computational mask : only points that were actually observed, and in 35d average
156 gforget 1.14 msk_point35d=cost_altimeter_read([dirEcco 'sladiff_raw'],listRecs);
157 gforget 1.11 msk_point35d=1*(msk_point35d~=0);
158     tmp1=sum(msk_point35d==0&msk_point~=0);
159     tmp2=sum(msk_point35d~=0&msk_point~=0);
160     fprintf('after masking : %d omitted, %d retained \n',tmp1,tmp2);
161     msk_point=msk_point.*msk_point35d;
162     msk_point(msk_point==0)=NaN;
163     %plotting mask : regions of less than 100 observations are omitted
164     msk_100pts=1*(nansum(msk_point,3)>=100);
165     msk_100pts(msk_100pts==0)=NaN;
166     %store
167     myflds.msk_100pts=msk_100pts;
168    
169     toc; tic;
170     %lsc cost function term:
171     if doDifObsOrMod==1;
172 gforget 1.14 sladiff_smooth=cost_altimeter_read([dirEcco 'sladiff_smooth'],listRecs);
173 gforget 1.11 elseif doDifObsOrMod==2;
174 gforget 1.14 sladiff_smooth=cost_altimeter_read([dirEcco 'slaobs_smooth'],listRecs);
175 gforget 1.11 else;
176 gforget 1.14 sladiff_smooth=cost_altimeter_read([dirEcco 'sladiff_smooth'],listRecs);
177     sladiff_smooth=sladiff_smooth+cost_altimeter_read([dirEcco 'slaobs_smooth'],listRecs);
178 gforget 1.11 end;
179     %mask missing points:
180     sladiff_smooth(sladiff_smooth==0)=NaN;
181     sladiff_smooth=sladiff_smooth.*msk_point;
182     count_tmp=sum(~isnan(sladiff_smooth),3); count_tmp(find(count_tmp<10))=NaN; msk_tmp=1+0*count_tmp;
183     %compute rms:
184     rms_sladiff_smooth=msk_tmp.*sqrt(nanmean(sladiff_smooth.^2,3));
185     std_sladiff_smooth=msk_tmp.*nanstd(sladiff_smooth,0,3);
186     %get weight:
187 gforget 1.12 sig_sladiff_smooth=read_bin([dirSigma nameSigma{2}],1,0);
188 gforget 1.11 sig_sladiff_smooth(find(sig_sladiff_smooth==0))=NaN;
189     %store:
190     myflds.rms_sladiff_smooth=rms_sladiff_smooth;
191     myflds.std_sladiff_smooth=std_sladiff_smooth;
192     myflds.sig_sladiff_smooth=sig_sladiff_smooth;
193    
194     toc; tic;
195     %pointwise/point35days cost function term:
196     if doDifObsOrMod==1;
197 gforget 1.14 sladiff_point35d=cost_altimeter_read([dirEcco 'sladiff_raw'],listRecs);
198 gforget 1.11 elseif doDifObsOrMod==2;
199 gforget 1.14 sladiff_point35d=cost_altimeter_read([dirEcco 'slaobs_raw'],listRecs);
200 gforget 1.11 else;
201 gforget 1.14 sladiff_point35d=cost_altimeter_read([dirEcco 'sladiff_raw'],listRecs);
202     sladiff_point35d=sladiff_point35d+cost_altimeter_read([dirEcco 'slaobs_raw'],listRecs);
203 gforget 1.11 end;
204     %mask missing points:
205     sladiff_point35d(sladiff_point35d==0)=NaN;
206     sladiff_point35d=sladiff_point35d.*msk_point;
207     count_tmp=sum(~isnan(sladiff_point35d),3); count_tmp(find(count_tmp<10))=NaN; msk_tmp=1+0*count_tmp;
208     %compute rms:
209     rms_sladiff_point35d=msk_tmp.*sqrt(nanmean(sladiff_point35d.^2,3));
210     std_sladiff_point35d=msk_tmp.*nanstd(sladiff_point35d,0,3);
211     %store:
212     myflds.rms_sladiff_point35d=rms_sladiff_point35d;
213     myflds.std_sladiff_point35d=std_sladiff_point35d;
214    
215     if 0;
216     %difference between scales
217     rms_sladiff_35dMsmooth=msk_tmp.*sqrt(nanmean((sladiff_point35d-sladiff_smooth).^2,3));
218     std_sladiff_35dMsmooth=msk_tmp.*nanstd((sladiff_point35d-sladiff_smooth),0,3);
219     %store:
220     myflds.rms_sladiff_35dMsmooth=rms_sladiff_35dMsmooth;
221     myflds.std_sladiff_35dMsmooth=std_sladiff_35dMsmooth;
222    
223     %difference between scales
224     rms_sladiff_pointMpoint35d=msk_tmp.*sqrt(nanmean((sladiff_point-sladiff_point35d).^2,3));
225     std_sladiff_pointMpoint35d=msk_tmp.*nanstd((sladiff_point-sladiff_point35d),0,3);
226     %store:
227     myflds.rms_sladiff_pointMpoint35d=rms_sladiff_pointMpoint35d;
228     myflds.std_sladiff_pointMpoint35d=std_sladiff_pointMpoint35d;
229     end;
230    
231     %compute zonal mean/median:
232     for ii=1:4;
233     switch ii;
234     case 1; tmp1='mdt'; cost_fld=(mygrid.mskC(:,:,1).*myflds.dif_mdt./myflds.sig_mdt).^2;
235     case 2; tmp1='lsc'; cost_fld=(mygrid.mskC(:,:,1).*myflds.rms_sladiff_smooth./myflds.sig_sladiff_smooth).^2;
236     case 3; tmp1='point35d'; cost_fld=(mygrid.mskC(:,:,1).*myflds.rms_sladiff_point35d./myflds.sig_sladiff_point).^2;
237     case 4; tmp1='point'; cost_fld=(mygrid.mskC(:,:,1).*myflds.rms_sladiff_point./myflds.sig_sladiff_point).^2;
238     end;
239     cost_zmean=calc_zonmean_T(cost_fld); eval(['mycosts_mean.' tmp1 '=cost_zmean;']);
240     cost_zmedian=calc_zonmedian_T(cost_fld); eval(['mycosts_median.' tmp1 '=cost_zmedian;']);
241     end;
242    
243     toc; %write to disk:
244     if doSave; eval(['save ' dirMat 'cost/cost_altimeter_' suf '.mat myflds mycosts_mean mycosts_median;']); end;
245    
246     end;%if doComp
247    
248     if doPlot;
249     cc=[-0.4:0.05:-0.25 -0.2:0.03:-0.05 -0.03:0.01:0.03 0.05:0.03:0.2 0.25:0.05:0.4];
250     figure; m_map_gcmfaces(myflds.dif_mdt,0,{'myCaxis',cc}); drawnow;
251     cc=[0:0.005:0.02 0.03:0.01:0.05 0.06:0.02:0.1 0.14:0.03:0.2 0.25:0.05:0.4];
252     figure; m_map_gcmfaces(myflds.rms_sladiff_smooth,0,{'myCaxis',cc}); drawnow;
253     figure; m_map_gcmfaces(myflds.rms_sladiff_point35d,0,{'myCaxis',cc}); drawnow;
254     figure; m_map_gcmfaces(myflds.rms_sladiff_point,0,{'myCaxis',cc}); drawnow;
255     end;
256    
257     if doPlot&doPrint;
258     dirFig='../figs/altimeter/'; ff0=gcf-4;
259     for ff=1:4;
260     figure(ff+ff0); saveas(gcf,[dirFig runName '_' suf num2str(ff)],'fig');
261     eval(['print -depsc ' dirFig runName '_' suf num2str(ff) '.eps;']);
262     eval(['print -djpeg90 ' dirFig runName '_' suf num2str(ff) '.jpg;']);
263     end;
264     end;
265    
266 gforget 1.4 end;%for doDifObsOrMod=1:3;
267    
268    
269    
270     function [fldOut]=cost_altimeter_read(fileIn,recIn);
271    
272     nrec=length(recIn);
273     global mygrid; siz=[size(convert2gcmfaces(mygrid.XC)) nrec];
274     lrec=siz(1)*siz(2)*4;
275     myprec='float32';
276    
277     fldOut=zeros(siz);
278     fid=fopen([fileIn '.data'],'r','b');
279     for irec=1:nrec;
280 gforget 1.11 status=fseek(fid,(recIn(irec)-1)*lrec,'bof');
281     fldOut(:,:,irec)=fread(fid,siz(1:2),myprec);
282 gforget 1.4 end;
283    
284     fldOut=convert2gcmfaces(fldOut);
285    
286    
287    

  ViewVC Help
Powered by ViewVC 1.1.22