interpleg_yp

PURPOSE ^

SYNOPSIS ^

function [XXf0,x_0] = interpleg_yp(XXf0_leg,mhu,momentumDKE)

DESCRIPTION ^

 Pitch-angle interpolation on a new mhu grid

 INPUT:

   - XXf0_leg: Legendre coefficients of the bounce-averaged distribution function (npn,n_pol_legendre + 1,nr)
   - mhu : new pitch-angle grid on which the distribution function will be interpoled (nmhu,1)
   - momentumDKE structure for calculating the interpolated normalization
     (optional)

 OUTPUT:

   - XXf0: interpolation distribution function (npn,nmhu,nr)
   - x_0: normalization at all radii (nr)

 By Yves Peysson (CEA-IRFM, yves.peysson@cea.fr) and Joan Decker (CEA-IRFM, joan.decker@cea.fr)

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 function [XXf0,x_0] = interpleg_yp(XXf0_leg,mhu,momentumDKE)
0002 %
0003 % Pitch-angle interpolation on a new mhu grid
0004 %
0005 % INPUT:
0006 %
0007 %   - XXf0_leg: Legendre coefficients of the bounce-averaged distribution function (npn,n_pol_legendre + 1,nr)
0008 %   - mhu : new pitch-angle grid on which the distribution function will be interpoled (nmhu,1)
0009 %   - momentumDKE structure for calculating the interpolated normalization
0010 %     (optional)
0011 %
0012 % OUTPUT:
0013 %
0014 %   - XXf0: interpolation distribution function (npn,nmhu,nr)
0015 %   - x_0: normalization at all radii (nr)
0016 %
0017 % By Yves Peysson (CEA-IRFM, yves.peysson@cea.fr) and Joan Decker (CEA-IRFM, joan.decker@cea.fr)
0018 %
0019 if nargin < 2,error('Wrong number of input arguments in interpleg_yp.m.');end
0020 %
0021 if abs(mhu) > 1,error('Inconsistent input in interpleg_yp.m: pitch-angle values must lie between -1 and 1.');end
0022 %
0023 if ~exist('momentumDKE'),x_0 = [];end
0024 %
0025 XXf0 = zeros(size(XXf0_leg,1),length(mhu),size(XXf0_leg,3));
0026 %
0027 [Pm,cL] = leg_dke_yp(mhu(:)',size(XXf0_leg,2) - 1);%Legendre polynomial calculation
0028 %
0029 %
0030 for ir = 1:size(XXf0_leg,3),%radial position
0031     XXf0(:,:,ir) = (Pm'*(XXf0_leg(:,:,ir).*(ones(size(XXf0_leg,1),1)*cL'))')';%Distribution function on the Legendre polynomials pitch-angle grid back calculated
0032     %
0033     if exist('momentumDKE'),
0034         x_0(ir) = 2*pi*integral_dke_jd(momentumDKE.dpn',momentumDKE.pn'.*momentumDKE.pn'.*trapz_dke_yp([mhu(:),XXf0(:,:,ir)'])',1);
0035     end
0036 end
0037 %
0038

Community support and wiki are available on Redmine. Last update: 18-Apr-2019.