R5-X2 - Calculation of the plasma bremsstrahlung emission Calculates the plasma bremsstrahlung emission. Call first the function bremchord_dke_yp.m INPUTS: - Zbremchord: HXR diagnostic chord arrangement data structure - equilHXR: HXR magnetic equilibrium structure - hxr: Bremsstrahlung structure - hxrparam : Bremsstrahlung calculation and display parameter structure - dkeparam: DKE data structure - radialDKE: DKE radial grid structure - mksa: MKSA structure - momentumDKE: DKE momentum structure - dke_out: DKE data output structure (distributions) - niewmode: (NaN) exact topology (0) perp (1) forward (-1) backward OUTPUTS: - Zbremplasma: Plasma bremsstrahlung data structure * Zbremplasma.kphot : photon energies (keV) [1,nk] * Zbremplasma.brem_bd_out : photon number per unit energy (keV) and time (s) for each chord [nchord,nk] * Zbremplasma.lchord : chord length (m) [1,nchord] by Joan Decker (CEA/IRFM,joan.decker@cea.fr) and Yves Peysson(CEA/IRFM,yves.peysson@cea.fr)
0001 function [Zbremplasma] = bremchordemission_dke_jd(Zbremplasma,Zbremchord,hxr,viewmode) 0002 %R5-X2 - Calculation of the plasma bremsstrahlung emission 0003 % 0004 % Calculates the plasma bremsstrahlung emission. Call first the function bremchord_dke_yp.m 0005 % 0006 % INPUTS: 0007 % 0008 % - Zbremchord: HXR diagnostic chord arrangement data structure 0009 % - equilHXR: HXR magnetic equilibrium structure 0010 % - hxr: Bremsstrahlung structure 0011 % - hxrparam : Bremsstrahlung calculation and display parameter structure 0012 % - dkeparam: DKE data structure 0013 % - radialDKE: DKE radial grid structure 0014 % - mksa: MKSA structure 0015 % - momentumDKE: DKE momentum structure 0016 % - dke_out: DKE data output structure (distributions) 0017 % - niewmode: (NaN) exact topology (0) perp (1) forward (-1) backward 0018 % 0019 % OUTPUTS: 0020 % 0021 % - Zbremplasma: Plasma bremsstrahlung data structure 0022 % * Zbremplasma.kphot : photon energies (keV) [1,nk] 0023 % * Zbremplasma.brem_bd_out : photon number per unit energy (keV) and time (s) for each chord [nchord,nk] 0024 % * Zbremplasma.lchord : chord length (m) [1,nchord] 0025 % 0026 % by Joan Decker (CEA/IRFM,joan.decker@cea.fr) and Yves Peysson(CEA/IRFM,yves.peysson@cea.fr) 0027 % 0028 if nargin < 4, 0029 viewmode = NaN; 0030 end 0031 if nargin < 3, 0032 error('Not enough input arguments in bremchordemission_dke_jd') 0033 end 0034 % 0035 vschords = Zbremchord.vschords; 0036 % 0037 nc = size(vschords,3); 0038 % 0039 zdatahxr_in = cell(1,nc); 0040 % 0041 for ic = 1:nc,%loop on each chord 0042 zdatahxr_in{ic}.vschords = vschords(:,:,ic); 0043 zdatahxr_in{ic}.EG_hxr = hxr.EG_hxr(ic); 0044 end 0045 % 0046 dkepath = load_structures_yp('dkepath','',''); 0047 % 0048 mdce_mode = 0; 0049 % 0050 clustermode.loop_brem_chord_yp.scheduler.mode = real(mdce_mode); 0051 if ~isreal(mdce_mode), 0052 clustermode.loop_brem_chord_yp.scheduler.gpu = imag(mdce_mode); 0053 else 0054 clustermode.loop_brem_chord_yp.scheduler.gpu = 0; 0055 end 0056 % 0057 dkecluster = clustermode_luke(clustermode,'loop_brem_chord_yp',dkepath); 0058 % 0059 [flag,zbrem_bd_out,zlchord] = mdce_run(@loop_brem_chord_yp,{viewmode,Zbremplasma.kphot,Zbremplasma.cL,Zbremplasma.ZZIB,''},5,zdatahxr_in,dkecluster); 0060 % 0061 for ic = 1:nc,%loop on each chord 0062 Zbremplasma.brem_bd_out(ic,:) = zbrem_bd_out{ic}; 0063 Zbremplasma.lchord(ic) = zlchord{ic}; 0064 end