ESPUTR
calculation of sputtering yields
esputr_cgi.f90
Go to the documentation of this file.
1 
55 
56 ! Copyright (c) 2016 Forschungszentrum Juelich GmbH
57 ! Markus Brenneis, Vladislav Kotov
58 !
59 ! This file is part of ESPUTR.
60 !
61 ! ESPUTR is free software: you can redistribute it and/or modify
62 ! it under the terms of the GNU General Public License as published by
63 ! the Free Software Foundation, either version 3 of the License, or
64 ! (at your option) any later version.
65 !
66 ! ESPUTR is distributed in the hope that it will be useful,
67 ! but WITHOUT ANY WARRANTY; without even the implied warranty of
68 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
69 ! GNU General Public License for more details.
70 !
71 ! You should have received a copy of the GNU General Public License
72 ! along with ESPUTR. If not, see <http://www.gnu.org/licenses/>.
73 !
74 
75 program esputr_cgi
76 
77  use cgi_protocol
78  use esputr
79  use esputr2001
80  use esputr1993
82  implicit none
84  integer :: unit
86  type(dict_struct), pointer :: dict => null()
88  character(*), parameter :: HEADER = "Content-type: text/csv; charset=utf-8"//char(10)//&
89  "Content-Disposition: attachment; filename=data.csv"//char(10)
91  integer, parameter :: MAX_POINTS = 1000000
92 
93  integer, parameter :: MAX_PARAM_LENGTH = 20
94  integer, parameter :: TMP_UNIT = 142
95  character(DICT_VALUE_LENGTH) :: input = ""
96  character(MAX_PARAM_LENGTH), dimension(:), allocatable :: inputs
97  character(MAX_PARAM_LENGTH), dimension(:), allocatable :: projs
98  integer :: n_projs
99  character(MAX_PARAM_LENGTH), dimension(:), allocatable :: targs
100  integer :: n_targs
101  character(MAX_PARAM_LENGTH), dimension(:), allocatable :: versions
102  integer :: n_versions = 0
103  integer :: n_combinations
104  integer :: err, i, c, p, t, v, th, e
105  real(kind=ESPUTR_DP) :: E0min = 10.
106  real(kind=ESPUTR_DP) :: E0max = 10000.
107  integer :: n_E0s = 1001
108  real(kind=ESPUTR_DP) :: deltaE0 = 10.
109  real(kind=ESPUTR_DP) :: thMin = 0
110  real(kind=ESPUTR_DP) :: thMax = 1.57
111  integer :: n_Ths = 0
112  real(kind=ESPUTR_DP) :: deltaTh = .1
113  real(kind=ESPUTR_DP), dimension(:), allocatable :: E0s
114  real(kind=ESPUTR_DP), dimension(:), allocatable :: Ths
115  real(kind=ESPUTR_DP), dimension(:,:), allocatable :: Ys
116  real(kind=ESPUTR_DP) :: epsfilter=0.
117  integer :: ilast
118  logical :: lskip
119  character(MAX_PARAM_LENGTH) :: ver
120 
121  intrinsic floor,transfer,count,max,index
122 
123  open(tmp_unit, iostat=err, status='scratch', action='readwrite')
124  esputr_unit = tmp_unit
125 
126  call get_environment_variable("QUERY_STRING", input)
127  if(index(input, '=') == 0) then
129  end if
130 
131  call cgi_begin(output_no_header, dict, unit)
132 
133  input = ""
134  call cgi_get(dict, "projectiles", input)
135  if(input == "") call cgi_error("esputr_cgi: no projectiles provided")
136  call split(input, projs)
137  n_projs = size(projs)
138 
139  input = ""
140  call cgi_get(dict, "targets", input)
141  if(input == "") call cgi_error("esputr_cgi: no targets provided")
142  call split(input, targs)
143  n_targs = size(targs)
144 
145  if ( dict_has_key( dict, "thdelta" ) .or. &
146  dict_has_key( dict, "thmin" ) .or. &
147  dict_has_key( dict, "thmax" ) ) then
148  call cgi_get(dict, "thdelta", deltath)
149  if(deltath < .001) deltath = .001d0
150  call cgi_get(dict, "thmin", thmin)
151  if(thmin < 0) thmin = 0
152  if(thmin > 1.57079d0) thmin = 1.57079d0
153  call cgi_get(dict, "thmax", thmax)
154  if(thmax > 1.57079d0) thmax = 1.57079d0
155  if(thmax < thmin) thmax = thmin + deltath
156  n_ths = floor((thmax-thmin) / deltath + 1.)
157  if(thmin >= 41) n_ths = 0
158  else
159  n_ths=0
160  end if
161 
162  input = ""
163  call cgi_get(dict, "versions", input)
164  if(input == "") call cgi_error("esputr_cgi: the models are not specified")
165  call split(input, versions)
166  n_versions = size(versions)
167 
168  call cgi_get(dict, "edelta", deltae0)
169  if(deltae0 < 0.) deltae0 = 0.
170  call cgi_get(dict, "emin", e0min)
171  if(e0min < 0.) e0min = 0.
172  call cgi_get(dict, "emax", e0max)
173  if(e0max < e0min) e0max = e0min + deltae0
174  call cgi_get(dict, "epsfilter", epsfilter)
175  input = ""
176  call cgi_get(dict, "e0s", input)
177  if(input /= "") then
178  call split(input, inputs)
179  call toreal8(inputs, e0s)
180  n_e0s = size(e0s)
181  else
182  if(deltae0>0.) then
183  n_e0s = floor( (e0max-e0min) / deltae0 + 1. )
184  else
185  n_e0s=1
186  end if
187  allocate(e0s(n_e0s), stat=err)
188  if(err /= 0) call cgi_error("esputr_cgi: memory allocation eror")
189  if (n_e0s>1) then
190  deltae0 = (e0max-e0min)/(n_e0s-1)
191  e0s(1)=e0min
192  do i=2,n_e0s-1
193  e0s(i)=e0s(i-1)+deltae0
194  end do
195  e0s(n_e0s)=e0max
196  else
197  e0s=e0min
198  end if
199  end if
200 
201  if(n_e0s*max(n_projs, n_targs, n_versions)*max(1, n_ths) > max_points) &
202  call cgi_error("esputr_cgi: too many data points are requested")
203 
204  allocate(ths(n_ths), stat=err)
205  if(err /= 0) call cgi_error("esputr_cgi: memory allocation eror")
206  ths = (/ (i*deltath+thmin, i = 0, n_ths-1) /)
207 
208  call initesputr()
209 
210  if(n_ths == 0) then
211  n_combinations = max(n_projs, n_targs, n_versions)
212  allocate(ys(n_combinations,n_e0s), stat=err)
213  if(err /= 0) call cgi_error("esputr_cgi: allocation error")
214 
215  foreachcombi: do c = 1, n_combinations
216  p = min(n_projs, c)
217  t = min(n_targs, c)
218  v = min(n_versions, c)
219  ver=versions(v)
220  if (index(ver,'2007')>0) ver='2001'
221  call esputr_yn(e0s, n_e0s, projs(p), targs(t), ver, ys(c,:), err)
222  if(err /= 0) call cgi_error("esputr_cgi: error while calculating the sputtering yield, check your input"&
223  //getlastmessage())
224  end do foreachcombi
225 
226  write(unit, '(A)') header
227 
228  write(unit, '(A)', advance='no') "Energy [eV]"
229  foreachcombi1: do c = 1, n_combinations
230  p = min(n_projs, c)
231  t = min(n_targs, c)
232  v = min(n_versions, c)
233  write(unit, '(6A)', advance='no') ", Y_",trim(projs(p)),"_on_",trim(targs(t)),"_",trim(versions(v))
234  end do foreachcombi1
235  write(unit, *) ""
236  ilast=1
237  do i = 1, n_e0s
238  !filter out nodes with almost equal (withing given accuracy epsfilter)
239  !sputtering yields
240  if(epsfilter>0.) then
241  if (i>1.and.i<n_e0s) then
242  lskip=.true.
243  do c = 1, n_combinations
244  if( abs(ys(c,i)-ys(c,ilast) ) > &
245  0.5*epsfilter*( abs(ys(c,i))+abs(ys(c,ilast)) ) ) then
246  lskip=.false.
247  exit
248  end if
249  end do
250  if(lskip) then
251  cycle
252  else
253  ilast=i
254  end if
255  end if !if (i>1.and.i<nE0s) then
256  end if !if(epsfilter>0.) then
257  write(unit, '(E12.5)', advance='no') e0s(i)
258  foreachcombi2: do c = 1, n_combinations
259  write(unit, '(A,E12.5)', advance='no') ",",ys(c,i)
260  end do foreachcombi2
261  write(unit, *) ""
262  end do
263  else
264  n_combinations = max(n_projs, n_targs, n_versions, n_e0s)
265  allocate(ys(n_combinations,n_ths), stat=err)
266  if(err /= 0) call cgi_error("esputr_cgi: allocation error")
267  foreachcombith: do c = 1, n_combinations
268  p = min(n_projs, c)
269  t = min(n_targs, c)
270  v = min(n_versions, c)
271  e = min(n_e0s, c)
272  ver=versions(v)
273  if (index(ver,'2007')>0) ver='2001'
274  call esputr_yth(ths, n_ths, (/e0s(e)/), 1, projs(p), targs(t), ver, ys(c,:), err)
275  if(err /= 0) call cgi_error("esputr_cgi: error while calculating the angular factor, check your input"//&
276  " (no data for this energy?)"//getlastmessage())
277  end do foreachcombith
278 
279  write(unit, '(A)') header
280 
281  write(unit, '(A)', advance='no') "theta [radian]"
282  foreachcombith1: do c = 1, n_combinations
283  p = min(n_projs, c)
284  t = min(n_targs, c)
285  v = min(n_versions, c)
286  e = min(n_e0s, c)
287  if (index(versions(v),'1993') > 0) then
288  write(unit, '(1A)', advance='no') ', Y(theta)/Y(0)_1993'
289  else
290  write(unit, '(5A,E12.5,2A)', advance='no') ", Y(theta)/Y(0)_",trim(projs(p)),"_on_",&
291  trim(targs(t)),"_E=_",e0s(e),"_eV_",trim(versions(v))
292  end if
293  end do foreachcombith1
294  write(unit, *) ""
295  do th = 1, n_ths
296  write(unit, '(E12.5)', advance='no') ths(th)
297  foreachcombith2: do c = 1, n_combinations
298  write(unit, '(A,E12.5)', advance='no') ",",ys(c, th)
299  end do foreachcombith2
300  write(unit, *) ""
301  end do
302  end if
303 
304  call esputr1993_deallocate(err)
305  call esputr2001_deallocate(err)
306  deallocate(e0s, ys, projs, targs, versions, ths)
307 
308  call cgi_end()
309 
310  contains
311 
313  subroutine split(csvString, array)
314  character(DICT_VALUE_LENGTH), intent(in) :: csvString
315  character(MAX_PARAM_LENGTH), dimension(:), allocatable, intent(inout) :: array
316  integer :: nCommas
317  ncommas = count(transfer(csvstring, 'a', dict_value_length) == ",")
318  allocate(array(ncommas + 1))
319  read(csvstring, *) array(1:ncommas+1)
320  end subroutine split
321 
323  subroutine toreal8(strArr, realArr)
324  character(MAX_PARAM_LENGTH), dimension(:), intent(in) :: strArr
325  real(kind=8), dimension(:), allocatable, intent(inout) :: realArr
326  integer :: i, err
327  allocate(realarr(size(strarr)))
328  do i = 1, size(strarr)
329  read(strarr(i), *, iostat=err) realarr(i)
330  if(err/=0) call cgi_error("esputr_cgi: invalid real value")
331  end do
332  end subroutine
333 
335  subroutine sendavailablecombinations()
336  character(ESPUTR2001_MAX_NAME_LENGTH), dimension(:,:), allocatable :: combinations_2001_N
337  character(ESPUTR2001_MAX_NAME_LENGTH*2), dimension(:,:), allocatable :: combinations_2001_TH
338  !TODO get rid of constants 11, 16
339  character(ESPUTR1993_FIRST_COL_WIDTH), dimension(:,:), allocatable :: combinations_1993
340  integer :: err
341  character(*), parameter :: JSON_HEADER = "Content-type: application/json; charset=utf-8"//char(10)//&
342  "Content-Disposition: attachment; filename=combinations.json"//char(10)
343 
344  unit = 6
345 
346  call initesputr()
347 
348  call esputr2001_availablecombinationsn(combinations_2001_n, err)
349  if(err /= 0) call cgi_error("esputr_cgi: error retreiving data")
350  call esputr2001_availablecombinationsth(combinations_2001_th, err)
351  if(err /= 0) call cgi_error("esputr_cgi: error retreiving data")
352  call esputr1993_availablecombinations(combinations_1993, err)
353  if(err /= 0) call cgi_error("esputr_cgi: error retreiving data")
354 
355  write(unit, '(A)') json_header
356  write(unit, *) '{'
357  call writejsonarray("1993", combinations_1993)
358  write(unit, *) ','
359  call writejsonarray("2007N", combinations_2001_n)
360  write(unit, *) ','
361  call writejsonarray("2007TH", combinations_2001_th)
362  write(unit, *) '}'
363 
364  call cgi_end()
365  end subroutine sendavailablecombinations
366 
372  subroutine writejsonarray(name, arr)
373  character(*), intent(in) :: name
374  character(*), dimension(:,:), intent(in) :: arr
375  integer :: i, iMax, err, indLo, indUp, pid, tid, ptid
376  real(kind=ESPUTR_DP) :: eMin, eMax, e_Th
377  write(unit, '(3A)') '"', name, '": ['
378  imax = size(arr, 1)
379  do i = 1, imax
380  write(unit, '(A)', advance='no') '{'
381  write(unit, '(3A)', advance='no') '"projectile":"', trim(arr(i, 1)), '",'
382  write(unit, '(3A)', advance='no') '"target":"', trim(arr(i, 2)), '"'
383  if(index(name, "TH") > 0) then
384  call esputr2001_getprojectiletargetidsth(arr(i, 1), arr(i, 2), indlo, indup, err)
385  if(err /= 0) call cgi_error("esputr_cgi: error while retrieving ids")
386  call esputr2001_getavailableenergyrange(indlo, indup, emin, emax, err)
387  if(err /= 0) call cgi_error("esputr_cgi: error while retrieving energy information")
388  write(unit, '(A,E12.5,A,E12.5,A)', advance='no') ',"energyRange":[',emin,',',emax,']'
389  else
390  if(index(name, "1993") > 0) then
391  pid = esputr1993_getprojectileid(arr(i, 1), err)
392  if(err /= 0) call cgi_error("esputr_cgi: error while retrieving projectile id")
393  tid = esputr1993_gettargetid(arr(i, 2), err)
394  if(err /= 0) call cgi_error("esputr_cgi: error while retrieving target id")
395  e_th = esputr1993_eth(pid, tid, err)
396  if(err /= 0) call cgi_error("esputr_cgi: error while retrieving E_th")
397  else
398  ptid = esputr2001_getprojectiletargetidn(arr(i, 1), arr(i, 2), err)
399  if(err /= 0) call cgi_error("esputr_cgi: error while retrieving projectile-target id")
400  e_th = esputr2001_eth(ptid, err)
401  if(err /= 0) call cgi_error("esputr_cgi: error while retrieving E_th")
402  end if
403  write(unit, '(A,E12.5)', advance='no') ',"E_th":',e_th
404  end if
405  write(unit, '(A)', advance='no') '}'
406  if(i < imax) write(unit, '(A)') ','
407  end do
408  write(unit, *) "]"
409  end subroutine writejsonarray
410 
412  subroutine initesputr()
413  call esputr1993_init("../data/SPUTER", err)
414  if(err /= 0) call cgi_error("esputr_cgi: initialization error")
415  call esputr2001_init("../data/ECKSTEIN2007N", "../data/ECKSTEIN2007TH", err)
416  if(err /= 0) call cgi_error("esputr_cgi: initialization error")
417  end subroutine initesputr
418 
420  character(DICT_VALUE_LENGTH*2) function getlastmessage()
422  backspace(tmp_unit, iostat=err)
423  backspace(tmp_unit, iostat=err)
424  read(tmp_unit, '(A)') getlastmessage(1:dict_value_length)
425  read(tmp_unit, '(A)') getlastmessage(len_trim(getlastmessage)+1:dict_value_length*2)
426  end function
427 
428 end program esputr_cgi
subroutine cgi_end
subroutine initesputr()
Initialize all types of models.
Definition: esputr_cgi.f90:413
integer function, public esputr1993_getprojectileid(proj, err)
Return ID of the projectile for the 1993-model.
Definition: esputr1993.f90:352
subroutine cgi_error(msg, template)
subroutine, public esputr1993_deallocate(err)
Deallocate dynamic arrays used by the module ESPUTR1993.
Definition: esputr1993.f90:329
subroutine toreal8(strArr, realArr)
Convert array of strings into array of reals.
Definition: esputr_cgi.f90:324
integer, parameter output_no_header
subroutine, public esputr2001_getavailableenergyrange(projTargStartId, projTargEndId, Emin, Emax, err)
Get the minimum and maximum energies for which the angular dependency factor is defined.
Definition: esputr2001.f90:776
subroutine, public esputr1993_init(constantsFile, err)
Initialization of the 1993-model.
Definition: esputr1993.f90:82
subroutine split(csvString, array)
Split one line into individual strings for each comma separated substring.
Definition: esputr_cgi.f90:314
character(dict_value_length *2) function getlastmessage()
Read two last lines printed in the output unit.
Definition: esputr_cgi.f90:421
subroutine, public esputr2001_deallocate(err)
Deallocate dynamic arrays used by this module.
Definition: esputr2001.f90:426
subroutine writejsonarray(name, arr)
Write the list of the projectile-target combinations in JSON (JavaScript Object Notation) format...
Definition: esputr_cgi.f90:373
integer function, public esputr2001_getprojectiletargetidn(proj, targ, err)
Return ID for a projectile-target combination for 2001-model for normal incidence.
Definition: esputr2001.f90:450
subroutine, public esputr2001_availablecombinationsth(combinations, err)
Return the list of available projectile-target combinations for angular dependence.
Definition: esputr2001.f90:740
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...
program esputr_cgi
Definition: esputr_cgi.f90:75
subroutine esputr_yn(E0s, n_E0s, proj, targ, version, YNs, err)
Calculate sputtering yield for normal incidence for selected incident energies and the specified proj...
integer function, public esputr1993_gettargetid(targ, err)
Return ID of the target for the 1993-model.
Definition: esputr1993.f90:382
real(esputr_dp) function, public esputr2001_eth(projectileTarget_id, err)
Return the threshold energy E_th for the given projectile-target combination.
Definition: esputr2001.f90:702
subroutine, public esputr1993_availablecombinations(combinations, err)
Return the list of available projectile-target combinations.
Definition: esputr1993.f90:532
subroutine cgi_begin(html, dict, luout)
subroutine sendavailablecombinations()
Send the list of available target-projectile combinations to the server.
Definition: esputr_cgi.f90:336
subroutine, public esputr2001_getprojectiletargetidsth(proj, targ, thLower, thUpper, err)
Return two IDs for a projectile-target combination for 2001-model for angular dependence.
Definition: esputr2001.f90:482
subroutine, public esputr2001_init(fileNName, fileThName, err)
Initialization of the 2001-model.
Definition: esputr2001.f90:128
subroutine, public esputr2001_availablecombinationsn(combinations, err)
Return the list of available projectile-target combinations for normal incidence. ...
Definition: esputr2001.f90:718
integer, save, public esputr_unit
Index of the unit for standard output, default value 6.
Definition: esputr.f90:73
real(esputr_dp) function, public esputr1993_eth(proj_id, targ_id, err)
Return threshold energy for the given projectile-target combination.
Definition: esputr1993.f90:505
integer, parameter dict_value_length