hotdisp_dke_jd

PURPOSE ^

This function solves the dispersion relation and the wave equation

SYNOPSIS ^

function [Nperp,e_xyz,e_pmz,e_pyk,D,X,X_s,X_sn,alphaphi_lin,phiT_xyz,phiP_xyz,phiT_pmz,phiP_pmz,phiT_pyk,phiP_pyk]= hotdisp_dke_jd(Npar,omega_rf,Te,ne,zTi,zni,zZi,zmi,Bfield,Nperp_in,herm_mode,display_mode,nmax,NZ,tol,method_type)

DESCRIPTION ^

 This function solves the dispersion relation and the wave equation
 in the linear kinetic plasma model. It calculates the perpendicular index of
 refraction and the polarization as a function of the wave
 frequency and the parallel index of refraction. 

 - (x,y,z): We assume that the magnetic field is in the z direction and the 
           perpendicular index of refraction in the x direction. 
 - (p,m,z): p and m refer to the left hand and right hand rotating fields 
           respectively.
 - (p,y,k): p and y are the transverse components of the polarization, with 
           respect to the wave vector. k is the longitudinal component.

   INPUTS: 

       - Npar: parallel index of refraction [1,1]       
       - omega_rf: wave frequency (rad/s) [1,1] 
       - Te: electron temperature (keV) [1,1]
       - ne: electron density (m-3) [1,1]
       - zTi: ion temperature (keV) [1,i]
       - zni: ion density (m-3) [1,i]
       - zZi: ion charges [1,i] 
       - zmi: ion masses (mp) [1,i]  
       - Bfield: magnetic field amplitude [1,1]
       - Nperp_in: perpendicular index of refraction (initial guess) [1,1] (default:1)    
       - herm_mode: mode for the calculation (0) full dielectric tensor, (1) hermitian part (default:0) 
           >>> In the case herm_mode = 1, the power flow and linear absorption
               are also calculated
       - display_mode: (0,1) none (2) convergence figures (default:0)  
       - nmax: maximum harmonic number (default:10) 
       - NZ: number of grid points in Z function integration calculation (default:100) 
       - tol: tolerance on the value of the dispersion relation (default:1e-6) 
       - method_type: 'newton', 'segments' (default:'newton') 

   OUTPUTS:

       - Nperp: perpendicular index of refraction [1,1]
       - e_xyz: polarization vector in the (x,y,z) coordinates [3,1]
       - e_pmz: polarization vector in the (p,m,z) coordinates [3,1] 
       - e_pyk: polarization vector in the (p,y,k) coordinates [3,1]  
       - X: kinetic plasma susceptibility tensor elements [3,3]
       - X_s: kinetic plasma susceptibility tensor elements for species s [3,3,1+i]
       - X_sn: kinetic plasma susceptibility tensor elements for species s and harmonic n [3,3,1+i,N]

 created 05/02/2003
 by Joan Decker <jodecker@alum.mit.edu> (MIT/RLE) and  Yves Peysson <yves.peysson@cea.fr> (CEA/DRFC)

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 function [Nperp,e_xyz,e_pmz,e_pyk,D,X,X_s,X_sn,alphaphi_lin,phiT_xyz,phiP_xyz,phiT_pmz,phiP_pmz,phiT_pyk,phiP_pyk]...
0002     = hotdisp_dke_jd(Npar,omega_rf,Te,ne,zTi,zni,zZi,zmi,Bfield,Nperp_in,herm_mode,display_mode,nmax,NZ,tol,method_type)
0003 % This function solves the dispersion relation and the wave equation
0004 % in the linear kinetic plasma model. It calculates the perpendicular index of
0005 % refraction and the polarization as a function of the wave
0006 % frequency and the parallel index of refraction.
0007 %
0008 % - (x,y,z): We assume that the magnetic field is in the z direction and the
0009 %           perpendicular index of refraction in the x direction.
0010 % - (p,m,z): p and m refer to the left hand and right hand rotating fields
0011 %           respectively.
0012 % - (p,y,k): p and y are the transverse components of the polarization, with
0013 %           respect to the wave vector. k is the longitudinal component.
0014 %
0015 %   INPUTS:
0016 %
0017 %       - Npar: parallel index of refraction [1,1]
0018 %       - omega_rf: wave frequency (rad/s) [1,1]
0019 %       - Te: electron temperature (keV) [1,1]
0020 %       - ne: electron density (m-3) [1,1]
0021 %       - zTi: ion temperature (keV) [1,i]
0022 %       - zni: ion density (m-3) [1,i]
0023 %       - zZi: ion charges [1,i]
0024 %       - zmi: ion masses (mp) [1,i]
0025 %       - Bfield: magnetic field amplitude [1,1]
0026 %       - Nperp_in: perpendicular index of refraction (initial guess) [1,1] (default:1)
0027 %       - herm_mode: mode for the calculation (0) full dielectric tensor, (1) hermitian part (default:0)
0028 %           >>> In the case herm_mode = 1, the power flow and linear absorption
0029 %               are also calculated
0030 %       - display_mode: (0,1) none (2) convergence figures (default:0)
0031 %       - nmax: maximum harmonic number (default:10)
0032 %       - NZ: number of grid points in Z function integration calculation (default:100)
0033 %       - tol: tolerance on the value of the dispersion relation (default:1e-6)
0034 %       - method_type: 'newton', 'segments' (default:'newton')
0035 %
0036 %   OUTPUTS:
0037 %
0038 %       - Nperp: perpendicular index of refraction [1,1]
0039 %       - e_xyz: polarization vector in the (x,y,z) coordinates [3,1]
0040 %       - e_pmz: polarization vector in the (p,m,z) coordinates [3,1]
0041 %       - e_pyk: polarization vector in the (p,y,k) coordinates [3,1]
0042 %       - X: kinetic plasma susceptibility tensor elements [3,3]
0043 %       - X_s: kinetic plasma susceptibility tensor elements for species s [3,3,1+i]
0044 %       - X_sn: kinetic plasma susceptibility tensor elements for species s and harmonic n [3,3,1+i,N]
0045 %
0046 % created 05/02/2003
0047 % by Joan Decker <jodecker@alum.mit.edu> (MIT/RLE) and  Yves Peysson <yves.peysson@cea.fr> (CEA/DRFC)
0048 %
0049 if nargin < 16,
0050    method_type = 'newton';
0051 end
0052 if nargin < 15,
0053    tol = 1e-5;
0054 end
0055 if nargin < 14,
0056    NZ = 1000;
0057 end
0058 if nargin < 13,
0059    nmax = 10;
0060 end
0061 if nargin < 12,
0062    display_mode = 0;
0063 end
0064 if nargin < 11,
0065    herm_mode = 0;
0066 end
0067 if nargin < 10,
0068    Nperp_in = 1.0;
0069 end
0070 %
0071 [qe,me,mp,mn,e0,mu0,re,mc2,clum] = pc_dke_yp;
0072 %
0073 ns = [ne;zni(:)];
0074 qs = qe*[-1;zZi(:)];
0075 ms = [me;mp*zmi(:)];
0076 %
0077 wps = sqrt(ns.*qs.^2/e0./ms);
0078 wcs = qs*Bfield./ms;
0079 %
0080 ps = wps/omega_rf;
0081 ys = wcs/omega_rf;
0082 %
0083 bTs = sqrt([Te*1e3*qe/me/clum^2;zTi(:)*1e3*qe./(mp*zmi(:))/clum^2]);
0084 %
0085 [Nperp,e_xyz,e_pmz,e_pyk,D,X,X_s,X_sn,...
0086           alphaphi_lin,phiT_xyz,phiP_xyz,phiT_pmz,phiP_pmz,phiT_pyk,phiP_pyk]...
0087     = kineticdisp_jd(ps.^2,ys,bTs,Npar,Nperp_in,herm_mode,display_mode,nmax,NZ,tol,method_type);
0088 %

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