/[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.3 - (hide annotations) (download)
Wed Nov 3 19:59:56 2010 UTC (14 years, 8 months ago) by gforget
Branch: MAIN
Changes since 1.2: +7 -0 lines
- take care of some print statements.
- introduce profiles_read_seals.m

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     %usage: []=MITprof_create(dirOut,fileOut,specs);
7     % ---> create an empty MITprof netcdf file that contains all T/S fields
8     % []=MITprof_create(dirOut,fileOut,specs,list_vars);
9     % ---> same but specifying the list of variables in list_vars cell
10     % array (e.g. list_vars={'prof_T','prof_Tweight'}).
11     %
12     %inputs: dirOut data file directory
13     % fileOut data file name
14     % specs {iProf,iLev,lTXT,fillval}
15     % nProf is the number of profiles
16     % nLev is the number of levels,
17     % nChar is the length of profiles description
18     % (list_vars optional variable list)
19     %
20     %outputs: (none)
21    
22     function []=MITprof_create(dirOut,fileOut,specs,varargin);
23    
24     %=============list of variables that will actually be in the file==============%
25     if nargin>4; list_vars=varargin{1};
26 gforget 1.2 else; list_vars={'prof_T','prof_Tweight','prof_Testim','prof_Terr','prof_Tflag',...
27     'prof_S','prof_Sweight','prof_Sestim','prof_Serr','prof_Sflag'};
28 gforget 1.1 end;
29    
30     list_vars_plus=[{'prof_depth','prof_date','prof_YYYYMMDD','prof_HHMMSS',...
31     'prof_lon','prof_lat','prof_basin','prof_point'}...
32     list_vars];
33    
34     %==========masters table of variables, units, names and dimensions=============%
35    
36     mt_v={'prof_depth'}; mt_u={'me'}; mt_n={'depth'}; mt_d={'iDEPTH'};
37     %mt_v=[mt_v '']; mt_u=[mt_u ' ']; mt_n=[mt_n '']; mt_d=[mt_d ''];
38     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'];
39     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'];
40     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'];
41     mt_v=[mt_v 'prof_lon']; mt_u=[mt_u ' ']; mt_n=[mt_n 'Longitude (degree East)']; mt_d=[mt_d 'iPROF'];
42     mt_v=[mt_v 'prof_lat']; mt_u=[mt_u ' ']; mt_n=[mt_n 'Latitude (degree North)']; mt_d=[mt_d 'iPROF'];
43     mt_v=[mt_v 'prof_basin']; mt_u=[mt_u ' ']; mt_n=[mt_n 'ocean basin index (ecco 4g)']; mt_d=[mt_d 'iPROF'];
44     mt_v=[mt_v 'prof_point']; mt_u=[mt_u ' ']; mt_n=[mt_n 'grid point index (ecco 4g)']; mt_d=[mt_d 'iPROF'];
45     %
46     mt_v=[mt_v 'prof_T']; mt_u=[mt_u 'degree C']; mt_n=[mt_n 'potential temperature']; mt_d=[mt_d 'iPROF,iDEPTH'];
47     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'];
48 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'];
49 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'];
50 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'];
51 gforget 1.1 %
52     mt_v=[mt_v 'prof_S']; mt_u=[mt_u 'psu']; mt_n=[mt_n 'salinity']; mt_d=[mt_d 'iPROF,iDEPTH'];
53     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'];
54 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'];
55 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'];
56 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'];
57 gforget 1.1 %
58     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'];
59     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'];
60     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'];
61     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'];
62     mt_v=[mt_v 'prof_ptr']; mt_u=[mt_u 'X']; mt_n=[mt_n 'passive tracer']; mt_d=[mt_d 'iPROF,iDEPTH'];
63     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'];
64     %
65     mt_v=[mt_v 'prof_bp']; mt_u=[mt_u 'cm']; mt_n=[mt_n 'bottom pressure']; mt_d=[mt_d 'iPROF'];
66     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'];
67     mt_v=[mt_v 'prof_ssh']; mt_u=[mt_u 'cm']; mt_n=[mt_n 'sea surface height']; mt_d=[mt_d 'iPROF'];
68     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'];
69    
70     %=============================create the file=================================%
71    
72     if ~strcmp(fileOut(end-2:end),'.nc'); suff='.nc'; else; suff=''; end;
73     nc=netcdf([dirOut fileOut suff],'clobber');
74    
75 gforget 1.3 %put a couple global attributes
76     aa=sprintf(['Format: MITprof. This file was created using \n' ...
77     'the matlab toolbox which can be obtained (see README) from \n'...
78     'http://mitgcm.org/viewvc/MITgcm/MITgcm_contrib/gael/profilesMatlabProcessing/']);
79     nc.description=aa;
80     nc.date = date;
81    
82 gforget 1.1 %define netcdf dimensions :
83     iPROF=specs{1}; iDEPTH=specs{2};
84     lTXT=specs{3}; fillval=specs{4};
85    
86     nc('iPROF') = iPROF;
87     nc('iDEPTH') = iDEPTH;
88     nc('lTXT') = lTXT;
89    
90     %mt_v=[mt_v 'prof_descr']; mt_u=[mt_u ' ']; mt_n=[mt_n 'profile description']; mt_d=[mt_d 'iPROF,lTXT'];
91     nc{'prof_descr'} = ncchar('iPROF','lTXT');
92     nc{'prof_descr'}.long_name='profile description';
93    
94     for ii=1:length(list_vars_plus);
95     jj=find(strcmp(mt_v,list_vars_plus{ii}));
96     if ~isempty(jj);
97     if strcmp(mt_d{jj},'iPROF,iDEPTH');
98     nc{mt_v{jj}} = ncdouble('iPROF','iDEPTH');
99     else;
100     nc{mt_v{jj}} = mt_d{jj};
101     end;
102     nc{mt_v{jj}}.long_name = mt_n{jj};
103     nc{mt_v{jj}}.units = mt_u{jj};
104     nc{mt_v{jj}}.missing_value=fillval;
105     end;
106     end;
107    
108     close(nc);
109    

  ViewVC Help
Powered by ViewVC 1.1.22