ETS  \$Id: Doxyfile 2162 2020-02-26 14:16:09Z g2dpc $
 All Classes Files Functions Variables Pages
generic_sources.f90
Go to the documentation of this file.
1 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
7 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
9 
10  USE euitm_schemas
11  USE itm_types
13  USE deallocate_structures
14  USE euitm_xml_parser
15  USE xml_file_reader
16 
17  implicit none
18 
19  public ecrh, icrh, nbi
20 
21  private
22 
23  integer, parameter :: buflen=256
24  REAL(R8), ALLOCATABLE :: rho(:), vprime(:), rhonrm(:)
25  TYPE (type_param), save :: code_parameters
26  integer :: nrho
27 !irena for impurity
28  INTEGER :: nnucl !number of nuclei species
29  INTEGER :: nion !number of ion species
30  INTEGER :: nimp, iimp !number of impurity species
31  INTEGER, ALLOCATABLE :: nzimp(:) !number of ionization states for each impurity
32  INTEGER :: nneut, ineut !number of neutrals species
33  INTEGER, ALLOCATABLE :: ncomp(:) !number of components for each neutral
34  INTEGER, ALLOCATABLE :: ntype(:) !number of types for each neutral
35 
36 contains
37 
38  subroutine ecrh(equilibrium, coreprof, coresource)
39 
40  TYPE (type_equilibrium), POINTER :: equilibrium(:)
41  TYPE (type_coreprof), POINTER :: coreprof(:)
42  TYPE (type_coresource), POINTER :: coresource(:)
43  character(len=BUFLEN), save :: j_src_f, qe_exp_f, qe_imp_f
44  character(len=BUFLEN), allocatable, save :: qi_exp_f(:), qi_imp_f(:), si_exp_f(:), si_imp_f(:)
45 !irena for impurity
46 ! character(len=BUFLEN), allocatable, save :: QZ_EXP_F(:,:), QZ_IMP_F(:,:), SZ_EXP_F(:,:), SZ_IMP_F(:,:)
47 
48  character(len=BUFLEN), allocatable, save :: ui_exp_f(:), ui_imp_f(:)
49  real(R8), save :: j_src_f_int, qe_exp_f_int, qe_imp_f_int
50  real(R8), allocatable, save :: qi_exp_f_int(:), qi_imp_f_int(:), si_exp_f_int(:), si_imp_f_int(:)
51 !irena for impurity
52 ! real(R8), allocatable, save :: QZ_EXP_F_INT(:), QZ_IMP_F_INT(:), SZ_EXP_F_INT(:), SZ_IMP_F_INT(:)
53 
54  real(R8), allocatable, save :: ui_exp_f_int(:), ui_imp_f_int(:)
55  integer :: return_status, irho, iion
56  logical, save :: first=.true.
57  REAL(R8) :: integral, integral_err
58 
59  write(*,*) 'ECRH module'
60 
61  if(associated(equilibrium(1)%profiles_1d%rho_vol)) then
62  nrho = size(equilibrium(1)%profiles_1d%rho_vol)
63  allocate(rhonrm(nrho))
64  rhonrm = equilibrium(1)%profiles_1d%rho_vol
65  write(*,*) 'RHONRM based on input equilibrium RHO_VOL'
66  elseif(associated(equilibrium(1)%profiles_1d%rho_tor)) then
67  nrho = size(equilibrium(1)%profiles_1d%rho_tor)
68  rhonrm = equilibrium(1)%profiles_1d%rho_tor / equilibrium(1)%profiles_1d%rho_tor(nrho)
69  write(*,*) 'RHONRM based on input equilibrium RHO_TOR'
70  elseif(associated(equilibrium(1)%profiles_1d%psi)) then
71  nrho = size(equilibrium(1)%profiles_1d%psi)
72  rhonrm = equilibrium(1)%profiles_1d%psi / equilibrium(1)%profiles_1d%psi(nrho)
73  write(*,*) 'RHONRM based on input equilibrium PSI'
74  else
75  write(*,*) 'Could not find a suitable radial coordinate'
76  stop
77  endif
78 
79 
80 
81 
82  if(first) then
83  CALL fill_param(code_parameters, 'XML/ecrh.xml', '', 'XML/sources.xsd')
84  call assign_code_parameters(code_parameters, return_status, &
85  j_src_f, qe_exp_f, qe_imp_f, qi_exp_f, qi_imp_f, &
86  si_exp_f, si_imp_f, ui_exp_f, ui_imp_f, &
87  j_src_f_int, qe_exp_f_int, qe_imp_f_int, qi_exp_f_int, qi_imp_f_int, &
88  si_exp_f_int, si_imp_f_int, ui_exp_f_int, ui_imp_f_int)
89  call deallocate_cpo(code_parameters)
90  first=.false.
91  endif
92 
93  CALL get_comp_dimensions(coreprof(1)%COMPOSITIONS, nnucl, nion, nimp, nzimp, nneut, ntype, ncomp)
94  CALL allocate_coresource_cpo(1, nrho, nnucl, nion, nimp, nzimp, nneut, ntype, ncomp, coresource)
95 
96 
97  allocate(rho(nrho), vprime(nrho))
98  rho = equilibrium(1)%profiles_1d%rho_tor
99  vprime = equilibrium(1)%profiles_1d%vprime
100 
101  coresource(1)%time = coreprof(1)%time
102  coresource(1)%VALUES(1)%rho_tor = equilibrium(1)%profiles_1d%rho_tor
103  coresource(1)%VALUES(1)%j(:) = profile(j_src_f, integral, integral_err)
104  coresource(1)%VALUES(1)%qe%exp(:) = profile(qe_exp_f, integral, integral_err)
105  if(qe_exp_f_int.ne.0.0_r8) then
106  coresource(1)%VALUES(1)%qe%exp(:) = coresource(1)%VALUES(1)%qe%exp(:) * qe_exp_f_int / integral
107  write(*,*) 'Integral Qe rescaled to ', qe_exp_f_int
108  ELSE
109  write(*,*) 'Integral Qe = ', integral, ' +- ', integral_err
110  ENDIF
111  coresource(1)%VALUES(1)%qe%imp(:) = profile(qe_imp_f, integral, integral_err)
112  do iion = 1, nion
113  coresource(1)%VALUES(1)%si%exp(:,iion) = profile(si_exp_f(iion), integral, integral_err)
114  if(si_exp_f_int(iion).ne.0.0_r8) then
115  coresource(1)%VALUES(1)%si%exp(:,iion) = coresource(1)%VALUES(1)%si%exp(:,iion) * si_exp_f_int(iion) / integral
116  write(*,*) 'Integral Si (',iion,') rescaled to ', si_exp_f_int(iion)
117  ELSE
118  write(*,*) 'Integral Si (',iion,') = ', integral, ' +- ', integral_err
119  ENDIF
120  coresource(1)%VALUES(1)%si%imp(:,iion) = profile(si_imp_f(iion), integral, integral_err)
121  coresource(1)%VALUES(1)%qi%exp(:,iion) = profile(qi_exp_f(iion), integral, integral_err)
122  if(qi_exp_f_int(iion).ne.0.0_r8) then
123  coresource(1)%VALUES(1)%qi%exp(:,iion) = coresource(1)%VALUES(1)%qi%exp(:,iion) * qi_exp_f_int(iion) / integral
124  write(*,*) 'Integral Qi (',iion,') rescaled to ', qi_exp_f_int(iion)
125  ELSE
126  write(*,*) 'Integral Qi (',iion,') = ', integral, ' +- ', integral_err
127  ENDIF
128  coresource(1)%VALUES(1)%qi%imp(:,iion) = profile(qi_imp_f(iion), integral, integral_err)
129  coresource(1)%VALUES(1)%ui%exp(:,iion) = profile(ui_exp_f(iion), integral, integral_err)
130  if(ui_exp_f_int(iion).ne.0.0_r8) then
131  coresource(1)%VALUES(1)%ui%exp(:,iion) = coresource(1)%VALUES(1)%ui%exp(:,iion) * ui_exp_f_int(iion) / integral
132  write(*,*) 'Integral Ui (',iion,') rescaled to ', ui_exp_f_int(iion)
133  ELSE
134  write(*,*) 'Integral Ui (',iion,') = ', integral, ' +- ', integral_err
135  ENDIF
136  coresource(1)%VALUES(1)%UI%IMP(:,iion) = profile(ui_imp_f(iion), integral, integral_err)
137  enddo
138 
139  deallocate(rhonrm, rho, vprime)
140 
141  end subroutine ecrh
142 
143  subroutine icrh(equilibrium, coreprof, coresource)
144 
145  TYPE (type_equilibrium), POINTER :: equilibrium(:)
146  TYPE (type_coreprof), POINTER :: coreprof(:)
147  TYPE (type_coresource), POINTER :: coresource(:)
148  character(len=BUFLEN), save :: j_src_f, qe_exp_f, qe_imp_f
149  character(len=BUFLEN), allocatable, save :: qi_exp_f(:), qi_imp_f(:), si_exp_f(:), si_imp_f(:)
150  character(len=BUFLEN), allocatable, save :: ui_exp_f(:), ui_imp_f(:)
151  real(R8), save :: j_src_f_int, qe_exp_f_int, qe_imp_f_int
152  real(R8), allocatable, save :: qi_exp_f_int(:), qi_imp_f_int(:), si_exp_f_int(:), si_imp_f_int(:)
153  real(R8), allocatable, save :: ui_exp_f_int(:), ui_imp_f_int(:)
154  integer :: return_status, irho, iion
155  logical, save :: first=.true.
156  REAL(R8) :: integral, integral_err
157 
158  write(*,*) 'ICRH module'
159 
160  if(associated(equilibrium(1)%profiles_1d%rho_vol)) then
161  nrho = size(equilibrium(1)%profiles_1d%rho_vol)
162  allocate(rhonrm(nrho))
163  rhonrm = equilibrium(1)%profiles_1d%rho_vol
164  write(*,*) 'RHONRM based on input equilibrium RHO_VOL'
165  elseif(associated(equilibrium(1)%profiles_1d%rho_tor)) then
166  nrho = size(equilibrium(1)%profiles_1d%rho_tor)
167  rhonrm = equilibrium(1)%profiles_1d%rho_tor / equilibrium(1)%profiles_1d%rho_tor(nrho)
168  write(*,*) 'RHONRM based on input equilibrium RHO_TOR'
169  elseif(associated(equilibrium(1)%profiles_1d%psi)) then
170  nrho = size(equilibrium(1)%profiles_1d%psi)
171  rhonrm = equilibrium(1)%profiles_1d%psi / equilibrium(1)%profiles_1d%psi(nrho)
172  write(*,*) 'RHONRM based on input equilibrium PSI'
173  else
174  write(*,*) 'Could not find a suitable radial coordinate'
175  stop
176  endif
177 
178  if(first) then
179  CALL fill_param(code_parameters, 'XML/icrh.xml', '', 'XML/sources.xsd')
180  call assign_code_parameters(code_parameters, return_status, &
181  j_src_f, qe_exp_f, qe_imp_f, qi_exp_f, qi_imp_f, &
182  si_exp_f, si_imp_f, ui_exp_f, ui_imp_f, &
183  j_src_f_int, qe_exp_f_int, qe_imp_f_int, qi_exp_f_int, qi_imp_f_int, &
184  si_exp_f_int, si_imp_f_int, ui_exp_f_int, ui_imp_f_int)
185  call deallocate_cpo(code_parameters)
186  first=.false.
187  endif
188 
189  CALL get_comp_dimensions(coreprof(1)%COMPOSITIONS, nnucl, nion, nimp, nzimp, nneut, ntype, ncomp)
190  CALL allocate_coresource_cpo(1, nrho, nnucl, nion, nimp, nzimp, nneut, ntype, ncomp, coresource)
191  allocate(rho(nrho), vprime(nrho))
192  rho = equilibrium(1)%profiles_1d%rho_tor
193  vprime = equilibrium(1)%profiles_1d%vprime
194 
195  coresource(1)%time = coreprof(1)%time
196  coresource(1)%VALUES(1)%rho_tor = equilibrium(1)%profiles_1d%rho_tor
197  coresource(1)%VALUES(1)%j(:) = profile(j_src_f, integral, integral_err)
198  coresource(1)%VALUES(1)%qe%exp(:) = profile(qe_exp_f, integral, integral_err)
199  if(qe_exp_f_int.ne.0.0_r8) then
200  coresource(1)%VALUES(1)%qe%exp(:) = coresource(1)%VALUES(1)%qe%exp(:) * qe_exp_f_int / integral
201  write(*,*) 'Integral Qe rescaled to ', qe_exp_f_int
202  ELSE
203  write(*,*) 'Integral Qe = ', integral, ' +- ', integral_err
204  ENDIF
205  coresource(1)%VALUES(1)%qe%imp(:) = profile(qe_imp_f, integral, integral_err)
206  do iion = 1, nion
207  coresource(1)%VALUES(1)%si%exp(:,iion) = profile(si_exp_f(iion), integral, integral_err)
208  if(si_exp_f_int(iion).ne.0.0_r8) then
209  coresource(1)%VALUES(1)%si%exp(:,iion) = coresource(1)%VALUES(1)%si%exp(:,iion) * si_exp_f_int(iion) / integral
210  write(*,*) 'Integral Si (',iion,') rescaled to ', si_exp_f_int(iion)
211  ELSE
212  write(*,*) 'Integral Si (',iion,') = ', integral, ' +- ', integral_err
213  ENDIF
214  coresource(1)%VALUES(1)%si%imp(:,iion) = profile(si_imp_f(iion), integral, integral_err)
215  coresource(1)%VALUES(1)%qi%exp(:,iion) = profile(qi_exp_f(iion), integral, integral_err)
216  if(qi_exp_f_int(iion).ne.0.0_r8) then
217  coresource(1)%VALUES(1)%qi%exp(:,iion) = coresource(1)%VALUES(1)%qi%exp(:,iion) * qi_exp_f_int(iion) / integral
218  write(*,*) 'Integral Qi (',iion,') rescaled to ', qi_exp_f_int(iion)
219  ELSE
220  write(*,*) 'Integral Qi (',iion,') = ', integral, ' +- ', integral_err
221  ENDIF
222  coresource(1)%VALUES(1)%qi%imp(:,iion) = profile(qi_imp_f(iion), integral, integral_err)
223  coresource(1)%VALUES(1)%ui%exp(:,iion) = profile(ui_exp_f(iion), integral, integral_err)
224  if(ui_exp_f_int(iion).ne.0.0_r8) then
225  coresource(1)%VALUES(1)%ui%exp(:,iion) = coresource(1)%VALUES(1)%ui%exp(:,iion) * ui_exp_f_int(iion) / integral
226  write(*,*) 'Integral Ui (',iion,') rescaled to ', ui_exp_f_int(iion)
227  ELSE
228  write(*,*) 'Integral Ui (',iion,') = ', integral, ' +- ', integral_err
229  ENDIF
230  coresource(1)%VALUES(1)%ui%imp(:,iion) = profile(ui_imp_f(iion), integral, integral_err)
231  enddo
232 
233  deallocate(rhonrm, rho, vprime)
234 
235  end subroutine icrh
236 
237  subroutine nbi(equilibrium, coreprof, coresource)
238 
239  TYPE (type_equilibrium), POINTER :: equilibrium(:)
240  TYPE (type_coreprof), POINTER :: coreprof(:)
241  TYPE (type_coresource), POINTER :: coresource(:)
242  character(len=BUFLEN), save :: j_src_f, qe_exp_f, qe_imp_f
243  character(len=BUFLEN), allocatable, save :: qi_exp_f(:), qi_imp_f(:), si_exp_f(:), si_imp_f(:)
244  character(len=BUFLEN), allocatable, save :: ui_exp_f(:), ui_imp_f(:)
245  real(R8), save :: j_src_f_int, qe_exp_f_int, qe_imp_f_int
246  real(R8), allocatable, save :: qi_exp_f_int(:), qi_imp_f_int(:), si_exp_f_int(:), si_imp_f_int(:)
247  real(R8), allocatable, save :: ui_exp_f_int(:), ui_imp_f_int(:)
248  integer :: return_status, irho, iion
249  logical, save :: first=.true.
250  REAL(R8) :: integral, integral_err
251 
252  write(*,*) 'NBI module'
253 
254  if(associated(equilibrium(1)%profiles_1d%rho_vol)) then
255  nrho = size(equilibrium(1)%profiles_1d%rho_vol)
256  allocate(rhonrm(nrho))
257  rhonrm = equilibrium(1)%profiles_1d%rho_vol
258  write(*,*) 'RHONRM based on input equilibrium RHO_VOL'
259  elseif(associated(equilibrium(1)%profiles_1d%rho_tor)) then
260  nrho = size(equilibrium(1)%profiles_1d%rho_tor)
261  rhonrm = equilibrium(1)%profiles_1d%rho_tor / equilibrium(1)%profiles_1d%rho_tor(nrho)
262  write(*,*) 'RHONRM based on input equilibrium RHO_TOR'
263  elseif(associated(equilibrium(1)%profiles_1d%psi)) then
264  nrho = size(equilibrium(1)%profiles_1d%psi)
265  rhonrm = equilibrium(1)%profiles_1d%psi / equilibrium(1)%profiles_1d%psi(nrho)
266  write(*,*) 'RHONRM based on input equilibrium PSI'
267  else
268  write(*,*) 'Could not find a suitable radial coordinate'
269  stop
270  endif
271 
272  if(first) then
273  CALL fill_param(code_parameters, 'XML/nbi.xml', '', 'XML/sources.xsd')
274  call assign_code_parameters(code_parameters, return_status, &
275  j_src_f, qe_exp_f, qe_imp_f, qi_exp_f, qi_imp_f, &
276  si_exp_f, si_imp_f, ui_exp_f, ui_imp_f, &
277  j_src_f_int, qe_exp_f_int, qe_imp_f_int, qi_exp_f_int, qi_imp_f_int, &
278  si_exp_f_int, si_imp_f_int, ui_exp_f_int, ui_imp_f_int)
279  call deallocate_cpo(code_parameters)
280  first=.false.
281  endif
282 
283  CALL get_comp_dimensions(coreprof(1)%COMPOSITIONS, nnucl, nion, nimp, nzimp, nneut, ntype, ncomp)
284  CALL allocate_coresource_cpo(1, nrho, nnucl, nion, nimp, nzimp, nneut, ntype, ncomp, coresource)
285 
286  allocate(rho(nrho), vprime(nrho))
287  rho = equilibrium(1)%profiles_1d%rho_tor
288  vprime = equilibrium(1)%profiles_1d%vprime
289 
290  coresource(1)%time = coreprof(1)%time
291  coresource(1)%VALUES(1)%rho_tor = equilibrium(1)%profiles_1d%rho_tor
292  coresource(1)%VALUES(1)%j(:) = profile(j_src_f, integral, integral_err)
293  coresource(1)%VALUES(1)%qe%exp(:) = profile(qe_exp_f, integral, integral_err)
294  if(qe_exp_f_int.ne.0.0_r8) then
295  coresource(1)%VALUES(1)%qe%exp(:) = coresource(1)%VALUES(1)%qe%exp(:) * qe_exp_f_int / integral
296  write(*,*) 'Integral Qe rescaled to ', qe_exp_f_int
297  ELSE
298  write(*,*) 'Integral Qe = ', integral, ' +- ', integral_err
299  ENDIF
300  coresource(1)%VALUES(1)%qe%imp(:) = profile(qe_imp_f, integral, integral_err)
301  do iion = 1, nion
302  coresource(1)%VALUES(1)%si%exp(:,iion) = profile(si_exp_f(iion), integral, integral_err)
303  if(si_exp_f_int(iion).ne.0.0_r8) then
304  coresource(1)%VALUES(1)%si%exp(:,iion) = coresource(1)%VALUES(1)%si%exp(:,iion) * si_exp_f_int(iion) / integral
305  write(*,*) 'Integral Si (',iion,') rescaled to ', si_exp_f_int(iion)
306  ELSE
307  write(*,*) 'Integral Si (',iion,') = ', integral, ' +- ', integral_err
308  ENDIF
309  coresource(1)%VALUES(1)%si%imp(:,iion) = profile(si_imp_f(iion), integral, integral_err)
310  coresource(1)%VALUES(1)%qi%exp(:,iion) = profile(qi_exp_f(iion), integral, integral_err)
311  if(qi_exp_f_int(iion).ne.0.0_r8) then
312  coresource(1)%VALUES(1)%qi%exp(:,iion) = coresource(1)%VALUES(1)%qi%exp(:,iion) * qi_exp_f_int(iion) / integral
313  write(*,*) 'Integral Qi (',iion,') rescaled to ', qi_exp_f_int(iion)
314  ELSE
315  write(*,*) 'Integral Qi (',iion,') = ', integral, ' +- ', integral_err
316  ENDIF
317  coresource(1)%VALUES(1)%qi%imp(:,iion) = profile(qi_imp_f(iion), integral, integral_err)
318  coresource(1)%VALUES(1)%ui%exp(:,iion) = profile(ui_exp_f(iion), integral, integral_err)
319  if(ui_exp_f_int(iion).ne.0.0_r8) then
320  coresource(1)%VALUES(1)%ui%exp(:,iion) = coresource(1)%VALUES(1)%ui%exp(:,iion) * ui_exp_f_int(iion) / integral
321  write(*,*) 'Integral Ui (',iion,') rescaled to ', ui_exp_f_int(iion)
322  ELSE
323  write(*,*) 'Integral Ui (',iion,') = ', integral, ' +- ', integral_err
324  ENDIF
325  coresource(1)%VALUES(1)%ui%imp(:,iion) = profile(ui_imp_f(iion), integral, integral_err)
326  enddo
327 
328  deallocate(rhonrm, rho, vprime)
329 
330  end subroutine nbi
331 
332  function profile(function_string, integral, integral_err)
333 
334  use fortranparser, only : equationparser
335 
336  real(R8) :: profile(1:nrho), integral, integral_err
337  character (len=BUFLEN) :: function_string
338 
339  type(equationparser) :: function_descriptor
340  character(len=10) :: variables(1) = ['x']
341 
342  integer :: i
343 
344  function_descriptor = equationparser(trim(function_string), variables)
345 ! if(.not.c_associated(function_descriptor)) then
346 ! write(*,*) 'Invalid function ', trim(function_string)
347 ! stop
348 ! endif
349 
350  do i = 1, nrho
351  profile(i) = function_descriptor%evaluate([rhonrm(i)])
352  enddo
353 
354 ! call evaluator_destroy(function_descriptor)
355 
356  call cubint(nrho, rho, vprime*profile, 1, nrho, integral, integral_err)
357 
358  end function profile
359 
360  subroutine assign_code_parameters(codeparameters, return_status, &
361  j_src_f, qe_exp_f, qe_imp_f, qi_exp_f, qi_imp_f, &
362  si_exp_f, si_imp_f, ui_exp_f, ui_imp_f, &
363  j_src_f_int, qe_exp_f_int, qe_imp_f_int, qi_exp_f_int, qi_imp_f_int, &
364  si_exp_f_int, si_imp_f_int, ui_exp_f_int, ui_imp_f_int)
365 
366  use mod_f90_kind
367 
368  implicit none
369 
370  character(len=BUFLEN) :: j_src_f, qe_exp_f, qe_imp_f
371  character(len=BUFLEN), allocatable :: qi_exp_f(:), qi_imp_f(:), si_exp_f(:), si_imp_f(:)
372  character(len=BUFLEN), allocatable :: ui_exp_f(:), ui_imp_f(:)
373  real(R8) :: j_src_f_int, qe_exp_f_int, qe_imp_f_int
374  real(R8), allocatable :: qi_exp_f_int(:), qi_imp_f_int(:), si_exp_f_int(:), si_imp_f_int(:)
375  real(R8), allocatable :: ui_exp_f_int(:), ui_imp_f_int(:)
376 
377  type (type_param), intent(in) :: codeparameters
378  integer(ikind), intent(out) :: return_status
379 
380  type(tree) :: parameter_list
381  type(element), pointer :: temp_pointer
382  integer(ikind) :: i, nparm, n_values
383  character(len = 132) :: cname
384  character (len=256), allocatable :: tmp_string(:)
385  integer :: lng
386 
387  j_src_f='0.0'
388  qe_exp_f='0.0'
389  qe_imp_f='0.0'
390  allocate(qi_exp_f(nion)) ; lng=nion ; qi_exp_f='0.0'
391  allocate(qi_imp_f(nion)) ; lng=nion ; qi_imp_f='0.0'
392  allocate(si_exp_f(nion)) ; lng=nion ; si_exp_f='0.0'
393  allocate(si_imp_f(nion)) ; lng=nion ; si_imp_f='0.0'
394  allocate(ui_exp_f(nion)) ; lng=nion ; ui_exp_f='0.0'
395  allocate(ui_imp_f(nion)) ; lng=nion ; ui_imp_f='0.0'
396  j_src_f_int=0.0_r8
397  qe_exp_f_int=0.0_r8
398  qe_imp_f_int=0.0_r8
399  allocate(qi_exp_f_int(nion)) ; lng=nion ; qi_exp_f_int=0.0_r8
400  allocate(qi_imp_f_int(nion)) ; lng=nion ; qi_imp_f_int=0.0_r8
401  allocate(si_exp_f_int(nion)) ; lng=nion ; si_exp_f_int=0.0_r8
402  allocate(si_imp_f_int(nion)) ; lng=nion ; si_imp_f_int=0.0_r8
403  allocate(ui_exp_f_int(nion)) ; lng=nion ; ui_exp_f_int=0.0_r8
404  allocate(ui_imp_f_int(nion)) ; lng=nion ; ui_imp_f_int=0.0_r8
405 
406  return_status = 0 ! no error
407 
408  !-- parse xml-string codeparameters%parameters
409 
410  write(*,*) 'Calling euitm_xml_parse'
411  call euitm_xml_parse(code_parameters, nparm, parameter_list)
412  write(*,*) 'Called euitm_xml_parse'
413 
414  !-- assign variables
415 
416  temp_pointer => parameter_list%first
417 
418  outer: do
419  cname = char2str(temp_pointer%cname) ! necessary for AIX
420  select case (cname)
421  case ('parameters')
422  temp_pointer => temp_pointer%child
423  cycle
424 
425 !-- coresource parameters
426  case ('coresource')
427  temp_pointer => temp_pointer%child
428  cycle
429  case ('j')
430  if (allocated(temp_pointer%cvalue)) &
431  j_src_f = char2str(temp_pointer%cvalue)
432  case ('qe_exp')
433  if (allocated(temp_pointer%cvalue)) &
434  qe_exp_f = char2str(temp_pointer%cvalue)
435  case ('qe_imp')
436  if (allocated(temp_pointer%cvalue)) &
437  qe_imp_f = char2str(temp_pointer%cvalue)
438  case ('qi_exp')
439  if (allocated(temp_pointer%cvalue)) &
440  call scan_str2str(char2str(temp_pointer%cvalue), 256, qi_exp_f, lng)
441  case ('qi_imp')
442  if (allocated(temp_pointer%cvalue)) &
443  call scan_str2str(char2str(temp_pointer%cvalue), 256, qi_imp_f, lng)
444  case ('si_exp')
445  if (allocated(temp_pointer%cvalue)) &
446  call scan_str2str(char2str(temp_pointer%cvalue), 256, si_exp_f, lng)
447  case ('si_imp')
448  if (allocated(temp_pointer%cvalue)) &
449  call scan_str2str(char2str(temp_pointer%cvalue), 256, si_imp_f, lng)
450  case ('ui_exp')
451  if (allocated(temp_pointer%cvalue)) &
452  call scan_str2str(char2str(temp_pointer%cvalue), 256, ui_exp_f, lng)
453  case ('ui_imp')
454  if (allocated(temp_pointer%cvalue)) &
455  call scan_str2str(char2str(temp_pointer%cvalue), 256, ui_imp_f, lng)
456  case ('j_int')
457  if (allocated(temp_pointer%cvalue)) &
458  call char2num(temp_pointer%cvalue, j_src_f_int)
459  case ('qe_exp_int')
460  if (allocated(temp_pointer%cvalue)) &
461  call char2num(temp_pointer%cvalue, qe_exp_f_int)
462  case ('qe_imp_int')
463  if (allocated(temp_pointer%cvalue)) &
464  call char2num(temp_pointer%cvalue, qe_imp_f_int)
465  case ('qi_exp_int')
466  if (allocated(temp_pointer%cvalue)) &
467  call scan_str2real(char2str(temp_pointer%cvalue), qi_exp_f_int, lng)
468  case ('qi_imp_int')
469  if (allocated(temp_pointer%cvalue)) &
470  call scan_str2real(char2str(temp_pointer%cvalue), qi_imp_f_int, lng)
471  case ('si_exp_int')
472  if (allocated(temp_pointer%cvalue)) &
473  call scan_str2real(char2str(temp_pointer%cvalue), si_exp_f_int, lng)
474  case ('si_imp_int')
475  if (allocated(temp_pointer%cvalue)) &
476  call scan_str2real(char2str(temp_pointer%cvalue), si_imp_f_int, lng)
477  case ('ui_exp_int')
478  if (allocated(temp_pointer%cvalue)) &
479  call scan_str2real(char2str(temp_pointer%cvalue), ui_exp_f_int, lng)
480  case ('ui_imp_int')
481  if (allocated(temp_pointer%cvalue)) &
482  call scan_str2real(char2str(temp_pointer%cvalue), ui_imp_f_int, lng)
483 
484 !-- default
485  case default
486  write(*, *) 'ERROR: invalid parameter', cname
487  return_status = 1
488  exit
489  end select
490  do
491  if (associated(temp_pointer%sibling)) then
492  temp_pointer => temp_pointer%sibling
493  exit
494  end if
495  if (associated(temp_pointer%parent, parameter_list%first )) &
496  exit outer
497  if (associated(temp_pointer%parent)) then
498  temp_pointer => temp_pointer%parent
499  else
500  write(*, *) 'ERROR: broken list.'
501  return
502  end if
503  end do
504  end do outer
505 
506  !-- destroy tree
507  call destroy_xml_tree(parameter_list)
508 
509  return
510 
511  end subroutine assign_code_parameters
512 
513 end module generic_sources
subroutine assign_code_parameters(codeparameters, return_status)
Definition: emeq.f90:671
subroutine, public icrh(equilibrium, coreprof, coresource)
subroutine get_comp_dimensions(COMPOSITIONS, NNUCL, NION, NIMP, NZIMP, NNEUT, NTYPE, NCOMP)
subroutine, public nbi(equilibrium, coreprof, coresource)
This module contains routines for allocation/deallocation if CPOs used in ETS.
subroutine allocate_coresource_cpo(NSLICE, NRHO, NNUCL, NION, NIMP, NZIMP, NNEUT, NTYPE, NCOMP, CORESOURCE)
This routine allocates CORESOURCE CPO.
subroutine, public ecrh(equilibrium, coreprof, coresource)
real(r8) function, dimension(1:size(x)) profile(function_string, x)
subroutine cubint(ntab, xtab, ftab, ia_in, ib_in, result, error)
Definition: cubint.f90:1
implements generic sources (ECRH, ICRH, NBI)
subroutine integral(n, h, r, f, int)
Definition: solution2.f90:527