/[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.7 - (hide annotations) (download)
Thu Apr 14 20:12:07 2011 UTC (14 years, 3 months ago) by roquet
Branch: MAIN
Changes since 1.6: +39 -28 lines
- Simplification of IO interface to MITcdf netcdf format.
  MITprof_create and MITprof_read are now low-level functions that should not be used by user.
  Instead, use MITprof_write to create/update a file, and MITprof_load to load data.
- Depth levels must now be provided as an argument of MITprof_create and MITprof_struct,
  and cannot be changed once the MITprof ncfile or struct have been created.
- profiles_prep_main and profiles_prep_write_nc have been changed accordingly.

1 roquet 1.7 function []=MITprof_create(fileOut,nProf,prof_depth,varargin)
2 gforget 1.1 %function: MITprof_create
3 roquet 1.7 %object: create a file in the "MIT format". Low-level function.
4 gforget 1.1 %author: Gael Forget (gforget@mit.edu)
5     %date: june 21th, 2006
6     %
7 roquet 1.7 %usage: [MITprof]=MITprof_create(fileOut,nProf,prof_depth);
8     % create an empty MITprof netcdf file
9     % vertical depth levels are set according to prof_depth
10     % create empty variables for all usual T/S fields
11     %
12     % [MITprof]=MITprof_create(fileOut,nProf,prof_depth,list_vars);
13     % same but specifying the list of variables in list_vars cell
14 gforget 1.1 % array (e.g. list_vars={'prof_T','prof_Tweight'}).
15     %
16 gforget 1.4 %inputs: fileOut data file name
17 roquet 1.7 % nProf number of profiles
18     % prof_depth vector of depth levels
19     % list_vars variable list (optional)
20 gforget 1.1 %
21    
22 roquet 1.7
23     % check that file exists and add prefix and suffix if necessary
24     [pathstr, name, ext] = fileparts(fileOut);
25     if isempty(pathstr) | strcmp(pathstr,'.'), pathstr=pwd; end
26     if isempty(ext) | ~strcmp(ext,'.nc'), ext='.nc'; end
27     fileOut=[pathstr '/' name ext];
28    
29     %define netcdf dimensions :
30     nLev=length(prof_depth);
31     prof_depth=reshape(prof_depth,length(prof_depth),1);
32     iPROF=nProf; iDEPTH=nLev;
33     lTXT=30; fillval=double(-9999);
34 gforget 1.1
35     %=============list of variables that will actually be in the file==============%
36 roquet 1.7 list_vars={'prof_T','prof_Tweight','prof_Testim','prof_Terr','prof_Tflag',...
37     'prof_S','prof_Sweight','prof_Sestim','prof_Serr','prof_Sflag'};
38     if nargin>3; list_vars=varargin{1}; end
39 gforget 1.1
40 roquet 1.7 list_vars_plus=[{'prof_depth','prof_descr','prof_date','prof_YYYYMMDD','prof_HHMMSS',...
41 gforget 1.5 'prof_lon','prof_lat','prof_basin','prof_point'}...
42     list_vars];
43 gforget 1.1
44 roquet 1.7 % eliminate doublons
45     [list,m]=unique(list_vars_plus);
46     list_vars_plus=list_vars_plus(sort(m));
47    
48 gforget 1.1 %==========masters table of variables, units, names and dimensions=============%
49    
50     mt_v={'prof_depth'}; mt_u={'me'}; mt_n={'depth'}; mt_d={'iDEPTH'};
51     %mt_v=[mt_v '']; mt_u=[mt_u ' ']; mt_n=[mt_n '']; mt_d=[mt_d ''];
52     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'];
53     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'];
54     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'];
55     mt_v=[mt_v 'prof_lon']; mt_u=[mt_u ' ']; mt_n=[mt_n 'Longitude (degree East)']; mt_d=[mt_d 'iPROF'];
56     mt_v=[mt_v 'prof_lat']; mt_u=[mt_u ' ']; mt_n=[mt_n 'Latitude (degree North)']; mt_d=[mt_d 'iPROF'];
57     mt_v=[mt_v 'prof_basin']; mt_u=[mt_u ' ']; mt_n=[mt_n 'ocean basin index (ecco 4g)']; mt_d=[mt_d 'iPROF'];
58     mt_v=[mt_v 'prof_point']; mt_u=[mt_u ' ']; mt_n=[mt_n 'grid point index (ecco 4g)']; mt_d=[mt_d 'iPROF'];
59     %
60     mt_v=[mt_v 'prof_T']; mt_u=[mt_u 'degree C']; mt_n=[mt_n 'potential temperature']; mt_d=[mt_d 'iPROF,iDEPTH'];
61     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'];
62 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'];
63 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'];
64 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'];
65 gforget 1.1 %
66     mt_v=[mt_v 'prof_S']; mt_u=[mt_u 'psu']; mt_n=[mt_n 'salinity']; mt_d=[mt_d 'iPROF,iDEPTH'];
67     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'];
68 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'];
69 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'];
70 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'];
71 gforget 1.1 %
72     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'];
73     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'];
74     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'];
75     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'];
76     mt_v=[mt_v 'prof_ptr']; mt_u=[mt_u 'X']; mt_n=[mt_n 'passive tracer']; mt_d=[mt_d 'iPROF,iDEPTH'];
77     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'];
78     %
79     mt_v=[mt_v 'prof_bp']; mt_u=[mt_u 'cm']; mt_n=[mt_n 'bottom pressure']; mt_d=[mt_d 'iPROF'];
80     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'];
81     mt_v=[mt_v 'prof_ssh']; mt_u=[mt_u 'cm']; mt_n=[mt_n 'sea surface height']; mt_d=[mt_d 'iPROF'];
82     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'];
83    
84     %=============================create the file=================================%
85    
86 roquet 1.6 % write the netcdf structure
87 roquet 1.7 ncid=nccreate(fileOut,'clobber');
88 roquet 1.6
89 roquet 1.7 aa=sprintf(['Format: MITprof netcdf. This file was created using \n' ...
90 roquet 1.6 'the matlab toolbox which can be obtained (see README) from \n'...
91     'http://mitgcm.org/viewvc/MITgcm/MITgcm_contrib/gael/profilesMatlabProcessing/']);
92     ncputAtt(ncid,'','description',aa);
93     ncputAtt(ncid,'','date',date);
94    
95     ncdefDim(ncid,'iPROF',iPROF);
96     ncdefDim(ncid,'iDEPTH',iDEPTH);
97     ncdefDim(ncid,'lTXT',lTXT);
98    
99     ncdefVar(ncid,'prof_descr','char',{'lTXT','iPROF'});
100     ncputAtt(ncid,'prof_descr','long_name','profile description');
101    
102     for ii=1:length(list_vars_plus);
103     jj=find(strcmp(mt_v,list_vars_plus{ii}));
104     if ~isempty(jj);
105     if strcmp(mt_d{jj},'iPROF,iDEPTH');
106     ncdefVar(ncid,mt_v{jj},'double',{'iDEPTH','iPROF'});%note the direction flip
107 gforget 1.5 else;
108 roquet 1.6 ncdefVar(ncid,mt_v{jj},'double',{mt_d{jj}});
109 gforget 1.5 end;
110 roquet 1.6 ncputAtt(ncid,mt_v{jj},'long_name',mt_n{jj});
111     ncputAtt(ncid,mt_v{jj},'units',mt_u{jj});
112     ncputAtt(ncid,mt_v{jj},'missing_value',fillval);
113     ncputAtt(ncid,mt_v{jj},'_FillValue',fillval);
114     else;
115     warning([list_vars_plus{ii} ' not included -- it is not a MITprof variable']);
116 gforget 1.1 end;
117     end;
118    
119 roquet 1.6 ncclose(ncid);
120    
121 roquet 1.7 %=============================set prof_depth=================================%
122    
123     ncid=ncopen(fileOut,'write');
124     ncputvar(ncid,'prof_depth',prof_depth);
125     ncclose(ncid);
126 gforget 1.5
127 gforget 1.1

  ViewVC Help
Powered by ViewVC 1.1.22