ESPUTR
calculation of sputtering yields
Data Types | Modules | Functions/Subroutines
esputr_interfaces.f90 File Reference

Wrappers which are used to call esputr from C/C++, as well as in Matlab and Octave interfaces. More...

Go to the source code of this file.

Data Types

interface  esputr_interfaces::esputr_wrapper
 

Modules

module  esputr_interfaces
 

Functions/Subroutines

subroutine esputr1993_init_wrapper (fileN, err)
 See include esputr_octave.hpp, esputr_YN.cpp, esputr_YTH.cpp as an example of how to use those subroutines in a C program. More...
 
subroutine esputr2001_init_wrapper (fileN, fileTH, err)
 module-less wraper for esputr2001_initN More...
 
subroutine esputr2001_initn_wrapper (fileN, err)
 module-less wraper for esputr2001_initN More...
 
subroutine esputr2001_initth_wrapper (fileTH, err)
 module-less wrapper for esputr2001_initTH More...
 
subroutine esputr_deallocate (err)
 call esputr2001_deallocate and esputr1993_deallocate More...
 
subroutine esputr_yn (E0s, n_E0s, proj, targ, version, YNs, err)
 Calculate sputtering yield for normal incidence for selected incident energies and the specified projectile-target combination. More...
 
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 for the specified projectile-target combination. More...
 
subroutine esputr2001_getavailableenergyrangebyname (proj, targ, Emin, Emax, err)
 Return the minimum and maximum energies for which the angular dependence data are defined. More...
 

Detailed Description

Wrappers which are used to call esputr from C/C++, as well as in Matlab and Octave interfaces.

Authors
Markus Brenneis Marku.nosp@m.s.Br.nosp@m.ennei.nosp@m.s@un.nosp@m.i-due.nosp@m.ssel.nosp@m.dorf..nosp@m.de, Vladislav Kotov v.kot.nosp@m.ov@f.nosp@m.z-jue.nosp@m.lich.nosp@m..de (contact person)

Definition in file esputr_interfaces.f90.

Function/Subroutine Documentation

subroutine esputr1993_init_wrapper ( character(*), intent(in)  fileN,
integer, intent(out)  err 
)

See include esputr_octave.hpp, esputr_YN.cpp, esputr_YTH.cpp as an example of how to use those subroutines in a C program.

Note
esputr_yn and esputr_yth are defined as subroutines rather than functions because in Octave the use of functions requires some extra work. Quotation from the Octave documentation: Note that we don't attempt to handle Fortran functions, we always use subroutine wrappers for them and pass the return value as an extra argument. The problem is that the results are written to a double* resultVar.fortran_vec(), but a function would return a double on the stack, so the return value gets invalidated when leaving the function. On could work around this problem by explicitly copying the result to the octave variable, but this is more complicated than using subroutines. module-less wrapper for esputr1993_init

Definition at line 43 of file esputr_interfaces.f90.

References esputr1993::esputr1993_init().

Referenced by DEFUN_DLD().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine esputr2001_init_wrapper ( character(*), intent(in)  fileN,
character(*), intent(in)  fileTH,
integer, intent(out)  err 
)

module-less wraper for esputr2001_initN

Definition at line 52 of file esputr_interfaces.f90.

References esputr2001::esputr2001_init().

Here is the call graph for this function:

subroutine esputr2001_initn_wrapper ( character(*), intent(in)  fileN,
integer, intent(out)  err 
)

module-less wraper for esputr2001_initN

Definition at line 61 of file esputr_interfaces.f90.

References esputr2001::esputr2001_initn().

Referenced by DEFUN_DLD().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine esputr2001_initth_wrapper ( character(*), intent(in)  fileTH,
integer, intent(out)  err 
)

module-less wrapper for esputr2001_initTH

Definition at line 70 of file esputr_interfaces.f90.

References esputr2001::esputr2001_initth().

Referenced by DEFUN_DLD().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine esputr_deallocate ( integer, intent(out)  err)

