/[MITgcm]/MITgcm/verification/tutorial_global_oce_latlon/diags_matlab/ncep2global_ocean.m
ViewVC logotype

Contents of /MITgcm/verification/tutorial_global_oce_latlon/diags_matlab/ncep2global_ocean.m

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


Revision 1.3 - (show annotations) (download)
Sat Aug 12 20:25:13 2006 UTC (17 years, 8 months ago) by jmc
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, checkpoint58u_post, checkpoint58w_post, checkpoint63p, checkpoint63q, checkpoint63r, checkpoint63s, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63g, checkpoint63a, checkpoint63b, checkpoint63c, checkpoint64, checkpoint65, checkpoint60, checkpoint61, checkpoint62, checkpoint63, checkpoint58r_post, checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint58n_post, checkpoint58x_post, checkpoint58t_post, checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65n, checkpoint65o, checkpoint65l, checkpoint65m, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65f, checkpoint65g, checkpoint65d, checkpoint65e, checkpoint58q_post, checkpoint59q, checkpoint59p, checkpoint59r, checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59a, checkpoint59c, checkpoint59b, checkpoint59m, checkpoint59l, checkpoint59o, checkpoint59n, checkpoint59i, checkpoint59h, checkpoint59k, checkpoint59j, checkpoint59, checkpoint58o_post, checkpoint62c, checkpoint62b, checkpoint62a, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62w, checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62z, checkpoint62y, checkpoint62x, checkpoint58y_post, checkpoint58v_post, checkpoint58s_post, checkpoint61f, checkpoint61g, checkpoint61d, checkpoint61e, checkpoint61b, checkpoint61c, checkpoint58p_post, checkpoint61a, checkpoint61n, checkpoint61o, checkpoint61l, checkpoint61m, checkpoint61j, checkpoint61k, checkpoint61h, checkpoint61i, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p, checkpoint61q, checkpoint61z, checkpoint61x, checkpoint61y, HEAD
Changes since 1.2: +1 -1 lines
accidentally removed ; put them back.

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

  ViewVC Help
Powered by ViewVC 1.1.22