/[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.13 - (hide annotations) (download)
Mon Feb 1 14:35:06 2016 UTC (9 years, 5 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint65t
Changes since 1.12: +12 -24 lines
- MITprof_create.m: replace nProf,prof_depth inputs by MITprofCur;
  determine iINTERP based upon MITprofCur.prof_interp_weights.
- MITprof_write.m: update call to MITprof_create.

1 gforget 1.13 function []=MITprof_create(fileOut,MITprofCur,varargin);
2     % MITPROF_CREATE(fileOut,MITprofCur,list_vars)
3     % Creates a file in the MITprof format. The list of variables
4     % can further be specified via list_vars (optional).
5 roquet 1.7
6     % check that file exists and add prefix and suffix if necessary
7     [pathstr, name, ext] = fileparts(fileOut);
8     if isempty(pathstr) | strcmp(pathstr,'.'), pathstr=pwd; end
9     if isempty(ext) | ~strcmp(ext,'.nc'), ext='.nc'; end
10     fileOut=[pathstr '/' name ext];
11    
12     %define netcdf dimensions :
13 gforget 1.13 prof_depth=MITprofCur.prof_depth;
14 roquet 1.7 prof_depth=reshape(prof_depth,length(prof_depth),1);
15 gforget 1.13 iPROF=MITprofCur.np; iDEPTH=MITprofCur.nr;
16     if isfield(MITprofCur,'prof_interp_weights');
17     iINTERP=size(MITprofCur.prof_interp_weights,2);
18     else;
19     iINTERP=1;
20     end;
21 roquet 1.7 lTXT=30; fillval=double(-9999);
22 gforget 1.1
23     %=============list of variables that will actually be in the file==============%
24 gforget 1.13 if nargin>2;
25 gforget 1.11 list_vars=varargin{1};
26     else;
27     list_vars={'prof_depth','prof_descr','prof_date','prof_YYYYMMDD','prof_HHMMSS',...
28     'prof_lon','prof_lat','prof_basin','prof_point','prof_flag','prof_T','prof_Tweight','prof_Testim','prof_Terr','prof_Tflag',...
29     'prof_S','prof_Sweight','prof_Sestim','prof_Serr','prof_Sflag','prof_D','prof_Destim'};
30     end;
31    
32 gforget 1.1
33 roquet 1.7 % eliminate doublons
34 gforget 1.11 [list,m]=unique(list_vars);
35     list_vars=list_vars(sort(m));
36    
37     list_vars_plus={};
38     for ii=1:length(list_vars);
39     if ~isempty(strfind(list_vars{ii},'estim'))&...
40     ~strcmp(list_vars{ii},'prof_Testim')&...
41     ~strcmp(list_vars{ii},'prof_Sestim');
42     list_vars_plus={list_vars_plus{:},list_vars{ii}(1:end-5)};
43     end;
44     end;
45 roquet 1.7
46 gforget 1.1 %==========masters table of variables, units, names and dimensions=============%
47    
48     mt_v={'prof_depth'}; mt_u={'me'}; mt_n={'depth'}; mt_d={'iDEPTH'};
49     %mt_v=[mt_v '']; mt_u=[mt_u ' ']; mt_n=[mt_n '']; mt_d=[mt_d ''];
50     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'];
51     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'];
52     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'];
53     mt_v=[mt_v 'prof_lon']; mt_u=[mt_u ' ']; mt_n=[mt_n 'Longitude (degree East)']; mt_d=[mt_d 'iPROF'];
54     mt_v=[mt_v 'prof_lat']; mt_u=[mt_u ' ']; mt_n=[mt_n 'Latitude (degree North)']; mt_d=[mt_d 'iPROF'];
55     mt_v=[mt_v 'prof_basin']; mt_u=[mt_u ' ']; mt_n=[mt_n 'ocean basin index (ecco 4g)']; mt_d=[mt_d 'iPROF'];
56     mt_v=[mt_v 'prof_point']; mt_u=[mt_u ' ']; mt_n=[mt_n 'grid point index (ecco 4g)']; mt_d=[mt_d 'iPROF'];
57 gforget 1.9 mt_v=[mt_v 'prof_flag']; mt_u=[mt_u ' ']; mt_n=[mt_n 'flag = i > 0 for suspicious profile ']; mt_d=[mt_d 'iPROF'];
58 gforget 1.1 %
59     mt_v=[mt_v 'prof_T']; mt_u=[mt_u 'degree C']; mt_n=[mt_n 'potential temperature']; mt_d=[mt_d 'iPROF,iDEPTH'];
60     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'];
61 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'];
62 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'];
63 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'];
64 gforget 1.1 %
65     mt_v=[mt_v 'prof_S']; mt_u=[mt_u 'psu']; mt_n=[mt_n 'salinity']; mt_d=[mt_d 'iPROF,iDEPTH'];
66     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'];
67 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'];
68 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'];
69 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'];
70 gforget 1.1 %
71 gforget 1.11 for ii=1:length(list_vars_plus);
72     mt_v=[mt_v list_vars_plus{ii}]; mt_u=[mt_u 'unknown']; mt_n=[mt_n 'unknown']; mt_d=[mt_d 'iPROF,iDEPTH'];
73     mt_v=[mt_v [list_vars_plus{ii} 'weight']]; mt_u=[mt_u '(unknown)^-2']; mt_n=[mt_n 'unknown least-square weight']; mt_d=[mt_d 'iPROF,iDEPTH'];
74     mt_v=[mt_v [list_vars_plus{ii} 'estim']]; mt_u=[mt_u 'unknown']; mt_n=[mt_n 'unknown estimate (e.g. from atlas)']; mt_d=[mt_d 'iPROF,iDEPTH'];
75     mt_v=[mt_v [list_vars_plus{ii} 'err']]; mt_u=[mt_u 'unknown']; mt_n=[mt_n 'unknown instrumental error']; mt_d=[mt_d 'iPROF,iDEPTH'];
76     mt_v=[mt_v [list_vars_plus{ii} 'flag']]; mt_u=[mt_u ' ']; mt_n=[mt_n 'flag = i > 0 means test i rejected data.']; mt_d=[mt_d 'iPROF,iDEPTH'];
77     end;
78 gforget 1.12 %
79     list_interpr={'prof_interp_XC11','prof_interp_YC11','prof_interp_XCNINJ','prof_interp_YCNINJ'};
80     for ii=1:length(list_interpr);
81     mt_v=[mt_v list_interpr{ii}]; mt_u=[mt_u 'unknown']; mt_n=[mt_n 'interpolation variable']; mt_d=[mt_d 'iPROF'];
82     end;
83     list_interpr={'prof_interp_i','prof_interp_j','prof_interp_lon','prof_interp_lat','prof_interp_weights'};
84     for ii=1:length(list_interpr);
85     mt_v=[mt_v list_interpr{ii}]; mt_u=[mt_u 'unknown']; mt_n=[mt_n 'interpolation variable']; mt_d=[mt_d 'iPROF,iINTERP'];
86     end;
87     %
88 gforget 1.11 if 0;
89     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'];
90     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'];
91     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'];
92     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'];
93     mt_v=[mt_v 'prof_ptr']; mt_u=[mt_u 'X']; mt_n=[mt_n 'passive tracer']; mt_d=[mt_d 'iPROF,iDEPTH'];
94     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'];
95     %
96     mt_v=[mt_v 'prof_D']; mt_u=[mt_u 'me']; mt_n=[mt_n 'variable depth']; mt_d=[mt_d 'iPROF,iDEPTH'];
97     mt_v=[mt_v 'prof_Destim']; mt_u=[mt_u 'me']; mt_n=[mt_n 'variable depth estimate (e.g. from atlas)']; mt_d=[mt_d 'iPROF,iDEPTH'];
98     %
99     mt_v=[mt_v 'prof_bp']; mt_u=[mt_u 'cm']; mt_n=[mt_n 'bottom pressure']; mt_d=[mt_d 'iPROF'];
100     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'];
101     mt_v=[mt_v 'prof_ssh']; mt_u=[mt_u 'cm']; mt_n=[mt_n 'sea surface height']; mt_d=[mt_d 'iPROF'];
102     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'];
103     end;
104 gforget 1.1
105     %=============================create the file=================================%
106    
107 roquet 1.6 % write the netcdf structure
108 roquet 1.7 ncid=nccreate(fileOut,'clobber');
109 roquet 1.6
110 roquet 1.7 aa=sprintf(['Format: MITprof netcdf. This file was created using \n' ...
111 roquet 1.6 'the matlab toolbox which can be obtained (see README) from \n'...
112     'http://mitgcm.org/viewvc/MITgcm/MITgcm_contrib/gael/profilesMatlabProcessing/']);
113     ncputAtt(ncid,'','description',aa);
114     ncputAtt(ncid,'','date',date);
115    
116     ncdefDim(ncid,'iPROF',iPROF);
117     ncdefDim(ncid,'iDEPTH',iDEPTH);
118 gforget 1.12 ncdefDim(ncid,'iINTERP',iINTERP);
119 roquet 1.6 ncdefDim(ncid,'lTXT',lTXT);
120    
121 gforget 1.11 for ii=1:length(list_vars);
122     jj=find(strcmp(mt_v,list_vars{ii}));
123 roquet 1.6 if ~isempty(jj);
124     if strcmp(mt_d{jj},'iPROF,iDEPTH');
125     ncdefVar(ncid,mt_v{jj},'double',{'iDEPTH','iPROF'});%note the direction flip
126 gforget 1.12 elseif strcmp(mt_d{jj},'iPROF,iINTERP');
127     ncdefVar(ncid,mt_v{jj},'double',{'iINTERP','iPROF'});%note the direction flip
128 gforget 1.5 else;
129 roquet 1.6 ncdefVar(ncid,mt_v{jj},'double',{mt_d{jj}});
130 gforget 1.5 end;
131 roquet 1.6 ncputAtt(ncid,mt_v{jj},'long_name',mt_n{jj});
132     ncputAtt(ncid,mt_v{jj},'units',mt_u{jj});
133     ncputAtt(ncid,mt_v{jj},'missing_value',fillval);
134     ncputAtt(ncid,mt_v{jj},'_FillValue',fillval);
135     else;
136 gforget 1.11 if strcmp(list_vars{ii},'prof_descr')
137 roquet 1.8 ncdefVar(ncid,'prof_descr','char',{'lTXT','iPROF'});
138     ncputAtt(ncid,'prof_descr','long_name','profile description');
139     else
140 gforget 1.11 warning([list_vars{ii} ' not included -- it is not a MITprof variable']);
141 roquet 1.8 end
142 gforget 1.1 end;
143     end;
144    
145 roquet 1.6 ncclose(ncid);
146    
147 roquet 1.7 %=============================set prof_depth=================================%
148    
149     ncid=ncopen(fileOut,'write');
150     ncputvar(ncid,'prof_depth',prof_depth);
151     ncclose(ncid);
152 gforget 1.5
153 gforget 1.1

  ViewVC Help
Powered by ViewVC 1.1.22