function [trace,error,varargout]=loadJETdata(shot,data_type,varargin) % % data_type: % 'Ip' = current % 'zmag' = vertical position of the center of the plasma (magnetic axis) % 'rmag' = radial position of the center of the plasma % 'rcont' = R of plama boundary vs time % 'zcont' = Z of plama boundary vs time % 'vol' = volume of flux surfaces vs rho=sqrt(psi) % 'qrho' = q profile on rho mesh % 'q95' = q95 vs time % 'kappa', 'elon' = edge elongation vs time % 'delta', 'triang' = edge averaged triangularity vs time % 'deltatop', 'triangtop' = edge upper (top) triangularity vs time % 'deltabot', 'triangbot' = edge lower (bottom) triangularity vs time % 'n1' or 'n2': n=1 or n=2 MHD signal % 'neint' = line-integrated electron density [m/m^3] % 'ne'= ne raw profile on (R,t). ADD error bars in .std % 'te'= Te raw profile on (R,t). ADD error bars in .std % 'nerho'= ne profile on (rho=sqrt(psi),time) mesh.Note rho is a 2D array as depends on time. ADD error bars in .std % 'terho'= Te profile on (rho=sqrt(psi),time) mesh.Note rho is a 2D array as depends on time. ADD error bars in .std % Now, use CHAIN2 lid2/neo and teo for nerho, terho % 'ece' = electron cyclotron emission % 'sxr' = soft x-ray emission % 'sxR' = soft x-ray emission with varargout{1} option (requires varargin{5}!) % 'halpha'= Dalpha signal % % Special case compatible with old gdat.m allows (JET related): % gdat(51994,'ppf','efit/xip',...) % omitting the 'JET' input as assumes JET if 3rd argument is a string % % examples: % aa=gdat(51994,'ppf/efit/xip',1,'JET'); % aa=gdat(55379,'jpf/di/c1f-chan8/131?type=lpf+diag=kc1f',1); % KC1F % aa=gdat(53290,'jpf/di/c3-cats<c:001?type=lpf+diag=cats1',0,'JET') % long magnetic 8s by 001, 002, 003, 004 blocks % aa=gdat(53290,'jpf/di/c3-cats<h:301?type=lpf+diag=cats1',0,'JET'); % fast data coils 301, 302, etc % sxr10=gdat(53290,'jpf/db/j3-sxr<v10/1',1,'JET'); % % INPUT: % shot: shot number % data_type: type of the required data. % % Allows extension for uid and seq number a la RDA: ?uid=jetthg+seq=110 % examples: % % data_type='Ip?uid=jetthg+seq=110' % data_type='ppf','efit/xip?uid=jetthg+seq=110' % % for EFIT traces, allows keyword extension '_m' to get data from ppf/efitm instead of ppf/efit % examples: % % data_type='Ip_m?uid=jetthg+seq=110' % data_type='ppf','efitm/xip?uid=jetthg+seq=110' % % Meaning of varargin depends on data_type: % % data_type=sxr or ece: % varargin{1}: [i1 i2] : if not empty, assumes need many chords from i1 to i2 % varargin{2}: channel status: 1=unread yet, 0=read % (for traces with many channel, enables to load additional channels, % like SXR, ECE, etc.) % varargin{3}: zmag for varargout{1} computation % varargin{4}, {5} (used by AUG, time interval and skip interval to reduce time points) % varargin{6}: sxr camera extra option: use to choose R projection: % empty or '0'=default (intersection with zmag with 5 degrees) % '1'=fixed to B. Alper R array (fixed in time) % % OUTPUT: % trace.data: data structure % trace.t: time of reference % trace.x: space of reference % error: error in loading signal (0=> OK, 1=> error) % % Additional Output arguments depending on data_type % % data_type=sxR: % varargout{1}: intersection of the view lines with magnetic axis % uses 5 degrees and zmag (varargin{3} if given) as default % if varargin{4}==1, uses B. Alper fixed R array for V camera % % functions needed: jetreaddata or mdsplus routines % % Example: % [zmag,error]=loadJETdata(shot,'zmag'); % [n2,error]=loadJETdata(shot,'jpf/da/c1-g102'); % [halpha,error]=loadJETdata(shot,'jpf/dd/s3-ad35'); % varargout={cell(1,1)}; error=1; % To allow multiple ways of writing a specific keyword, use data_type_eff within this routine data_type_eff=data_type; if size(data_type_eff,1)==1 i=findstr('/',data_type_eff); if length(i)>1 % assumes given a la 'ppf/efit/xip' data_type_eff=[{data_type_eff(1:i(1)-1)} ; {data_type_eff(i(1)+1:end)}]; end end i_efitm=0; i_ext=length(data_type_eff)+1; name_ext=''; if size(data_type_eff,1)==1 i=findstr('?',data_type_eff); if ~isempty(i) i_ext=i; name_ext=data_type_eff(i_ext:end); end data_type_eff_noext=data_type_eff(1:i_ext-1); i=findstr('_m',data_type_eff_noext); if ~isempty(i) i_efitm=1; data_type_eff_noext=data_type_eff(1:i-1); end if ~isempty(strmatch(data_type_eff_noext,[{'ip'} {'i_p'} {'xip'}],'exact')) data_type_eff_noext='Ip'; end if ~isempty(strmatch(data_type_eff_noext,[{'Te'} {'t_e'} {'TE'} {'T_e'}],'exact')) data_type_eff_noext='te'; end if ~isempty(strmatch(data_type_eff_noext,[{'Ne'} {'n_e'} {'NE'} {'N_e'}],'exact')) data_type_eff_noext='ne'; end if ~isempty(strmatch(data_type_eff_noext,[{'Terho'}],'exact')) data_type_eff_noext='terho'; end if ~isempty(strmatch(data_type_eff_noext,[{'SXR'}],'exact')) data_type_eff_noext='sxr'; end if ~isempty(strmatch(data_type_eff_noext,[{'ECE'}],'exact')) data_type_eff_noext='ece'; end if ~isempty(strmatch(data_type_eff_noext,[{'VOL'} {'volume'}],'exact')) data_type_eff_noext='vol'; end if ~isempty(strmatch(data_type_eff_noext,[{'q_95'} {'Q95'}],'exact')) data_type_eff_noext='q95'; end if ~isempty(strmatch(data_type_eff_noext,[{'elongation'} {'elon'}],'exact')) data_type_eff_noext='kappa'; end if ~isempty(strmatch(data_type_eff_noext,[{'triangularity'} {'triang'}],'exact')) data_type_eff_noext='delta'; end if ~isempty(strmatch(data_type_eff_noext,[{'deltaup'} {'deltau'} {'triangtop'} {'triangu'} {'triangup'}],'exact')) data_type_eff_noext='deltatop'; end if ~isempty(strmatch(data_type_eff_noext,[{'deltalow'} {'deltal'} {'triangbot'} {'triangl'} {'trianglow'}],'exact')) data_type_eff_noext='deltabot'; end if ~isempty(strmatch(data_type_eff_noext,[{'halpha'} {'Halpha'}],'exact')) data_type_eff_noext='halpha'; end if ~isempty(strmatch(data_type_eff_noext,[{'n1'} {'N1'}],'exact')) data_type_eff_noext='n1'; end if ~isempty(strmatch(data_type_eff_noext,[{'n2'} {'N2'}],'exact')) data_type_eff_noext='n2'; end if ~isempty(strmatch(data_type_eff_noext,[{'Rmag'}],'exact')) data_type_eff_noext='rmag'; end if ~isempty(strmatch(data_type_eff_noext,[{'Zmag'}],'exact')) data_type_eff_noext='zmag'; end if ~isempty(strmatch(data_type_eff_noext,[{'Rcont'}],'exact')) data_type_eff_noext='rcont'; end if ~isempty(strmatch(data_type_eff_noext,[{'Zcont'}],'exact')) data_type_eff_noext='zcont'; end else i_ext=length(data_type_eff{2})+1; name_ext=''; i=findstr('?',data_type_eff{2}); if ~isempty(i) i_ext=i; name_ext=data_type_eff{2}(i_ext:end); end data_type_eff_noext=data_type_eff{2}(1:i_ext-1); end % all keywords and corresponding case to run below JETkeywrdall=[{'Ip'} {'zmag'} {'rmag'} {'rcont'} {'zcont'} {'vol'} {'qrho'} {'q95'} {'kappa'} ... {'delta'} {'deltatop'} {'deltabot'} {'halpha'} {'n1'} {'n2'} {'neint'} ... {'ne'} {'te'} {'nerho'} {'terho'} ... {'sxr'} {'sxR'} {'ece'}]; JETsig.iip=strmatch('Ip',JETkeywrdall,'exact'); JETsig.izmag=strmatch('zmag',JETkeywrdall,'exact'); JETsig.irmag=strmatch('rmag',JETkeywrdall,'exact'); JETsig.ircont=strmatch('rcont',JETkeywrdall,'exact'); JETsig.izcont=strmatch('zcont',JETkeywrdall,'exact'); JETsig.ivol=strmatch('vol',JETkeywrdall,'exact'); JETsig.iqrho=strmatch('qrho',JETkeywrdall,'exact'); JETsig.iq95=strmatch('q95',JETkeywrdall,'exact'); JETsig.ikappa=strmatch('kappa',JETkeywrdall,'exact'); JETsig.idelta=strmatch('delta',JETkeywrdall,'exact'); JETsig.ideltatop=strmatch('deltatop',JETkeywrdall,'exact'); JETsig.ideltabot=strmatch('deltabot',JETkeywrdall,'exact'); JETsig.ihalpha=strmatch('halpha',JETkeywrdall,'exact'); JETsig.in1=strmatch('n1',JETkeywrdall,'exact'); JETsig.in2=strmatch('n2',JETkeywrdall,'exact'); JETsig.ineint=strmatch('neint',JETkeywrdall,'exact'); JETsig.ine=strmatch('ne',JETkeywrdall,'exact'); JETsig.ite=strmatch('te',JETkeywrdall,'exact'); JETsig.inerho=strmatch('nerho',JETkeywrdall,'exact'); JETsig.iterho=strmatch('terho',JETkeywrdall,'exact'); JETsig.isxr=strmatch('sxr',JETkeywrdall,'exact'); JETsig.isxR=strmatch('sxR',JETkeywrdall,'exact'); JETsig.iece=strmatch('ece',JETkeywrdall,'exact'); % For each keyword, specify which case to use. As most common is 'simplereaddata', fill in with this and change % only indices needed. Usually use name of case same as keyword name JETkeywrdcase=cell(size(JETkeywrdall)); JETkeywrdcase(:)={'simplereaddata'}; JETkeywrdcase(JETsig.iqrho)=JETkeywrdall(JETsig.iqrho); % special as efit q on psi JETkeywrdcase(JETsig.idelta)=JETkeywrdall(JETsig.idelta); % special as average of triu and tril JETkeywrdcase(JETsig.ine)=JETkeywrdall(JETsig.ine); % special as adds error bars JETkeywrdcase(JETsig.ite)=JETkeywrdall(JETsig.ite); % idem JETkeywrdcase(JETsig.inerho)=JETkeywrdall(JETsig.inerho); % idem JETkeywrdcase(JETsig.iterho)=JETkeywrdall(JETsig.iterho); % idem JETkeywrdcase(JETsig.isxr)=JETkeywrdall(JETsig.isxr); JETkeywrdcase(JETsig.isxR)=JETkeywrdall(JETsig.isxR); JETkeywrdcase(JETsig.iece)=JETkeywrdall(JETsig.iece); % Information about which dimension has time, always return 2D data as (x,t) array % as most are 1D arrays with time as first index, fill in with ones and change only those needed JETsigtimeindx=ones(size(JETkeywrdall)); % For the 'simplereaddata' cases, we need the full node in case gdat was called with full location directly % for the other cases, leave this location empty JETsiglocation=cell(2,size(JETkeywrdall,2)); JETsiglocation(:)={''}; JETsiglocation(:,JETsig.iip)={'ppf'; 'efit/xip'}; JETsiglocation(:,JETsig.izmag)={'ppf'; 'efit/zmag'}; JETsiglocation(:,JETsig.irmag)={'ppf'; 'efit/rmag'}; JETsiglocation(:,JETsig.ircont)={'ppf' ; 'efit/rbnd'}; JETsigtimeindx(JETsig.ircont)=2; JETsiglocation(:,JETsig.izcont)={'ppf' ; 'efit/zbnd'}; JETsigtimeindx(JETsig.izcont)=2; JETsiglocation(:,JETsig.ivol)={'ppf'; 'equi/vol'}; JETsiglocation(:,JETsig.iq95)={'ppf'; 'efit/q95'}; JETsiglocation(:,JETsig.ikappa)={'ppf'; 'efit/elon'}; JETsiglocation(:,JETsig.ideltatop)={'ppf'; 'efit/triu'}; JETsiglocation(:,JETsig.ideltabot)={'ppf'; 'efit/tril'}; JETsiglocation(:,JETsig.ihalpha)={'jpf'; 'dd/s3-ad35'}; JETsiglocation(:,JETsig.in1)={'jpf'; 'da/c1-g101'}; JETsiglocation(:,JETsig.in2)={'jpf'; 'da/c1-g102'}; JETsiglocation(:,JETsig.ineint)={'ppf'; 'kg1v/lid3'}; % initialize order of substructures and allows just a "return" if data empty trace.data=[]; trace.x=[]; trace.t=[]; trace.dim=[]; trace.dimunits=[]; trace.name=[]; % find index of signal called upon if size(data_type_eff,1)==2 % in case node name was given in 2 parts directly (as old way) ii1=strmatch(data_type_eff(1),JETsiglocation(1,:),'exact'); iiindex=strmatch(data_type_eff_noext,JETsiglocation(2,ii1),'exact'); if ~isempty(iiindex) index=ii1(iiindex); else index=[]; end if isempty(index) % $$$ disp('********************') % $$$ disp('trace not yet registered.') % $$$ disp('If standard data, ask andrea.scarabosio@epfl.ch or olivier.sauter@epfl.ch to create a keyqord entry for this data') % eval(['!mail -s ''' data_type_eff{1} ' ' data_type_eff{2} ' ' num2str(shot) ' ' ... % getenv('USER') ' JET'' olivier.sauter@epfl.ch < /dev/null']) disp('********************') % temporarily add entry in arrays, so can work below index=length(JETkeywrdall)+1; JETkeywrdall(end+1)={'new'}; JETkeywrdcase(end+1)={'simplereaddata'}; JETsiglocation(1:2,end+1)=[data_type_eff(1) ; {data_type_eff_noext}]; JETsigtimeindx(end+1)=0; elseif ~strcmp(JETkeywrdcase{index},'simplereaddata') msgbox(['Problem in loadJETdata with data_type_eff = ' char(data_type_eff(end)) ... '. Full paths of nodes should only be for case simplereaddata'],'in loadJETdata','error') error('in loadJETdata') end else index=strmatch(data_type_eff_noext,JETkeywrdall,'exact'); if isempty(index) disp(' ') disp('********************') if iscell(data_type_eff) disp(['no such keyword: ' data_type_eff{1} '/' data_type_eff{2}]) else disp(['no such keyword: ' data_type_eff]) end disp(' ') disp('Available keywords:') JETkeywrdall(:) disp('********************') return end end disp(' ') if iscell(data_type_eff) disp(['loading' ' ' data_type_eff{1} '/' data_type_eff{2} ' from JET shot #' num2str(shot)]); else disp(['loading' ' ' data_type_eff ' from JET shot #' num2str(shot)]); end disp(['case ' JETkeywrdcase{index}]) disp(' ') switch JETkeywrdcase{index} %&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&& case 'simplereaddata' ppftype=JETsiglocation{1,index}; if i_efitm; tracename=['eftm' JETsiglocation{2,index}(5:end) name_ext]; else tracename=[JETsiglocation{2,index} name_ext]; end ij=find(tracename~=''''); tracename=tracename(ij); keyboard [a,x,t,d,e]=rda_eff(shot,ppftype,tracename); switch tracename case {'efit/btpd','efit/btpd?uid=jetppf+seq=0'} if isempty(a) | isempty(t); disp('data or t empty, assumes means btpd not defined'); [xip,x,t,d,e]=rda_eff(shot,'ppf','efit/xip'); shot_mg3_list=[47274 47275 47276 47280 47281 47282 47283 47284 47285 47286 47287 47290 47295 47296 47301]; if isempty(find(shot_mg3_list==shot)) [wdia,x1,t1,d,e]=rda_eff(shot,'ppf','efit/wdia'); else [wdia,x1,t1,d,e]=rda_eff(shot,'ppf','mg3/wpd'); wdia=interp1(t1,wdia,t); end [rgeo,x3,t3,d,e]=rda_eff(shot,'ppf','efit/rgeo'); a=2.122e6 .* wdia ./xip.^2 ./ rgeo; end case {'efit/btnd','efit/btnd?uid=jetppf+seq=0'} if isempty(a) | isempty(t); disp('data or t empty, assumes means btnd not defined'); [xip,x,t,d,e]=rda_eff(shot,'ppf','efit/xip'); shot_mg3_list=[47274 47275 47276 47280 47281 47282 47283 47284 47285 47286 47287 47290 47295 47296 47301]; if isempty(find(shot_mg3_list==shot)) [wdia,x1,t1,d,e]=rda_eff(shot,'ppf','efit/wdia'); else [wdia,x1,t1,d,e]=rda_eff(shot,'ppf','mg3/wpd'); wdia=interp1(t1,wdia,t); end [rgeo,x3,t3,d,e]=rda_eff(shot,'ppf','efit/rgeo'); [cr0,x3,t3,d,e]=rda_eff(shot,'ppf','efit/cr0'); [bvac,x3,t3,d,e]=rda_eff(shot,'ppf','efit/bvac'); [volm,x3,t3,d,e]=rda_eff(shot,'ppf','efit/volm'); a=56.605.*wdia.*cr0.*rgeo./xip./volm./bvac; end case {'LIDR/Z','lidr/z','Lidr/Z','LIDR/Z?uid=jetppf+seq=0'} % 1D but vs R instead of t x=t; t=[]; end trace.data=a; trace.x=x; trace.t=t; clear error error=e; if length(size(trace.data))==1 | (length(size(trace.data))==2 & size(trace.data,2)==1) trace.dim=[{trace.t}]; trace.dimunits={'time [s]'}; elseif length(size(trace.data))==2 trace.dim=[{trace.x} ; {trace.t}]; trace.dimunits=[{'R [m] or rho=sqrt(psi_norm)'} ; {'time [s]'}]; else disp('how to deal with 3D arrays?') trace.dim=[{trace.x} ; {trace.t} ; {d}]; trace.dimunits=[{'R [m] or rho=sqrt(psi_norm)'} ; {'time [s]'} ; {'d'}]; trace.d=d; end trace.name=[ppftype '/' num2str(shot) '/' tracename]; %&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&& case JETkeywrdall{JETsig.iqrho} % q profile on sqrt(psi_norm) ppftype='ppf'; if i_efitm tracename=['eftm/q' name_ext]; else tracename=['efit/q' name_ext]; end [a,x,t,d,e]=rda_eff(shot,ppftype,tracename); trace.data=a; trace.x=sqrt(x); % x is psi (? to test) trace.t=t; trace.dim=[{trace.x} ; {trace.t}]; trace.dimunits=[{'sqrt(\psi)'} ; {'time [s]'}]; trace.name=[ppftype '/' num2str(shot) '/' tracename]; error=e; %&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&& case JETkeywrdall{JETsig.idelta} % average of delatop and deltabot ppftype='ppf'; tracename1=['efit/triu' name_ext]; tracename2=['efit/tril' name_ext]; [a1,x,t,d,e]=rda_eff(shot,ppftype,tracename1); [a2,x,t,d,e]=rda_eff(shot,ppftype,tracename2); trace.data=0.5.*(a1+a2); trace.x=x; trace.t=t; trace.dim=[{trace.t}]; trace.dimunits=[{'time [s]'}]; trace.name=[ppftype '/' num2str(shot) '/efit/0.5(triu+tril)']; error=e; %&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&& case {JETkeywrdall{JETsig.ine} , JETkeywrdall{JETsig.ite}} % ne, te raw data from LIDR vs R,t. Add error bars ppftype='ppf'; if strcmp(JETkeywrdcase{index},JETkeywrdall{JETsig.ine}) tracename=['LIDR/NE' name_ext]; else tracename=['LIDR/TE' name_ext]; end [a,x,t,d,e]=rda_eff(shot,ppftype,tracename); trace.data=a; trace.x=x; trace.t=t; trace.dim=[{trace.x} ; {trace.t}]; trace.dimunits=[{'R [m]'} ; {'time [s]'}]; trace.std=[]; trace.name=[ppftype '/' num2str(shot) '/' tracename]; error=e; %&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&& case {JETkeywrdall{JETsig.inerho} , JETkeywrdall{JETsig.iterho}} % ne, te on rho mesh. use lid2, thus need chain2 to have been run. Add error bars ppftype='ppf'; if strcmp(JETkeywrdcase{index},JETkeywrdall{JETsig.inerho}) tracename=['LID2/NEO' name_ext]; else tracename=['LID2/TEO' name_ext]; end [a,x,t,d,e]=rda_eff(shot,ppftype,tracename); trace.data=a; trace.x=x; trace.t=t; trace.dim=[{trace.x} ; {trace.t}]; trace.dimunits=[{'rho=sqrt(psi)'} ; {'time [s]'}]; trace.std=[]; trace.name=[ppftype '/' num2str(shot) '/' tracename]; error=e; %&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&& case {'sxr','sxR'} % LOAD MULTI CHANNEL DATA % load JET soft x-ray data % parameters needed for correct convertion of JET Sxr data vconvert= [1.379 1.311 1.249 1.191 1.139 1.093 1.049 ... 1.011 0.975 0.945 0.917 0.893 0.873 0.856 ... 0.842 0.829 0.821 0.815 0.821 0.829 0.842 ... 0.856 0.873 0.894 0.918 0.946 0.976 1.012 ... 1.050 1.094 1.141 1.193 1.251 1.313 1.382]; rconvert= [3.45 3.43 3.41 3.37 3.33 3.28 3.23 3.18 3.14 ... 3.09 3.05 3.00 2.94 2.89 2.83 2.77 2.72 2.68 2.63 ... 2.59 2.55 2.49 2.44 2.40 2.37 2.33 2.29 2.26 2.23 ... 2.19 2.14 2.12 2.10 2.08 2.06]; if nargin>=3 & ~isempty(varargin{1}) starti=varargin{1}(1); endi=varargin{1}(2); else starti=1; endi=24; end if nargin>=4 & ~isempty(varargin{2}) status=varargin{2}; else status=ones(endi-starti+1,1); end trace.t=[]; trace.x=[]; iloaded_data=0; for i=starti:endi % Read channels from lowchannel to upchannel if necessary if status(i)==1 iloaded_data=iloaded_data+1; % Status=1 => Not Read Yet % vertical SXR chords ppftype='jpf'; tracename=['db/j3-sxr<v' num2str(i) '/1' name_ext]; if shot<48000 tracename=['db/j3-sxr<t' num2str(i) '/1' name_ext]; disp('Using T camera: Radius data probably wrong') end a=which('jpfdat'); if isempty(a) [a,x,t,d,e]=rda_eff(shot,ppftype,tracename); % Convert from raw sxr data to W/m^2 % EDIT: on JET a is not empty on error if ~isempty(a) if ~isempty(t) trace.data(i,:) = a * vconvert(i); trace.t=t; if ~isempty(x) trace.x(i,:)=x; end end error=e; end else disp(['Reading channel ' tracename]); [a, t, nwds, title, unit, ier] = jpfdat(tracename, shot); if ~ier % Convert from raw sxr data to W/m^2 trace.data(i,:) = a' * vconvert(i); trace.t=t'; end end trace.dim=[{[starti:endi]'} ; {trace.t}]; trace.dimunits=[{'channels'} ; {'time [s]'}]; trace.name=[ 'jpf/' num2str(shot) '/' 'db/j3-sxr<vXX' '/1' name_ext]; end end if isempty(trace.t) disp(['no data in ' trace.name]) return end % calculating intersection of the view lines with magnetics axis if strcmp(data_type_eff_noext,'sxR') if iloaded_data>0 if nargin>=5 & ~isempty(varargin{3}) zmag=varargin{3}; else zmag=loadJETdata(shot,'zmag'); end zmageff=interp1(zmag.t,zmag.data,trace.t); for i=starti:endi radius.data(i,:)=2.848 + (2.172-zmageff') .* tan(-5.0/180.*3.14159 - atan2(0.99.*(i-18),35.31)); end iii=0; if nargin>=8 & ~isempty(varargin{6}) iii=str2num(varargin{6}); end if iii==1 for i=starti:endi radius.data(i,:)=rconvert(i); end disp('uses B. Alper fixed R array for SXR R intersection') end radius.t=t; varargout{1}={radius}; trace.R=radius.data; else varargout{1}={struct([])}; trace.R=[]; end end %&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&& case 'ece' if nargin>=3 & ~isempty(varargin{1}) starti=varargin{1}(1); endi=varargin{1}(2); else starti=1; endi=24; end if nargin>=4 & ~isempty(varargin{2}) status=varargin{2}; else status=ones(endi,1); end % Read channels from lowchannel to upchannel if necessary for i=starti:endi if status(i)==1 % ECE, te0 % Status=1 => Not Read Yet ppftype='ppf'; tracename=['kk3/te' num2str(i,'%2.2d') name_ext]; disp(tracename) a=which('ppfread'); if isempty(a) | ~isempty(name_ext) [a,x,t,d,e]=rda_eff(shot,ppftype,tracename); if isempty(a) & size(trace.data,2)>1 trace.data(i,:)=NaN; else trace.data(i,:)=a; trace.t=t; end if ~isempty(x); trace.x(i,:)=x'; end; error=e; else [a,x,t,unitd,unitx,unitt,comment,sequence,e]= ... ppfread(shot,'KK3',['TE' num2str(i,'%2.2d')]); if isempty(a) & size(trace.data,2)>1 trace.data(i,:)=NaN; else trace.data(i,:)=a'; trace.t=t'; end if ~isempty(x); trace.x(i,:)=x; end; error=e; end ppftypeR='ppf'; tracenameR=['kk3/rc' num2str(i,'%2.2d') name_ext]; a=which('ppfdat'); if isempty(a) | ~isempty(name_ext) [a,x,t,d,e]=rda_eff(shot,ppftypeR,tracenameR); if isempty(a) & size(trace.data,2)>1 radius.data(i,:)=NaN; else radius.data(i,:)=a; radius.t=t; radius.x=x; end else [a,x,t,unitd,unitx,unitt,comment,sequence,e]= ... ppfread(shot,'KK3',['RC' num2str(i,'%2.2d')]); if isempty(a) & size(trace.data,2)>1 radius.data(i,:)=NaN; else radius.data(i,:)=a'; radius.t=t'; radius.x=x'; end end end end trace.dim=[{[starti:endi]'} ; {trace.t}]; trace.dimunits=[{'channels'} ; {'time [s]'}]; trace.name=[ 'ppf/' num2str(shot) '/' 'kk3/teXX' name_ext ]; if exist('radius') varargout={{radius}}; for i=starti:endi trace.R(i,:)=interp1(radius.t,radius.data(i,:),trace.t); end else varargout={{struct([])}}; trace.R=[]; end otherwise disp('case not yet defined') end