call esputr2001_deallocate and esputr1993_deallocate

Definition at line 79 of file esputr_interfaces.f90.

References esputr1993::esputr1993_deallocate(), and esputr2001::esputr2001_deallocate().

Referenced by DEFUN_DLD().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine esputr_yn ( real(kind=esputr_dp), dimension(n_e0s), intent(in)  E0s,
integer, intent(in)  n_E0s,
character(*), intent(in)  proj,
character(*), intent(in)  targ,
character(*), intent(in)  version,
real(kind=esputr_dp), dimension(size(e0s)), intent(out)  YNs,
integer, intent(out)  err 
)

Calculate sputtering yield for normal incidence for selected incident energies and the specified projectile-target combination.

Warning
In case of an error this subroutine deallocates the modules
Parameters
[in]e0sArray of incident energies, eV
[in]n_e0sNumber of energies in E0s
[in]projString which defines the projectile (H, He etc.)
[in]targString which defines the target (Be, W, Tungsten etc.)
[in]versionVersion of the model ('1993', '2001' etc.)
[out]errError code
[out]ynsSputtering yields for normal incidence

Definition at line 125 of file esputr_interfaces.f90.

References esputr1993::esputr1993_deallocate(), esputr1993::esputr1993_getprojectileid(), esputr1993::esputr1993_gettargetid(), esputr1993::esputr1993_yn(), esputr2001::esputr2001_deallocate(), esputr2001::esputr2001_getprojectiletargetidn(), esputr2001::esputr2001_yn(), and esputr::esputr_unit.

Referenced by DEFUN_DLD(), esputr_cgi(), and mexfunction().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine esputr_yth ( real(kind=esputr_dp), dimension(n_ths), intent(in)  THs,
integer, intent(in)  n_THs,
real(kind=esputr_dp), dimension(n_e0s), intent(in)  E0s,
integer, intent(in)  n_E0s,
character(*), intent(in)  proj,
character(*), intent(in)  targ,
character(*), intent(in)  version,
real(kind=esputr_dp), dimension(size(ths), size(e0s)), intent(out)  YTHs,
integer, intent(out)  err 
)

Calculate the angular dependence factor Y(E,theta)/Y(E,0) for selected incident angles and energies for the specified projectile-target combination.

Warning
In case of an error this subroutine deallocates the modules
Parameters
[in]thsArray of incident angles, radian
[in]n_thsNumber of elements in THs
[in]e0sArray of incident energies, eV
[in]n_e0sNumber of elements in E0s
[in]projString which defines the projectile (H, He etc.)
[in]targString which defines the target (Be, W, Tungsten etc.)
[in]versionVersion of the model ('1993', '2001' etc.)
[out]errError code
[out]ythsAngular dependence factor Y(E,theta)/Y(E,0) for each angle-energy combination first index: angles; second index: energy

Definition at line 184 of file esputr_interfaces.f90.

References esputr1993::esputr1993_deallocate(), esputr1993::esputr1993_yth(), esputr2001::esputr2001_deallocate(), esputr2001::esputr2001_getprojectiletargetidsth(), esputr2001::esputr2001_yth(), and esputr::esputr_unit.

Referenced by DEFUN_DLD(), esputr_cgi(), and mexfunction().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine esputr2001_getavailableenergyrangebyname ( character(*), intent(in)  proj,
character(*), intent(in)  targ,
real(kind=esputr_dp), intent(out)  Emin,
real(kind=esputr_dp), intent(out)  Emax,
integer, intent(out)  err 
)

Return the minimum and maximum energies for which the angular dependence data are defined.

Parameters
[in]projName of the projectile
[in]targName of the target
[out]eminMinimal energy
[out]emaxMaximal energy

Definition at line 247 of file esputr_interfaces.f90.

References esputr2001::esputr2001_getavailableenergyrange(), and esputr2001::esputr2001_getprojectiletargetidsth().

Referenced by mexfunction().

Here is the call graph for this function:

Here is the caller graph for this function: