/[MITgcm]/MITgcm_contrib/ESMF/global_ocean.128x64x15/diags_matlab/ncep2global_ocean.m
ViewVC logotype

Annotation of /MITgcm_contrib/ESMF/global_ocean.128x64x15/diags_matlab/ncep2global_ocean.m

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


Revision 1.1.1.1 - (hide annotations) (download) (vendor branch)
Sun Feb 15 22:28:30 2004 UTC (21 years, 5 months ago) by cnh
Branch: MAIN, Initial
CVS Tags: adoption_1_0_pre_A, Baseline, HEAD
Changes since 1.1: +0 -0 lines
Initial checkin

1 cnh 1.1 function ncep2global_ocean
2     % function ncep2global_ocean
3    
4     % read various fluxes, adds them up to have net heat flux and net
5     % freshwater flux (E-P only)
6     % interpolates the fluxes onto 90x40 grid of global_ocean.90x40x15
7    
8     % constants
9     ql_evap = 2.5e6; % latent heat due to evaporation
10     rho_fresh = 1000; % density of fresh water
11     onemonth = 60*60*24*30;
12     % grid parameter
13     grd = mit_loadgrid(DIRECTORY WHERE THE GRID INFORMATION IS HELD ...
14     (XC*, YC*, HFAC* etc.));
15     lonc = grd.lonc;
16     latc = grd.latc;
17     cmask = grd.cmask; % 1=wet point, NaN = dry point
18     rac = grd.rac;
19    
20     bdir = DIRECTORY WHERE THE NCEP DATA FILES ARE STORED;
21     % taux
22     ncload(fullfile(bdir,'ncep_taux_monthly.cdf'),'taux','X','Y','T')
23     % tauy
24     ncload(fullfile(bdir,'ncep_tauy_monthly.cdf'),'tauy')
25     % there appears to be a different east west convention in the data
26     taux = - taux;
27     tauy = - tauy;
28    
29     % qnet
30     ncload(fullfile(bdir,'ncep_netsolar_monthly.cdf'),'solr')
31     ncload(fullfile(bdir,'ncep_netlw_monthly.cdf'),'lwflx') % long wave radiation
32     ncload(fullfile(bdir,'ncep_latent_monthly.cdf'),'heat_flux');
33     lhflx = heat_flux;
34     ncload(fullfile(bdir,'ncep_sensible_monthly.cdf'),'heat_flux');
35     shflx = heat_flux;
36    
37     % emp
38     ncload(fullfile(bdir,'ncep_precip_monthly.cdf'),'prate')
39     ncload(fullfile(bdir,'ncep_runoff_monthly.cdf'),'runoff')
40     evap = lhflx/ql_evap/rho_fresh; % evaporation rate estimated from latent
41     % heat flux
42     precip = prate/rho_fresh; % change the units
43    
44     % add the fields
45     % net heat flux: sum of solar radiation, long wave flux, latent heat flux,
46     % and sensible heat flux
47     qnet = solr + lwflx + lhflx + shflx;
48     % net fresh water flux: difference between evaporation and precipitation
49     emp = evap - precip;
50     % substract the runoff (doesn't lead anywhere, unfortunately)
51     empmr = emp - change(runoff,'==',NaN,0)/rho_fresh/onemonth;
52     runoff(find(isnan(runoff)))=0;
53     empmr = emp - runoff/rho_fresh/onemonth;
54    
55     % interpolate the three fields onto 90x40 grid
56     ir = interp2global(X,Y,lonc,latc); % this sets up the
57     % interpolator ir!
58     for k = 1:length(T);
59     tauxg(:,:,k) = interp2global(X,Y,squeeze(taux(k,:,:)),lonc,latc,ir);
60     tauyg(:,:,k) = interp2global(X,Y,squeeze(tauy(k,:,:)),lonc,latc,ir);
61     qnetg(:,:,k) = interp2global(X,Y,squeeze(qnet(k,:,:)),lonc,latc,ir);
62     empg(:,:,k) = interp2global(X,Y,squeeze(emp(k,:,:)),lonc,latc,ir);
63     empmrg(:,:,k) = interp2global(X,Y,squeeze(empmr(k,:,:)),lonc,latc,ir);
64     end
65    
66     % apply landmask
67     for k = 1:length(T);
68     tauxg(:,:,k) = tauxg(:,:,k).*cmask(:,:,1);
69     tauyg(:,:,k) = tauyg(:,:,k).*cmask(:,:,1);
70     qnetg(:,:,k) = qnetg(:,:,k).*cmask(:,:,1);
71     empg(:,:,k) = empg(:,:,k).*cmask(:,:,1);
72     empmrg(:,:,k) = empmrg(:,:,k).*cmask(:,:,1);
73     end
74     % balance the fluxes
75     effrac = rac.*cmask(:,:,1);
76     mqnet = mit_mean(qnetg,effrac);
77     qnetb = qnetg-mean(mqnet);
78    
79     memp = mit_mean(empg,effrac);
80     empb = empg-mean(memp);
81     mempmr = mit_mean(empmrg,effrac);
82     empmrb = empmrg-mean(mempmr);
83    
84     % apply the landmasks one more time (because we have added something
85     % to what should have been zero == land)
86     for k = 1:length(T);
87     qnetb(:,:,k) = qnetb(:,:,k).*smask;
88     empb(:,:,k) = empb(:,:,k).*smask;
89     empmrb(:,:,k) = empmrb(:,:,k).*smask;
90     end
91    
92     % save the fields
93     acc = 'real*4';
94     mit_writefield('ncep_taux.bin',change(tauxg,'==',NaN,0),acc);
95     mit_writefield('ncep_tauy.bin',change(tauyg,'==',NaN,0),acc);
96     mit_writefield('ncep_qnet.bin',change(qnetb,'==',NaN,0),acc);
97     mit_writefield('ncep_emp.bin',change(empb,'==',NaN,0),acc);
98     mit_writefield('ncep_empmr.bin',change(empmrb,'==',NaN,0),acc);
99    
100     return
101    
102     %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
103    
104     function Y = mit_mean(X,area);
105     %function Y = mit_mean(X,area);
106    
107     na = ndims(area);
108     n = ndims(X);
109     full_area = nansum(area(:));
110     if n == 2 & na == 2
111     [nx ny] = size(X);
112     Xa = X.*area;
113     Y = nansum(Xa(:))/full_area;
114     elseif n == 3 & na == 2
115     [nx ny nt] = size(X);
116     Xa = X.*repmat(area,[1 1 nt]);
117     for k=1:nt
118     Y(k) = nansum(nansum(Xa(:,:,k)))/full_area;
119     end
120     elseif n == 3 & na == 3
121     [nx ny nz] = size(X);
122     Xa = X.*area;
123     Y = nansum(Xa(:))/full_area;
124     else
125     error('X and area have to have consistent dimensions')
126     end
127    
128     return
129    
130     %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
131    
132     function count = mit_writefield(fname,h,accuracy)
133     %function count = mit_writefield(fname,h,accuracy)
134    
135     ieee='ieee-be';
136    
137     [fid message] = fopen(fname,'w',ieee);
138     if fid <= 0
139     error([message ', filename: ', [fname]])
140     end
141    
142     clear count
143     dims = size(h);
144     if length(dims) == 2
145     count = fwrite(fid,h,accuracy);
146     elseif length(dims) == 3
147     for k=1:dims(3);
148     count(k) = fwrite(fid,squeeze(h(:,:,k)),accuracy);
149     end
150     elseif length(dims) == 4
151     for kt=1:dims(4)
152     for k=1:dims(3)
153     count(k,kt) = fwrite(fid,squeeze(h(:,:,k,kt)),accuracy);
154     end
155     end
156     end
157    
158     fclose(fid);
159    
160     return
161    
162     %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
163    
164     function [zi] = interp2global(varargin);
165     %function [zi] = interp2global(x,y,z,xi,yi,ir);
166    
167     if nargin == 4
168     [xx yy] = meshgrid(varargin{1},varargin{2});
169     xx = reshape(xx,prod(size(xx)),1);
170     yy = reshape(yy,prod(size(yy)),1);
171     [xxi yyi] = meshgrid(varargin{3},varargin{4});
172     xxv = reshape(xxi,prod(size(xxi)),1);
173     yyv = reshape(yyi,prod(size(yyi)),1);
174     % create map
175     if ~exist('ir','var');
176     r = 230e3; %m
177     for k = 1:length(xxv);
178     % find all points within radius r
179     radius = lonlat2dist(xxv(k),yyv(k),xx,yy);
180     ir{k} = find(radius<=r);
181     end
182     end
183     zi = ir;
184     else
185     % zi = interp2_global(x,y,z,ix,iy','spline')';
186     z = varargin{3};
187     xi = varargin{4};
188     yi = varargin{5};
189     [xxi yyi] = meshgrid(varargin{4},varargin{5});
190     ir = varargin{6};
191     zzv = repmat(NaN,length(ir),1);
192     % interpolate
193     for k = 1:length(ir);
194     % find all points within radius r
195     zzv(k) = mean(z(ir{k}));
196     end
197     zi = reshape(zzv,size(xxi))';
198     end
199    
200     return
201    
202     %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
203    
204     function r = lonlat2dist(lon,lat,lons,lats)
205     %function r = lonlat2dist(lon,lat,lons,lats)
206     %
207     % returns distance (in meters) between postion (lon,lat) and positions
208     % (lons,lats) on the earth (sphere). length(r) = length(lons)
209    
210     earth=6371000;
211     deg2rad=pi/180;
212     nx=length(lon);
213     lt = deg2rad*lat;
214     ln = deg2rad*lon;
215     lts = deg2rad*lats;
216     lns = deg2rad*lons;
217     alpha = ...
218     acos( ( cos(lt)*cos(lts) ).*cos(ln-lns) ...
219     + ( sin(lt)*sin(lts) ) );
220     r = earth*abs(alpha');
221    
222     return
223    

  ViewVC Help
Powered by ViewVC 1.1.22