/[MITgcm]/MITgcm/utils/matlab/rdmnc.m
ViewVC logotype

Diff of /MITgcm/utils/matlab/rdmnc.m

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

revision 1.21 by mlosch, Sun May 16 08:56:43 2010 UTC revision 1.22 by mlosch, Fri May 28 07:16:43 2010 UTC
# Line 299  function [S] = rdmnc_local_matlabAPI(fna Line 299  function [S] = rdmnc_local_matlabAPI(fna
299    % now open the file for reading    % now open the file for reading
300    nc = netcdf.open(fname,'NC_NOWRITE');    nc = netcdf.open(fname,'NC_NOWRITE');
301    % get basic information about netcdf file    % get basic information about netcdf file
302    [ndims nvars natts dimm] = netcdf.inq(nc);    [ndims nvars natts recdim] = netcdf.inq(nc);
303    
304    % No variables specified? Default to all    % No variables specified? Default to all
305    if isempty(varlist),    if isempty(varlist),
# Line 336  function [S] = rdmnc_local_matlabAPI(fna Line 336  function [S] = rdmnc_local_matlabAPI(fna
336      end      end
337            
338      [varname,xtype,dimids,natts] = netcdf.inqVar(nc,varid);      [varname,xtype,dimids,natts] = netcdf.inqVar(nc,varid);
339        % does this variable contain a record (unlimited) dimension?
340        [isrecvar,recpos] = ismember(recdim,dimids);
341            
342      % Dimensions      % Dimensions
343      clear sizVar dims      clear sizVar dims
344      for k=1:length(dimids)      for k=1:length(dimids)
345        [dims{k}, sizVar(k)] = netcdf.inqDim(nc,dimids(k));        [dims{k}, sizVar(k)] = netcdf.inqDim(nc,dimids(k));
346      end      end
     if length(sizVar)==1; sizVar(2) = 1; end;  
347      nDims=length(sizVar);      nDims=length(sizVar);
348      if dims{1} == 'T'      if isrecvar
349        if isempty(find(fii)), error('Iters not found'); end        if isempty(find(fii)), error('Iters not found'); end
350        it = length(dims);        it = length(dims);
351        tmpdata = netcdf.getVar(nc,varid,'double');        if length(dimids) == 1
352        tmpdata = tmpdata(fii,:);          % just a time or iteration variable, this will result in a vector
353        % leading unity dimensions get lost; add them back:          % and requires special treatment
354        tmpdata=reshape(tmpdata,[length(fii) sizVar(2:end)]);          icount=1;
355            tmpdata = zeros(length(fii),1);
356            for k=1:length(fii)
357              istart = fii(k)-1;
358              tmpdata(k) = netcdf.getVar(nc,varid,istart,icount,'double');
359            end
360          else
361            % from now on we assume that the record dimension is always the
362            % last dimension. This may not always be the case
363            if recpos ~= nDims
364              error(sprintf('%s\n%s\n%s%s%i%s%i', ...
365                            ['The current code assumes that the record ' ...
366                             'dimension is the last dimension,'], ...
367                            'this is not the case for variable', cvar, ...
368                            ': nDims = ', nDims,  ...
369                            ', position of recDim = ', recpos))
370            end
371            istart = zeros(1,it); % indexing starts a 0
372            icount = sizVar;
373            % we always want to get only on time slice at a time
374            icount(recpos) = 1;
375            % make your life simpler by putting the time dimension first
376            tmpdata = zeros([length(fii) sizVar(1:end-1)]);
377            for k=1:length(fii)
378              istart(recpos) = fii(k)-1; % indexing starts at 0
379              tmp = netcdf.getVar(nc,varid,istart,icount,'double');
380              tmpdata(k,:) = tmp(:);
381            end
382            % move time dimension to the end ...
383            tmpdata = shiftdim(tmpdata,1);
384            % ... and restore original shape
385            tmpdata = reshape(tmpdata,[sizVar(1:end-1) length(fii)]);
386          end
387      else      else
388        it = 0;        it = 0;
389        tmpdata = netcdf.getVar(nc,varid,'double');        tmpdata = netcdf.getVar(nc,varid,'double');
       % leading unity dimensions get lost; add them back:  
       tmpdata=reshape(tmpdata,sizVar);  
390      end      end
391            %
392      if dBug > 1,      if dBug > 1,
393        fprintf(['  var:',cvar,': nDims=%i ('],nDims);fprintf(' %i',sizVar);        fprintf(['  var:',cvar,': nDims=%i ('],nDims);fprintf(' %i',sizVar);
394        fprintf('):%i,nD=%i,it=%i ;',length(size(tmpdata)),length(dims),it);        fprintf('):%i,nD=%i,it=%i ;',length(size(tmpdata)),length(dims),it);
395      end      end
396      [ni nj nk nm nn no np]=size(tmpdata);      [ni nj nk nm nn no np]=size(tmpdata);
397              %
398      [i0,j0,fn]=findTileOffset(S);      [i0,j0,fn]=findTileOffset(S);
399      cdim=dims{1}; if cdim(1)~='X'; i0=0; end      cdim=dims{1}; if cdim(1)~='X'; i0=0; end
400      cdim=dims{1}; if cdim(1)=='Y'; i0=j0; j0=0; end      cdim=dims{1}; if cdim(1)=='Y'; i0=j0; j0=0; end
# Line 375  function [S] = rdmnc_local_matlabAPI(fna Line 406  function [S] = rdmnc_local_matlabAPI(fna
406      if dBug > 1,      if dBug > 1,
407        fprintf(' i0,ni= %i %i; j,nj= %i %i; nk=%i :',i0,ni,j0,nj,nk);        fprintf(' i0,ni= %i %i; j,nj= %i %i; nk=%i :',i0,ni,j0,nj,nk);
408      end      end
409            %
410      Sstr = '';      Sstr = '';
411      for istr = 1:max(nDims,length(dims)),      for istr = 1:max(nDims,length(dims)),
412        if     istr == it,  Sstr = [Sstr,'dii,'];        if     istr == it,  Sstr = [Sstr,'dii,'];
# Line 392  function [S] = rdmnc_local_matlabAPI(fna Line 423  function [S] = rdmnc_local_matlabAPI(fna
423      eval(['S.(cvar)(',Sstr(1:end-1),')=tmpdata;'])      eval(['S.(cvar)(',Sstr(1:end-1),')=tmpdata;'])
424      %S.(cvar)(i0+(1:ni),j0+(1:nj),(1:nk),(1:nm),(1:nn),(1:no),(1:np))=tmpdata;      %S.(cvar)(i0+(1:ni),j0+(1:nj),(1:nk),(1:nm),(1:nn),(1:no),(1:np))=tmpdata;
425      if dBug > 1, fprintf(' %i',size(S.(cvar))); fprintf('\n'); end      if dBug > 1, fprintf(' %i',size(S.(cvar))); fprintf('\n'); end
426            %
427      S.attributes.(cvar)=ncgetatt(nc,cvar);      S.attributes.(cvar)=ncgetatt(nc,cvar);
428      % replace missing or FillValues with NaN      % replace missing or FillValues with NaN
429      if isfield(S.attributes.(cvar),'missing_value');      if isfield(S.attributes.(cvar),'missing_value');

Legend:
Removed from v.1.21  
changed lines
  Added in v.1.22

  ViewVC Help
Powered by ViewVC 1.1.22