/[MITgcm]/MITgcm_contrib/gael/profilesMatlabProcessing/profiles_IO_v2/MITprof_create.m
ViewVC logotype

Annotation of /MITgcm_contrib/gael/profilesMatlabProcessing/profiles_IO_v2/MITprof_create.m

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


Revision 1.6 - (hide annotations) (download)
Mon Apr 11 20:56:08 2011 UTC (14 years, 3 months ago) by roquet
Branch: MAIN
Changes since 1.5: +36 -80 lines
use new interface for old/new netcdf toolbox

1 gforget 1.1 %function: MITprof_create
2     %object: read netcdf data files in the "MIT format"
3     %author: Gael Forget (gforget@mit.edu)
4     %date: june 21th, 2006
5     %
6 gforget 1.5 %usage: [MITprof]=MITprof_create(fileOut,nProf,nLev);
7 gforget 1.1 % ---> create an empty MITprof netcdf file that contains all T/S fields
8 gforget 1.5 % [MITprof]=MITprof_create(fileOut,nProf,nLev,list_vars);
9 gforget 1.1 % ---> same but specifying the list of variables in list_vars cell
10     % array (e.g. list_vars={'prof_T','prof_Tweight'}).
11     %
12 gforget 1.4 %inputs: fileOut data file name
13     % nProf number of profiles
14     % nLev number of levels
15     % list_vars variable list (optional)
16 gforget 1.1 %
17 gforget 1.4 %outputs: (none)
18 gforget 1.1
19 gforget 1.4 function []=MITprof_create(fileOut,nProf,nLev,varargin);
20 gforget 1.1
21     %=============list of variables that will actually be in the file==============%
22 roquet 1.6 if nargin>3;
23     list_vars=varargin{1};
24     else;
25     list_vars={'prof_T','prof_Tweight','prof_Testim','prof_Terr','prof_Tflag',...
26 gforget 1.5 'prof_S','prof_Sweight','prof_Sestim','prof_Serr','prof_Sflag'};
27 gforget 1.1 end;
28    
29     list_vars_plus=[{'prof_depth','prof_date','prof_YYYYMMDD','prof_HHMMSS',...
30 gforget 1.5 'prof_lon','prof_lat','prof_basin','prof_point'}...
31     list_vars];
32 gforget 1.1
33     %==========masters table of variables, units, names and dimensions=============%
34    
35     mt_v={'prof_depth'}; mt_u={'me'}; mt_n={'depth'}; mt_d={'iDEPTH'};
36     %mt_v=[mt_v '']; mt_u=[mt_u ' ']; mt_n=[mt_n '']; mt_d=[mt_d ''];
37     mt_v=[mt_v 'prof_date']; mt_u=[mt_u ' ']; mt_n=[mt_n 'Julian day since Jan-1-0000']; mt_d=[mt_d 'iPROF'];
38     mt_v=[mt_v 'prof_YYYYMMDD']; mt_u=[mt_u ' ']; mt_n=[mt_n 'year (4 digits), month (2 digits), day (2 digits)']; mt_d=[mt_d 'iPROF'];
39     mt_v=[mt_v 'prof_HHMMSS']; mt_u=[mt_u ' ']; mt_n=[mt_n 'hour (2 digits), minute (2 digits), second (2 digits)']; mt_d=[mt_d 'iPROF'];
40     mt_v=[mt_v 'prof_lon']; mt_u=[mt_u ' ']; mt_n=[mt_n 'Longitude (degree East)']; mt_d=[mt_d 'iPROF'];
41     mt_v=[mt_v 'prof_lat']; mt_u=[mt_u ' ']; mt_n=[mt_n 'Latitude (degree North)']; mt_d=[mt_d 'iPROF'];
42     mt_v=[mt_v 'prof_basin']; mt_u=[mt_u ' ']; mt_n=[mt_n 'ocean basin index (ecco 4g)']; mt_d=[mt_d 'iPROF'];
43     mt_v=[mt_v 'prof_point']; mt_u=[mt_u ' ']; mt_n=[mt_n 'grid point index (ecco 4g)']; mt_d=[mt_d 'iPROF'];
44     %
45     mt_v=[mt_v 'prof_T']; mt_u=[mt_u 'degree C']; mt_n=[mt_n 'potential temperature']; mt_d=[mt_d 'iPROF,iDEPTH'];
46     mt_v=[mt_v 'prof_Tweight']; mt_u=[mt_u '(degree C)^-2']; mt_n=[mt_n 'pot. temp. least-square weight']; mt_d=[mt_d 'iPROF,iDEPTH'];
47 gforget 1.2 mt_v=[mt_v 'prof_Testim']; mt_u=[mt_u 'degree C']; mt_n=[mt_n 'pot. temp. estimate (e.g. from atlas)']; mt_d=[mt_d 'iPROF,iDEPTH'];
48 gforget 1.1 mt_v=[mt_v 'prof_Terr']; mt_u=[mt_u 'degree C']; mt_n=[mt_n 'pot. temp. instrumental error']; mt_d=[mt_d 'iPROF,iDEPTH'];
49 gforget 1.2 mt_v=[mt_v 'prof_Tflag']; mt_u=[mt_u ' ']; mt_n=[mt_n 'flag = i > 0 means test i rejected data.']; mt_d=[mt_d 'iPROF,iDEPTH'];
50 gforget 1.1 %
51     mt_v=[mt_v 'prof_S']; mt_u=[mt_u 'psu']; mt_n=[mt_n 'salinity']; mt_d=[mt_d 'iPROF,iDEPTH'];
52     mt_v=[mt_v 'prof_Sweight']; mt_u=[mt_u '(psu)^-2']; mt_n=[mt_n 'salinity least-square weight']; mt_d=[mt_d 'iPROF,iDEPTH'];
53 gforget 1.2 mt_v=[mt_v 'prof_Sestim']; mt_u=[mt_u 'psu']; mt_n=[mt_n 'salinity estimate (e.g. from atlas)']; mt_d=[mt_d 'iPROF,iDEPTH'];
54 gforget 1.1 mt_v=[mt_v 'prof_Serr']; mt_u=[mt_u 'psu']; mt_n=[mt_n 'salinity instrumental error']; mt_d=[mt_d 'iPROF,iDEPTH'];
55 gforget 1.2 mt_v=[mt_v 'prof_Sflag']; mt_u=[mt_u ' ']; mt_n=[mt_n 'flag = i > 0 means test i rejected data.']; mt_d=[mt_d 'iPROF,iDEPTH'];
56 gforget 1.1 %
57     mt_v=[mt_v 'prof_U']; mt_u=[mt_u 'm/s']; mt_n=[mt_n 'eastward velocity comp.']; mt_d=[mt_d 'iPROF,iDEPTH'];
58     mt_v=[mt_v 'prof_Uweight']; mt_u=[mt_u '(m/s)^-2']; mt_n=[mt_n 'east. v. least-square weight']; mt_d=[mt_d 'iPROF,iDEPTH'];
59     mt_v=[mt_v 'prof_V']; mt_u=[mt_u 'm/s']; mt_n=[mt_n 'northward velocity comp.']; mt_d=[mt_d 'iPROF,iDEPTH'];
60     mt_v=[mt_v 'prof_Vweight']; mt_u=[mt_u '(m/s)^-2']; mt_n=[mt_n 'north. v. least-square weight']; mt_d=[mt_d 'iPROF,iDEPTH'];
61     mt_v=[mt_v 'prof_ptr']; mt_u=[mt_u 'X']; mt_n=[mt_n 'passive tracer']; mt_d=[mt_d 'iPROF,iDEPTH'];
62     mt_v=[mt_v 'prof_ptrweight']; mt_u=[mt_u '(X)^-2']; mt_n=[mt_n 'pass. tracer least-square weight']; mt_d=[mt_d 'iPROF,iDEPTH'];
63     %
64     mt_v=[mt_v 'prof_bp']; mt_u=[mt_u 'cm']; mt_n=[mt_n 'bottom pressure']; mt_d=[mt_d 'iPROF'];
65     mt_v=[mt_v 'prof_bpweight']; mt_u=[mt_u '(cm)^-2']; mt_n=[mt_n 'bot. pres. least-square weight']; mt_d=[mt_d 'iPROF'];
66     mt_v=[mt_v 'prof_ssh']; mt_u=[mt_u 'cm']; mt_n=[mt_n 'sea surface height']; mt_d=[mt_d 'iPROF'];
67     mt_v=[mt_v 'prof_sshweight']; mt_u=[mt_u '(cm)^-2']; mt_n=[mt_n 'ssh least-square weight']; mt_d=[mt_d 'iPROF'];
68    
69     %=============================create the file=================================%
70    
71     if ~strcmp(fileOut(end-2:end),'.nc'); suff='.nc'; else; suff=''; end;
72 gforget 1.3
73 gforget 1.1 %define netcdf dimensions :
74 gforget 1.4 iPROF=nProf; iDEPTH=nLev;
75 roquet 1.6 lTXT=30; fillval=double(-9999);
76 gforget 1.1
77 gforget 1.5
78    
79 roquet 1.6
80     % write the netcdf structure
81     ncid=nccreate([fileOut suff],'clobber');
82    
83     aa=sprintf(['Format: MITprof. This file was created using \n' ...
84     'the matlab toolbox which can be obtained (see README) from \n'...
85     'http://mitgcm.org/viewvc/MITgcm/MITgcm_contrib/gael/profilesMatlabProcessing/']);
86     ncputAtt(ncid,'','description',aa);
87     ncputAtt(ncid,'','date',date);
88    
89     ncdefDim(ncid,'iPROF',iPROF);
90     ncdefDim(ncid,'iDEPTH',iDEPTH);
91     ncdefDim(ncid,'lTXT',lTXT);
92    
93     ncdefVar(ncid,'prof_descr','char',{'lTXT','iPROF'});
94     ncputAtt(ncid,'prof_descr','long_name','profile description');
95    
96     for ii=1:length(list_vars_plus);
97     jj=find(strcmp(mt_v,list_vars_plus{ii}));
98     if ~isempty(jj);
99     if strcmp(mt_d{jj},'iPROF,iDEPTH');
100     ncdefVar(ncid,mt_v{jj},'double',{'iDEPTH','iPROF'});%note the direction flip
101 gforget 1.5 else;
102 roquet 1.6 ncdefVar(ncid,mt_v{jj},'double',{mt_d{jj}});
103 gforget 1.5 end;
104 roquet 1.6 ncputAtt(ncid,mt_v{jj},'long_name',mt_n{jj});
105     ncputAtt(ncid,mt_v{jj},'units',mt_u{jj});
106     ncputAtt(ncid,mt_v{jj},'missing_value',fillval);
107     ncputAtt(ncid,mt_v{jj},'_FillValue',fillval);
108     else;
109     warning([list_vars_plus{ii} ' not included -- it is not a MITprof variable']);
110 gforget 1.1 end;
111     end;
112    
113 roquet 1.6 ncclose(ncid);
114    
115 gforget 1.5
116 gforget 1.1

  ViewVC Help
Powered by ViewVC 1.1.22