ESPUTR
calculation of sputtering yields
esputr_YTH.cpp
Go to the documentation of this file.
1 
6 // Copyright (c) 2016 Forschungszentrum Juelich GmbH
7 // Markus Brenneis, Vladislav Kotov
8 //
9 // This file is part of ESPUTR.
10 //
11 // ESPUTR is free software: you can redistribute it and/or modify
12 // it under the terms of the GNU General Public License as published by
13 // the Free Software Foundation, either version 3 of the License, or
14 // (at your option) any later version.
15 //
16 // ESPUTR is distributed in the hope that it will be useful,
17 // but WITHOUT ANY WARRANTY; without even the implied warranty of
18 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19 // GNU General Public License for more details.
20 //
21 // You should have received a copy of the GNU General Public License
22 // along with ESPUTR. If not, see <http://www.gnu.org/licenses/>.
23 //
24 #include <octave/oct.h>
25 #include <cstring>
26 #include <string>
27 #include "esputr_octave.hpp"
28 
29 // Copyright (c) 2016 Forschungszentrum Juelich GmbH
30 // Markus Brenneis, Vladislav Kotov
31 //
32 // This file is part of ESPUTR.
33 //
34 // ESPUTR is free software: you can redistribute it and/or modify
35 // it under the terms of the GNU General Public License as published by
36 // the Free Software Foundation, either version 3 of the License, or
37 // (at your option) any later version.
38 //
39 // ESPUTR is distributed in the hope that it will be useful,
40 // but WITHOUT ANY WARRANTY; without even the implied warranty of
41 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
42 // GNU General Public License for more details.
43 //
44 // You should have received a copy of the GNU General Public License
45 // along with ESPUTR. If not, see <http://www.gnu.org/licenses/>.
46 //
47 
48 DEFUN_DLD (esputr_YTH, args, nargout,
49 " \n\
50  Calculate the angular dependency factor of the sputtering yield Y(E,tehta)/Y(E,0) (interface to ESPUTR) \n\
51 \n\
52  [Y1 Y2 ...] = esputr_YN(TH, E0, proj, targ, version, fileTH) \n\
53 \n\
54  Y1, Y2 ... : matrices with calculated angular factors - \n\
55  - one matrix for each projectile-target combination, \n\
56  rows: over angle THs, columns: over eneregy E0s \n\
57  TH : column vector of the incident agles, radian \n\
58  Incident angle is the angle between surface normal \n\
59  and velocity of the incident particle \n\
60  E0 : row vector of the projectile energies, eV \n\
61  proj : column vector of the projectile names \n\
62  targ : column vector of the target names \n\
63  version : designation of the model ('1993' or '2001') \n\
64  fileN : name of the file (full path) with fitting parameters \n\
65 \n\
66  If the lengths of proj and targ differ, the last element of the shorter vector is repeated \n\
67 \n\
68  Examples:\n\
69 \n\
70  1993-model\n\
71  f = esputr_YTH([0:5:90]*pi/180, 100, 'E', 'l', '1993')\n\
72 \n\
73  D on Be, 2007 model (2001 with coefficients from [Eckstein 2007]) \n\
74  f = esputr_YTH([0:5:90]*pi/180, [500,700,1000], 'D', 'Be', '2001','../data/ECKSTEIN2007TH')\n\
75 \n\
76  D on W, D on Mo, He on W, 2007 model\n\
77  [D_W,D_Mo,He_W] = esputr_YTH([0:1:90]*pi/180, [500,1000], strvcat('D','D','He'),strvcat('W','Mo','W'), '2001', ... \n\
78  '/home/v.kotov/JuNiPlib/esputr/data/ECKSTEIN2007TH')") {
79  octave_value_list retval;
80  int nargin = args.length();
81 
82  if(nargin < 5) {
83  error("at least 5 input arguments required");
84 //vk print_usage();
85  } else {
86  NDArray THs = args(0).array_value();
87  NDArray E0s = args(1).array_value();
88  if(!error_state) {
89  double* THs_v = THs.fortran_vec();
90  octave_idx_type n_THs = THs.numel();
91 
92  double* E0s_v = E0s.fortran_vec();
93  octave_idx_type n_E0s = E0s.numel();
94 
95  charMatrix projs = args(2).char_matrix_value();
96  octave_idx_type n_projs = projs.rows();
97  charMatrix targs = args(3).char_matrix_value();
98  octave_idx_type n_targs = targs.rows();
99 
100  octave_idx_type n_combinations = std::max(n_E0s, std::max(n_projs, n_targs));
101 
102  std::string version_str = args(4).char_matrix_value().row_as_string(0);
103  const char* version = version_str.c_str();
104  octave_idx_type version_length = version_str.length();
105 
106  int err = 0;
107 
108  std::string fileTH = "";
109 
110  if(version_str != "2001" && version_str !="1993") {
111  error("unknown model");
112  F77_XFCN(esputr_deallocate, FORTSUB, (err));
113  return retval;
114  }
115 
116  if(version_str == "2001" && nargin == 6) {
117  fileTH = args(5).char_matrix_value().row_as_string(0);
118  } else if(version_str == "2001" && nargin < 6) {
119  error("6 arguments required");
120  F77_XFCN(esputr_deallocate, FORTSUB, (err));
121  return retval;
122  }
123 
124  if(error_state) {
125  error("sth went wrong");
126  return retval;
127  }
128 
129  if(version_str == "2001") {
130  F77_XFCN(esputr2001_initth_wrapper, FORTSUB, (fileTH.c_str(), err, fileTH.length()));
131  if(err != 0) {
132  error("error in esputr_matlab_init");
133  F77_XFCN(esputr_deallocate, FORTSUB, (err));
134  return retval;
135  }
136  }
137 
138  for(int c = 0; c < n_combinations; c++) {
139  std::string proj = projs.row_as_string(std::min(c, n_projs-1));
140  std::string targ = targs.row_as_string(std::min(c, n_targs-1));
141  Matrix YNs (n_THs, n_E0s);
142  F77_XFCN(esputr_yth, FORTSUB,
143  (THs_v, n_THs,
144  E0s_v, n_E0s,
145  proj.c_str(),
146  targ.c_str(),
147  version,
148  YNs.fortran_vec(),
149  err
150  F77_CHAR_ARG_LEN(proj.length())
151  F77_CHAR_ARG_LEN(targ.length())
152  F77_CHAR_ARG_LEN(version_length)
153  ));
154  if(err != 0) {
155  error("error in esputr_matlab_yth");
156  F77_XFCN(esputr_deallocate, FORTSUB, (err));
157  return retval;
158  }
159  retval(c) = YNs;
160  }
161 
162  F77_XFCN(esputr_deallocate, FORTSUB, (err));
163  }
164  }
165  return retval;
166 }
subroutine esputr_deallocate(err)
call esputr2001_deallocate and esputr1993_deallocate
subroutine esputr2001_initth_wrapper(fileTH, err)
module-less wrapper for esputr2001_initTH
integer, parameter n_combinations
Number of tested projectile-target combinations.
subroutine esputr_yth(THs, n_THs, E0s, n_E0s, proj, targ, version, YTHs, err)
Calculate the angular dependence factor Y(E,theta)/Y(E,0) for selected incident angles and energies f...
DEFUN_DLD(esputr_YTH, args, nargout," \n\ Calculate the angular dependency factor of the sputtering yield Y(E,tehta)/Y(E,0) (interface to ESPUTR) \n\ \n\ [Y1 Y2 ...] = esputr_YN(TH, E0, proj, targ, version, fileTH) \n\ \n\ Y1, Y2 ... : matrices with calculated angular factors - \n\ - one matrix for each projectile-target combination, \n\ rows: over angle THs, columns: over eneregy E0s \n\ TH : column vector of the incident agles, radian \n\ Incident angle is the angle between surface normal \n\ and velocity of the incident particle \n\ E0 : row vector of the projectile energies, eV \n\ proj : column vector of the projectile names \n\ targ : column vector of the target names \n\ version : designation of the model ('1993' or '2001') \n\ fileN : name of the file (full path) with fitting parameters \n\ \n\ If the lengths of proj and targ differ, the last element of the shorter vector is repeated \n\ \n\ Examples:\n\ \n\ 1993-model\n\ f = esputr_YTH([0:5:90]*pi/180, 100, 'E', 'l', '1993')\n\ \n\ D on Be, 2007 model (2001 with coefficients from [Eckstein 2007]) \n\ f = esputr_YTH([0:5:90]*pi/180, [500,700,1000], 'D', 'Be', '2001','../data/ECKSTEIN2007TH')\n\ \n\ D on W, D on Mo, He on W, 2007 model\n\ [D_W,D_Mo,He_W] = esputr_YTH([0:1:90]*pi/180, [500,1000], strvcat('D','D','He'),strvcat('W','Mo','W'), '2001', ... \n\ '/home/v.kotov/JuNiPlib/esputr/data/ECKSTEIN2007TH')")
Definition: esputr_YTH.cpp:48