rundke

PURPOSE ^

Script for running the DKE solver (can be modified by the user for specific simulations)

SYNOPSIS ^

This is a script file.

DESCRIPTION ^

Script for running the DKE solver (can be modified by the user for specific simulations)
by Y.Peysson CEA-DRFC <yves.peysson@cea.fr> and Joan Decker MIT-RLE (jodecker@mit.edu)

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 %Script for running the DKE solver (can be modified by the user for specific simulations)
0002 %by Y.Peysson CEA-DRFC <yves.peysson@cea.fr> and Joan Decker MIT-RLE (jodecker@mit.edu)
0003 %
0004 clc
0005 clear all
0006 clear mex
0007 clear functions
0008 close all
0009 warning off
0010 global nfig
0011 %
0012 p_opt = 2;
0013 %
0014 permission = test_permissions_yp;
0015 %
0016 if ~permission 
0017     disp('Please move the script to a local folder where you have write permission before to run it')
0018     return;
0019 end
0020 %
0021 % ***********************This part must be specified by the user, run make files if necessary) *****************************
0022 %
0023 id_simul = 'LH_karney_3D_num';%Simulation ID
0024 path_simul = '';%if nothing is specified, the working directory is first used and then MatLab is looking in all the path
0025 %
0026 psin_S = [];%Normalized poloidal flux grid where calculations are performed (0 < psin_S < 1) (If one value: local calculation only, not used if empty)
0027 rho_S = [0:0.05:1];%Normalized radial flux grid where calculations are performed (0 < rho_S < 1) (If one value: local calculation only, not used if empty)
0028 %
0029 id_path = '';%For all paths used by DKE solver
0030 path_path = '';%if nothing is specified, the working directory is first used and then MatLab is looking in all the path
0031 %
0032 id_equil = 'TScyl';%For plasma equilibrium
0033 path_equil = '';%if nothing is specified, the working directory is first used and then MatLab is looking in all the path
0034 %
0035 id_dkeparam = 'UNIFORM10010020';%For DKE code parameters
0036 path_dkeparam = '';%if nothing is specified, the working directory is first used and then MatLab is looking in all the path
0037 %
0038 id_display = 'NO_DISPLAY';%For output code display
0039 path_display = '';%if nothing is specified, the working directory is first used and then MatLab is looking in all the path
0040 %
0041 id_ohm = '';%For Ohmic electric contribution
0042 path_ohm = '';%if nothing is specified, the working directory is first used and then MatLab is looking in all the path
0043 %
0044 ids_wave = {''};%For RF waves contribution (put all the type of waves needed)
0045 paths_wave = {''};%if nothing is specified, the working directory is first used and then MatLab is looking in all the path
0046 %
0047 id_transpfaste = 'nomagneticturbulence';%For fast electron radial transport
0048 path_transpfaste = '';%if nothing is specified, the working directory is first used and then MatLab is looking in all the path
0049 %
0050 id_ripple = '';%For fast electron magnetic ripple losses
0051 path_ripple = '';%if nothing is specified, the working directory is first used and then MatLab is looking in all the path
0052 %
0053 %************************************************************************************************************************************
0054 %************************************************************************************************************************************
0055 %************************************************************************************************************************************
0056 %
0057 [dkepath,equil,dkeparam,dkedisplay,ohm,waves,transpfaste,ripple] = load_structures_yp('dkepath',id_path,path_path,'equil',id_equil,path_equil,'dkeparam',id_dkeparam,path_dkeparam,'dkedisplay',id_display,path_display,'ohm',id_ohm,path_ohm,'waves',ids_wave,paths_wave,'transpfaste',id_transpfaste,path_transpfaste,'ripple',id_ripple,path_ripple);
0058 %
0059 %************************************************************************************************************************************
0060 %
0061 wavestruct.omega_lh = [4]*2*pi*1e9; %(GHz -> rad/s). Wave frequency [1,1] Indicative, no effect in small FLR limit opt_lh > 0
0062 %Option parameter for cross-comparison between old LH code:
0063 %    - (1): 1/vpar dependence
0064 %    - (2): no 1/vpar dependence and old grid technique for Dlh calculations (Karney, Shoucri) (see rfdiff_dke_jd)
0065 wavestruct.opt_lh = 2; % [1,1]
0066 %
0067 % Choose (vparmin_lh,vparmax_lh) or (Nparmin_lh,Nparmax_lh) for square n// LH wave power spectrum,
0068 % or (Npar_lh,dNpar_lh) for Gaussian shape
0069 %
0070 wavestruct.norm_ref = 1;%Normalization procedure for the LH quasilinear diffusion coefficient and spectrum boundaries
0071 %
0072 wavestruct.yvparmin_lh = [4];%LH wave square N// Spectrum: Lower limit of the plateau (vth_ref or vth) [1,n_scenario_lh]
0073 wavestruct.yvparmax_lh = [7];%LH wave square N// Spectrum: Upper limit of the plateau (vth_ref or vth) [1,n_scenario_lh]
0074 %
0075 wavestruct.yNparmin_lh = [NaN];%LH wave square N// Spectrum: Lower limit [1,n_scenario_lh]
0076 wavestruct.yNparmax_lh = [NaN];%LH wave square N// Spectrum: Upper limit [1,n_scenario_lh]
0077 wavestruct.yNpar_lh = [NaN];%LH wave Gaussian N// Spectrum: peak [1,n_scenario_lh]
0078 wavestruct.ydNpar_lh = [NaN];%LH wave Gaussian N// Spectrum: width [1,n_scenario_lh]
0079 %
0080 %   Note: this diffusion coefficient is different from the general QL D0. It has a benchmarking purpose only
0081 wavestruct.yD0_in_c_lh = [1];%Central LH QL diffusion coefficient (nhuth_ref*pth_ref^2 or nhuth*pth^2) [1,n_scenario_lh]
0082 %
0083 wavestruct.yD0_in_lh_prof = [0];%Quasilinear diffusion coefficient radial profile: (0) uniform, (1) gaussian radial profile [1,n_scenario_lh]
0084 wavestruct.ypeak_lh = [NaN];%Radial peak position of the LH quasi-linear diffusion coefficient (r/a on midplane) [1,n_scenario_lh]
0085 wavestruct.ywidth_lh = [NaN];%Radial width of the LH quasi-linear diffusion coefficient (r/a on midplane) [1,n_scenario_lh]
0086 %
0087 wavestruct.ythetab_lh = [0]*pi/180;%(deg -> rad). Poloidal location of LH beam [0..2pi] [1,n_scenario_lh]
0088 %               (0) from local values Te and ne, (1) from central values Te0 and ne0
0089 %
0090 %************************************************************************************************************************************
0091 %
0092 dkeparam.boundary_mode_f = 0;%Number of points where the Maxwellian distribution is enforced from p = 0 (p=0, free conservative mode but param_inv(1) must be less than 1e-4, otherwise 1e-3 is OK most of the time. Sensitive to the number of points in p)
0093 dkeparam.norm_mode_f = 1;%Local normalization of f0 at each iteration: (0) no, the default value when the numerical conservative scheme is correct, (1) yes
0094 dkeparam.coll_mode = 0;% Relativistic Maxwellian background
0095 %
0096 dkeparam.ludroptol = 1e-6;%Drop tolerance for the built-in MatLab partial LU decomposition (recipes: droptol = max(1/dtn/nr,param_inv(1)))
0097 dkeparam.mlmaxitinv = 20;%Maximum number of iterations for iterative built-in MatLab matrix inversion method (between 1 and 20 usualy, if empty, mlmaxitinv = min(matrixdim,20). See matlab documentation.
0098 dkeparam.mlprecinv = 1e-5;%Accuracy for the MatLab build-in matrix inversion.
0099 %
0100 dkeparam.psin_S = psin_S;
0101 dkeparam.rho_S = rho_S;
0102 %
0103 waves{1} = make_idealLHwave_jd(equil,wavestruct);
0104 %
0105 np_list = 25:25:250;
0106 nmhu_list = 25:25:250;
0107 nr_list = [1,5:5:45];
0108 Dr0_list = [0,0.001,0.01,0.1,0.3,1,3,10,30,100];
0109 %
0110 % Incomplete LU factorization
0111 %
0112 dkeparam.invproc  = 1;%Numerical method (CGS, BICG, QMR)
0113 transpfaste.Dr0 = 0;
0114 %
0115 nnp = length(np_list);
0116 memory_LUi_np = NaN(1,nnp);
0117 time_LUi_np = NaN(1,nnp);
0118 %
0119 for inp = 1:nnp,
0120     %
0121     dkeparam.np_S = np_list(inp) + 1;
0122     %
0123     [Znorm,Zcurr,ZP0,dke_out] = main_dke_yp(id_simul,dkepath,equil,dkeparam,dkedisplay,ohm,waves,transpfaste,ripple,[],[]);
0124     %
0125     memory_LUi_np(inp) = dke_out.memoryLU_f_out;
0126     time_LUi_np(inp) = dke_out.time_DKE;
0127     %
0128 end
0129 %
0130 dkeparam.np_S = 101;
0131 %
0132 nnmhu = length(nmhu_list);
0133 memory_LUi_nmhu = NaN(1,nnmhu);
0134 time_LUi_nmhu = NaN(1,nnmhu);
0135 %
0136 for inmhu = 1:nnmhu,
0137     %
0138     dkeparam.nmhu_S = nmhu_list(inmhu) + 1;
0139     %
0140     [Znorm,Zcurr,ZP0,dke_out] = main_dke_yp(id_simul,dkepath,equil,dkeparam,dkedisplay,ohm,waves,transpfaste,ripple,[],[]);
0141     %
0142     memory_LUi_nmhu(inmhu) = dke_out.memoryLU_f_out;
0143     time_LUi_nmhu(inmhu) = dke_out.time_DKE;
0144     %
0145 end
0146 %
0147 dkeparam.nmhu_S = 101;
0148 %
0149 nnr = length(nr_list);
0150 memory_LUi_nr = NaN(1,nnr);
0151 time_LUi_nr = NaN(1,nnr);
0152 %
0153 for inr = 1:nnr,
0154     %
0155     rho_S = linspace(0,1,nr_list(inr)+1);
0156     %
0157     [Znorm,Zcurr,ZP0,dke_out] = main_dke_yp(id_simul,dkepath,equil,dkeparam,dkedisplay,ohm,waves,transpfaste,ripple,[],[]);
0158     %
0159     memory_LUi_nr(inr) = dke_out.memoryLU_f_out;
0160     time_LUi_nr(inr) = dke_out.time_DKE;
0161     %
0162 end
0163 %
0164 rho_S = 0:0.05:1;
0165 %
0166 nDr0 = length(Dr0_list);
0167 memory_LUi_Dr0 = NaN(1,nDr0);
0168 time_LUi_Dr0 = NaN(1,nDr0);
0169 %
0170 for iDr0 = 1:nDr0,
0171     %
0172     transpfaste.Dr0 = Dr0_list(iDr0);
0173     %
0174     [Znorm_MATLAB{iDr0},Zcurr_MATLAB{iDr0},ZP0,dke_out] = main_dke_yp(id_simul,dkepath,equil,dkeparam,dkedisplay,ohm,waves,transpfaste,ripple,[],[]);
0175     %
0176     memory_LUi_Dr0(iDr0) = dke_out.memoryLU_f_out;
0177     time_LUi_Dr0(iDr0) = dke_out.time_DKE;
0178     %
0179 end
0180 %
0181 transpfaste.Dr0 =  0;
0182 %
0183 disp('Incomplete LU calculation done')
0184 %
0185 % MUMPSMEX package
0186 %
0187 if exist('dmumpsmex'),
0188     dkeparam.invproc = -2;
0189     %
0190     nnp = length(np_list);
0191     memory_MX_np = NaN(1,nnp);
0192     time_MX_np = NaN(1,nnp);
0193     %
0194     for inp = 1:nnp,
0195         %
0196         disp(['inp = ',int2str(inp),'/',int2str(nnp)])
0197         %
0198         dkeparam.np_S = np_list(inp) + 1;
0199         %
0200         [Znorm,Zcurr,ZP0,dke_out] = main_dke_yp(id_simul,dkepath,equil,dkeparam,dkedisplay,ohm,waves,transpfaste,ripple,[],[]);
0201         %
0202         memory_MX_np(inp) = dke_out.memoryLU_f_out;
0203         time_MX_np(inp) = dke_out.time_DKE;
0204         %
0205     end
0206     %
0207     dkeparam.np_S = 101;
0208     %
0209     nnmhu = length(nmhu_list);
0210     memory_MX_nmhu = NaN(1,nnmhu);
0211     time_MX_nmhu = NaN(1,nnmhu);
0212     %
0213     for inmhu = 1:nnmhu,
0214         %
0215         disp(['inmhu = ',int2str(inmhu),'/',int2str(nnmhu)])
0216         %
0217         dkeparam.nmhu_S = nmhu_list(inmhu) + 1;
0218         %
0219         [Znorm,Zcurr,ZP0,dke_out] = main_dke_yp(id_simul,dkepath,equil,dkeparam,dkedisplay,ohm,waves,transpfaste,ripple,[],[]);
0220         %
0221         memory_MX_nmhu(inmhu) = dke_out.memoryLU_f_out;
0222         time_MX_nmhu(inmhu) = dke_out.time_DKE;
0223         %
0224     end
0225     %
0226     dkeparam.nmhu_S = 101;
0227     %
0228     nnr = length(nr_list);
0229     memory_MX_nr = NaN(1,nnr);
0230     time_MX_nr = NaN(1,nnr);
0231     %
0232     for inr = 1:nnr,
0233         %
0234         disp(['inr = ',int2str(inr),'/',int2str(nnr)])
0235         %
0236         rho_S = linspace(0,1,nr_list(inr)+1);
0237         %
0238         [Znorm,Zcurr,ZP0,dke_out] = main_dke_yp(id_simul,dkepath,equil,dkeparam,dkedisplay,ohm,waves,transpfaste,ripple,[],[]);
0239         %
0240         memory_MX_nr(inr) = dke_out.memoryLU_f_out;
0241         time_MX_nr(inr) = dke_out.time_DKE;
0242         %
0243     end
0244     %
0245     rho_S = 0:0.05:1;
0246     %
0247     nDr0 = length(Dr0_list);
0248     memory_MX_Dr0 = NaN(1,nDr0);
0249     time_MX_Dr0 = NaN(1,nDr0);
0250     %
0251     for iDr0 = 1:nDr0,
0252         %
0253         disp(['iDr0 = ',int2str(iDr0),'/',int2str(nDr0)])
0254         %
0255         transpfaste.Dr0 = Dr0_list(iDr0);
0256         %
0257         [Znorm_MUMPS{iDr0},Zcurr_MUMPS{iDr0},ZP0,dke_out] = main_dke_yp(id_simul,dkepath,equil,dkeparam,dkedisplay,ohm,waves,transpfaste,ripple,[],[]);
0258         %
0259         memory_MX_Dr0(iDr0) = dke_out.memoryLU_f_out;
0260         time_MX_Dr0(iDr0) = dke_out.time_DKE;
0261         %
0262     end
0263     %
0264     transpfaste.Dr0 =  0;
0265     %
0266     disp(['MUMPSMEX calculation done'])
0267 end
0268 %
0269 ratio_LUi_np = time_LUi_np./memory_LUi_np;
0270 ratio_LUi_nmhu = time_LUi_nmhu./memory_LUi_nmhu;
0271 ratio_LUi_nr = time_LUi_nr./memory_LUi_nr;
0272 ratio_LUi_Dr0 = time_LUi_Dr0./memory_LUi_Dr0;
0273 ratio_MX_np = time_MX_np./memory_MX_np;
0274 ratio_MX_nmhu = time_MX_nmhu./memory_MX_nmhu;
0275 ratio_MX_nr = time_MX_nr./memory_MX_nr;
0276 ratio_MX_Dr0 = time_MX_Dr0./memory_MX_Dr0;
0277 %
0278 %************************************************************************************************************************************
0279 %
0280 figure(1),clf
0281 %
0282 leg = {'MATLAB LUi','MUMPSMEX'};
0283 xlim = [0,300];
0284 ylim = [0,50];
0285 xlab = 'n_p';
0286 ylab = 'FPE calculation time (s)';
0287 tit = '';
0288 siz = 20+14i;
0289 %
0290 graph1D_jd(np_list,time_LUi_np,0,0,xlab,ylab,tit,NaN,xlim,ylim,'-','+','r',2,siz,gca,0.9,0.7,0.7);
0291 graph1D_jd(np_list,time_MX_np,0,0,'','','',leg,xlim,ylim,'-','+','b',2,siz,gca);
0292 %
0293 set(gca,'ytick',[0:0.2:1]*ylim(2))
0294 %set(gca,'xtick',[1e-7,1e-6,1e-5,1e-4,1e-3,1e-2,1e-1])
0295 set(gca,'xtick',[0:0.2:1]*xlim(2))
0296 %set(gca,'XTickLabel',{'10^{-5}','10^{-4}','10^{-3}','10^{-2}','10^{-1}','1'})
0297 %set(gca,'XMinorGrid','off')
0298 %set(gca,'XMinorTick','on')
0299 %
0300 figure(2),clf
0301 %
0302 xlim = [0,300];
0303 ylim = [0,200];
0304 xlab = 'n_{\xi}';
0305 %
0306 graph1D_jd(nmhu_list,time_LUi_nmhu,0,0,xlab,ylab,tit,NaN,xlim,ylim,'-','+','r',2,siz,gca,0.9,0.7,0.7);
0307 graph1D_jd(nmhu_list,time_MX_nmhu,0,0,'','','',leg,xlim,ylim,'-','+','b',2,siz,gca);
0308 set(gca,'ytick',[0:0.2:1]*ylim(2))
0309 %set(gca,'xtick',[1e-7,1e-6,1e-5,1e-4,1e-3,1e-2,1e-1])
0310 set(gca,'xtick',[0:0.2:1]*xlim(2))
0311 %set(gca,'XTickLabel',{'10^{-5}','10^{-4}','10^{-3}','10^{-2}','10^{-1}','1'})
0312 %set(gca,'XMinorGrid','off')
0313 %set(gca,'XMinorTick','on')
0314 %
0315 figure(3),clf
0316 %
0317 xlim = [0,30];
0318 ylim = [0,50];
0319 xlab = 'n_r';
0320 %
0321 graph1D_jd(nr_list,time_LUi_nr,0,0,xlab,ylab,tit,NaN,xlim,ylim,'-','+','r',2,siz,gca,0.9,0.7,0.7);
0322 graph1D_jd(nr_list,time_MX_nr,0,0,'','','',leg,xlim,ylim,'-','+','b',2,siz,gca);
0323 set(gca,'ytick',[0:0.2:1]*ylim(2))
0324 %set(gca,'xtick',[1e-7,1e-6,1e-5,1e-4,1e-3,1e-2,1e-1])
0325 set(gca,'xtick',[0:0.2:1]*xlim(2))
0326 %set(gca,'XTickLabel',{'10^{-5}','10^{-4}','10^{-3}','10^{-2}','10^{-1}','1'})
0327 %set(gca,'XMinorGrid','off')
0328 %set(gca,'XMinorTick','on')
0329 %
0330 figure(4),clf
0331 %
0332 xlim = [0.001,100];
0333 ylim = 10.^[0,4];%[0,5000];%
0334 xlab = 'D_{\psi\psi} (m^{2}/s)';
0335 %
0336 graph1D_jd(Dr0_list(2:end),time_LUi_Dr0(2:end),1,1,xlab,ylab,tit,NaN,xlim,ylim,'-','+','r',2,siz,gca,0.9,0.7,0.7);
0337 graph1D_jd(Dr0_list(2:end),time_MX_Dr0(2:end),1,1,'','','',leg,xlim,ylim,'-','+','b',2,siz,gca);
0338 graph1D_jd(xlim,[time_LUi_Dr0(1),time_LUi_Dr0(1)],1,1,'','','',leg,xlim,ylim,'--','none','r',2,siz,gca);
0339 graph1D_jd(xlim,[time_MX_Dr0(1),time_MX_Dr0(1)],1,1,'','','',leg,xlim,ylim,'--','none','b',2,siz,gca);
0340 set(gca,'ytick',10.^[0:4]);%set(gca,'ytick',[0:0.2:1]*ylim(2))%
0341 set(gca,'xtick',[0.001,0.01,0.1,1,10,100])
0342 %set(gca,'xtick',[0:0.2:1]*xlim(2))
0343 %set(gca,'XTickLabel',{'10^{-5}','10^{-4}','10^{-3}','10^{-2}','10^{-1}','1'})
0344 set(gca,'XMinorGrid','off')
0345 set(gca,'XMinorTick','on')
0346 %
0347 %
0348 figure(5),clf
0349 %
0350 xlim = [0,300];
0351 ylim = [0,1000];
0352 xlab = 'n_p';
0353 ylab = 'LU matrix size (MBytes)';
0354 %
0355 graph1D_jd(np_list,memory_LUi_np,0,0,xlab,ylab,tit,NaN,xlim,ylim,'-','+','r',2,siz,gca,0.9,0.7,0.7);
0356 graph1D_jd(np_list,memory_MX_np,0,0,'','','',leg,xlim,ylim,'-','+','b',2,siz,gca);
0357 %set(gca,'ytick',[0:0.2:1]*ylim(2))
0358 %set(gca,'xtick',[1e-7,1e-6,1e-5,1e-4,1e-3,1e-2,1e-1])
0359 set(gca,'xtick',[0:0.2:1]*xlim(2))
0360 %set(gca,'XTickLabel',{'10^{-5}','10^{-4}','10^{-3}','10^{-2}','10^{-1}','1'})
0361 %set(gca,'XMinorGrid','off')
0362 %set(gca,'XMinorTick','on')
0363 %
0364 figure(6),clf
0365 %
0366 xlim = [0,300];
0367 xlab = 'n_{\xi}';
0368 %
0369 graph1D_jd(nmhu_list,memory_LUi_nmhu,0,0,xlab,ylab,tit,NaN,xlim,ylim,'-','+','r',2,siz,gca,0.9,0.7,0.7);
0370 graph1D_jd(nmhu_list,memory_MX_nmhu,0,0,'','','',leg,xlim,ylim,'-','+','b',2,siz,gca);
0371 %set(gca,'ytick',[0:0.2:1]*ylim(2))
0372 %set(gca,'xtick',[1e-7,1e-6,1e-5,1e-4,1e-3,1e-2,1e-1])
0373 set(gca,'xtick',[0:0.2:1]*xlim(2))
0374 %set(gca,'XTickLabel',{'10^{-5}','10^{-4}','10^{-3}','10^{-2}','10^{-1}','1'})
0375 %set(gca,'XMinorGrid','off')
0376 %set(gca,'XMinorTick','on')
0377 %
0378 figure(7),clf
0379 %
0380 xlim = [0,30];
0381 xlab = 'n_r';
0382 %
0383 graph1D_jd(nr_list,memory_LUi_nr,0,0,xlab,ylab,tit,NaN,xlim,ylim,'-','+','r',2,siz,gca,0.9,0.7,0.7);
0384 graph1D_jd(nr_list,memory_MX_nr,0,0,'','','',leg,xlim,ylim,'-','+','b',2,siz,gca);
0385 %set(gca,'ytick',[0:0.2:1]*ylim(2))
0386 %set(gca,'xtick',[1e-7,1e-6,1e-5,1e-4,1e-3,1e-2,1e-1])
0387 set(gca,'xtick',[0:0.2:1]*xlim(2))
0388 %set(gca,'XTickLabel',{'10^{-5}','10^{-4}','10^{-3}','10^{-2}','10^{-1}','1'})
0389 %set(gca,'XMinorGrid','off')
0390 %set(gca,'XMinorTick','on')
0391 %
0392 figure(8),clf
0393 %
0394 xlim = [0.001,100];
0395 ylim = 10.^[0,5];%[0,1000];%
0396 xlab = 'D_{\psi\psi} (m^{2}/s)';
0397 %
0398 graph1D_jd(Dr0_list(2:end),memory_LUi_Dr0(2:end),1,1,xlab,ylab,tit,NaN,xlim,ylim,'-','+','r',2,siz,gca,0.9,0.7,0.7);
0399 graph1D_jd(Dr0_list(2:end),memory_MX_Dr0(2:end),1,1,'','','',leg,xlim,ylim,'-','+','b',2,siz,gca);
0400 graph1D_jd(xlim,[memory_LUi_Dr0(1),memory_LUi_Dr0(1)],1,1,'','','',leg,xlim,ylim,'--','none','r',2,siz,gca);
0401 graph1D_jd(xlim,[memory_MX_Dr0(1),memory_MX_Dr0(1)],1,1,'','','',leg,xlim,ylim,'--','none','b',2,siz,gca);
0402 set(gca,'ytick',10.^[0:5]);%set(gca,'ytick',[0:0.2:1]*ylim(2));%
0403 set(gca,'xtick',[0.001,0.01,0.1,1,10,100])
0404 %set(gca,'xtick',[0:0.2:1]*xlim(2))
0405 %set(gca,'XTickLabel',{'10^{-5}','10^{-4}','10^{-3}','10^{-2}','10^{-1}','1'})
0406 set(gca,'XMinorGrid','off')
0407 set(gca,'XMinorTick','on')
0408 %
0409 figure(9),clf
0410 %
0411 xlim = 10.^[-2,4];%[0.1,2000];%
0412 ylim = 10.^[-2,1];%[0.05,5];%
0413 xlab = 'FPE time (s)';
0414 ylab = 'FPE time/LU size  (s/MBytes)';
0415 %
0416 graph1D_jd(time_LUi_np,ratio_LUi_np,1,1,xlab,ylab,tit,NaN,xlim,ylim,'none','+','r',2,siz,gca,0.9,0.7,0.7);
0417 graph1D_jd(time_MX_np,ratio_MX_np,1,1,'','','',leg,xlim,ylim,'none','o','r',2,siz,gca);
0418 graph1D_jd(time_LUi_nmhu,ratio_LUi_nmhu,1,1,'','','',NaN,xlim,ylim,'none','+','b',2,siz,gca);
0419 graph1D_jd(time_LUi_nr,ratio_LUi_nr,1,1,'','','',NaN,xlim,ylim,'none','+','g',2,siz,gca);
0420 graph1D_jd(time_LUi_Dr0,ratio_LUi_Dr0,1,1,'','','',NaN,xlim,ylim,'none','+','m',2,siz,gca);
0421 graph1D_jd(time_MX_nmhu,ratio_MX_nmhu,1,1,'','','',NaN,xlim,ylim,'none','o','b',2,siz,gca);
0422 graph1D_jd(time_MX_nr,ratio_MX_nr,1,1,'','','',NaN,xlim,ylim,'none','o','g',2,siz,gca);
0423 graph1D_jd(time_MX_Dr0,ratio_MX_Dr0,1,1,'','','',NaN,xlim,ylim,'none','o','m',2,siz,gca);
0424 %set(gca,'ytick',[0:0.2:1]*ylim(2))
0425 set(gca,'ytick',10.^[-2:1]);%set(gca,'ytick',[1e-1,1]);%
0426 set(gca,'xtick',10.^[-2:4]);%set(gca,'xtick',[1e-1,1,10,100,1000]);%
0427 %set(gca,'xtick',[0:0.2:1]*xlim(2))
0428 %set(gca,'XTickLabel',{'10^{-5}','10^{-4}','10^{-3}','10^{-2}','10^{-1}','1'})
0429 set(gca,'XMinorGrid','off')
0430 set(gca,'XMinorTick','on')
0431 set(gca,'YMinorGrid','off')
0432 set(gca,'YMinorTick','on')
0433 %
0434 %
0435 print_jd(2,'fig_num_3D_np_time','./figures',1)
0436 print_jd(2,'fig_num_3D_nmhu_time','./figures',2)
0437 print_jd(2,'fig_num_3D_nr_time','./figures',3)
0438 print_jd(2,'fig_num_3D_Dr0_time','./figures',4)
0439 print_jd(2,'fig_num_3D_np_memory','./figures',5)
0440 print_jd(2,'fig_num_3D_nmhu_memory','./figures',6)
0441 print_jd(2,'fig_num_3D_nr_memory','./figures',7)
0442 print_jd(2,'fig_num_3D_Dr0_memory','./figures',8)
0443 print_jd(2,'fig_num_3D_ratio','./figures',9)
0444 %
0445 %************************************************************************************************************************************
0446 eval(['save ',path_simul,'DKE_RESULTS_',id_equil,'_',id_simul,'.mat']);
0447 info_dke_yp(2,['Data saved in ',path_simul,'DKE_RESULTS_',id_equil,'_',id_simul,'.mat']);

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