ETS-Core  version:0.0.4-46-ge2d8
Core actors for the ETS-6
 All Classes Files Functions Variables Pages
convert.f90
Go to the documentation of this file.
1 MODULE convert
2 
3  IMPLICIT NONE
4 
5 
6 
7 
8 CONTAINS
9 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
10 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
12  (equilibrium_old, equilibrium_iter, coreprofiles_old, coreprofiles_iter, &
13  coretransport_iter, coresource_iter, coresolver_iter, &
14 !
15  geometry, profiles, transport, sources, impurity, numerics, &
16  composition_iter, diag)
17 
18 !-------------------------------------------------------!
19 ! This routine converts IDSs into the ETS !
20 ! derived types. !
21 !-------------------------------------------------------!
22 ! Source: --- !
23 ! Developers: D.Kalupin !
24 ! Contacts: Denis.Kalupin@euro-fusion.org !
25 ! !
26 ! Comments: --- !
27 ! !
28 !-------------------------------------------------------!
29 
30 
31  USE imas_constants_module
32  USE ids_schemas
33  USE l3interps, only: l3interp, l3deriv
35 
36  USE ets_plasma
38  USE get_composition
40  USE stringtools, only: num2str
42 
43  IMPLICIT NONE
44 
45 
46 ! +++ IDS derived types:
47  TYPE (ids_core_profiles) :: coreprofiles_old
48  TYPE (ids_core_profiles) :: coreprofiles_iter
49  TYPE (ids_equilibrium) :: equilibrium_old
50  TYPE (ids_equilibrium) :: equilibrium_iter
51  TYPE (ids_core_transport) :: coretransport_iter
52  TYPE (ids_core_sources) :: coresource_iter
53  TYPE (ids_transport_solver_numerics) :: coresolver_iter
54 
55 
56 ! +++ ETS derived types:
57  TYPE (magnetic_geometry) :: geometry !contains all geometry quantities
58  TYPE (plasma_profiles) :: profiles !contains profiles of plasma parameters
59  TYPE (transport_coefficients) :: transport !contains profiles of trasport coefficients
60  TYPE (sources_and_sinks) :: sources !contains profiles of sources
61  TYPE (impurity_profiles) :: impurity !contains profiles of impurities
62  TYPE (transport_solver_numerics) :: numerics !contains all parameters required by numerical part
63  TYPE (diagnostic) :: diag !contains warnings and error messages
64  TYPE (plasma_composition) :: composition_iter
65  TYPE (plasma_composition) :: composition_old
66  TYPE (plasma_composition) :: composition_transp
67  TYPE (plasma_composition) :: composition_source
68 
69  CHARACTER (len=132) :: ref_combiner_transp = "combined"
70  CHARACTER (len=132) :: ref_combiner_src = "total"
71 
72 ! +++ Dimensions & indexes:
73  INTEGER :: neq, neq_old
74  INTEGER :: nprof, nprof_old
75  INTEGER :: ndiff, nconv, nflux
76  INTEGER :: nsrc
77  INTEGER :: nsol
78 
79 
80 ! +++ Local variables:
81  REAL (IDS_REAL), ALLOCATABLE :: fun1(:), fun2(:), fun3(:), fun4(:)
82  REAL (IDS_REAL), ALLOCATABLE :: dfun1(:), dfun2(:)
83 
84  INTEGER (IDS_INT) :: ifound
85  INTEGER (IDS_INT) :: ids_iion
86  INTEGER (IDS_INT) :: nion, nimp, max_nzimp
87  INTEGER (IDS_INT) :: iion, iimp, izimp
88  INTEGER (IDS_INT) :: i_index
89  INTEGER (IDS_INT), ALLOCATABLE :: ion_index(:), imp_index(:)
90  INTEGER (IDS_INT) :: ind, isol, imod, isrc, ibnd, ib0, i, ival
91  INTEGER (IDS_INT) :: max_solver_number
92  INTEGER (IDS_INT) :: ifail
93 
94 
95 
96 ! +++ Check if input IDSs follow the convention of homogeneous time:
97 ! IF(COREPROFILES_OLD%ids_properties%homogeneous_time.NE.1) CALL REPORT_ERROR (-1, "inhomogeneous time in COREPROFILES_OLD", DIAG)
98 ! IF(COREPROFILES_ITER%ids_properties%homogeneous_time.NE.1) CALL REPORT_ERROR (-1, "inhomogeneous time in COREPROFILES_ITER", DIAG)
99 ! IF(EQUILIBRIUM_OLD%ids_properties%homogeneous_time.NE.1) CALL REPORT_ERROR (-1, "inhomogeneous time in EQUILIBRIUM_OLD", DIAG)
100 ! IF(EQUILIBRIUM_ITER%ids_properties%homogeneous_time.NE.1) CALL REPORT_ERROR (-1, "inhomogeneous time in EQUILIBRIUM_ITER", DIAG)
101 ! IF(CORETRANSPORT_ITER%ids_properties%homogeneous_time.NE.1) CALL REPORT_ERROR (-1, "inhomogeneous time in CORETRANSPORT_ITER", DIAG)
102 ! IF(CORESOURCE_ITER%ids_properties%homogeneous_time.NE.1) CALL REPORT_ERROR (-1, "inhomogeneous time in CORESOURCE_ITER", DIAG)
103 ! IF(CORESOLVER_ITER%ids_properties%homogeneous_time.NE.1) CALL REPORT_ERROR (-1, "inhomogeneous time in CORESOLVER_ITER", DIAG)
104 ! IF(DIAG%ERROR_INDEX.LT.0) GOTO 666
105 
106 
107 
108 !!! assumption: all arrays of structures have only one element, or only the first element is read
109 ! IF(SIZE(EQUILIBRIUM_ITER%time_slice).GT.1) CALL REPORT_ERROR (-1, "> 1 array elements in EQUILIBRIUM_ITER", DIAG)
110 ! IF(SIZE(EQUILIBRIUM_OLD%time_slice).GT.1) CALL REPORT_ERROR (-1, "> 1 array elements in EQUILIBRIUM_OLD", DIAG)
111 ! IF(SIZE(COREPROFILES_ITER%profiles_1d).GT.1) CALL REPORT_ERROR (-1, "> 1 array elements in COREPROFILES_ITER", DIAG)
112 ! IF(SIZE(COREPROFILES_OLD%profiles_1d).GT.1) CALL REPORT_ERROR (-1, "> 1 array elements in COREPROFILES_OLD", DIAG)
113 ! IF(SIZE(CORETRANSPORT_ITER%MODEL(1)%profiles_1d).GT.1) CALL REPORT_ERROR (-1, "> 1 array elements in CORETRANSPORT_ITER", DIAG)
114 ! IF(SIZE(CORESOURCE_ITER%SOURCE(1)%profiles_1d).GT.1) CALL REPORT_ERROR (-1, "> 1 array elements in CORESOURCE_ITER", DIAG)
115 ! IF(SIZE(CORESOLVER_ITER%solver_1d).GT.1) CALL REPORT_ERROR (-1, "> 1 array elements in CORESOLVER_ITER", DIAG)
116 ! IF(DIAG%ERROR_INDEX.LT.0) GOTO 666
117 
118 
119 ! +++ Dimensions & Grid:
120  !Check grids
121  CALL check_coreprof_grids(coreprofiles_iter)
122  CALL check_coreprof_grids(coreprofiles_old)
123  CALL check_coretransp_grids(coretransport_iter)
124  CALL check_coresource_grids(coresource_iter)
125 
126  neq = SIZE(equilibrium_iter%time_slice(1)%profiles_1d%psi)
127  neq_old = SIZE(equilibrium_old%time_slice(1)%profiles_1d%psi)
128 
129  nprof = SIZE(coreprofiles_iter%profiles_1d(1)%grid%rho_tor_norm)
130  nprof_old = SIZE(coreprofiles_old%profiles_1d(1)%grid%rho_tor_norm)
131 
132  geometry%NRHO = nprof
133  profiles%NRHO = nprof
134  transport%NRHO = nprof
135  sources%NRHO = nprof
136  impurity%NRHO = nprof
137 
138 
139 
140 ! +++ Get composition from input IDSs:
141  CALL get_prof_composition(coreprofiles_iter, composition_iter)
142  CALL get_prof_composition(coreprofiles_old, composition_old)
143  nion = composition_iter%NION
144  nimp = composition_iter%NIMP
145  max_nzimp = composition_iter%MAX_NZIMP
146 
147 
148 ! +++ Allocation of transport solver derived types:
149  CALL allocate_magnetic_geometry(nprof, geometry, ifail)
150  CALL allocate_plasma_profiles(nprof, nion, profiles, ifail)
151  CALL allocate_transport_coefficients(nprof, nion, transport, ifail)
152  CALL allocate_sources_and_sinks(nprof, nion, sources, ifail)
153  CALL allocate_impurity_profiles(nprof, nimp, max_nzimp, impurity, ifail)
154 
155 
156 ! +++ Allocation of numerics, MAX: psi + ne + Te + NION * 3(ni+Ti+Vtor) + NIMP * MAX_NZIMP:
157  max_solver_number = 1 + 2 + nion*3 + nimp*max_nzimp
158 
159  ALLOCATE (numerics%SOLVER(max_solver_number))
160 
161 ! +++ Grids:
162  geometry%RHO_NORM = coreprofiles_iter%profiles_1d(1)%grid%rho_tor_norm
163  geometry%RHO = geometry%RHO_NORM*equilibrium_iter%time_slice(1)%profiles_1d%rho_tor(neq)
164 
165 
166 ! +++ Convert geometry:
167  geometry%R0 = equilibrium_iter%vacuum_toroidal_field%r0
168  geometry%B0 = equilibrium_iter%vacuum_toroidal_field%b0(1)
169  geometry%RGEO = equilibrium_iter%time_slice(1)%boundary%geometric_axis%r
170  geometry%BGEO = geometry%B0*geometry%R0/geometry%RGEO
171 
172  geometry%PHI_BND = equilibrium_iter%time_slice(1)%profiles_1d%phi(neq)
173  geometry%RHO_BND = equilibrium_iter%time_slice(1)%profiles_1d%rho_tor(neq)
174  IF (equilibrium_iter%time(1).GT.equilibrium_old%time(1)) THEN
175  geometry%PHI_BND_PRIME= (equilibrium_iter%time_slice(1)%profiles_1d%phi(neq)-equilibrium_old%time_slice(1)%profiles_1d%phi(neq_old))&
176  /(equilibrium_iter%time(1)-equilibrium_old%time(1))
177  geometry%RHO_BND_PRIME= (equilibrium_iter%time_slice(1)%profiles_1d%rho_tor(neq)-equilibrium_old%time_slice(1)%profiles_1d%rho_tor(neq_old))&
178  /(equilibrium_iter%time(1)-equilibrium_old%time(1))
179  END IF
180 
181  CALL l3interp_pointer( &
182  equilibrium_iter%time_slice(1)%profiles_1d%volume, &
183  equilibrium_iter%time_slice(1)%profiles_1d%rho_tor_norm, &
184  geometry%VOL, geometry%RHO_NORM)
185  CALL l3interp_pointer( &
186  equilibrium_iter%time_slice(1)%profiles_1d%area, &
187  equilibrium_iter%time_slice(1)%profiles_1d%rho_tor_norm, &
188  geometry%AREA, geometry%RHO_NORM)
189  CALL l3interp_pointer( &
190  equilibrium_iter%time_slice(1)%profiles_1d%surface, &
191  equilibrium_iter%time_slice(1)%profiles_1d%rho_tor_norm, &
192  geometry%SURFACE, geometry%RHO_NORM)
193  CALL l3deriv_pointer( &
194  equilibrium_iter%time_slice(1)%profiles_1d%volume, &
195  equilibrium_iter%time_slice(1)%profiles_1d%rho_tor_norm, &
196  geometry%VPR, geometry%RHO_NORM)
197  geometry%VPR = geometry%VPR / equilibrium_iter%time_slice(1)%profiles_1d%rho_tor(neq)
198  CALL l3deriv_pointer( &
199  equilibrium_old%time_slice(1)%profiles_1d%volume, &
200  equilibrium_old%time_slice(1)%profiles_1d%rho_tor_norm, &
201  geometry%VPRM,geometry%RHO_NORM)
202  geometry%VPRM = geometry%VPRM / equilibrium_old%time_slice(1)%profiles_1d%rho_tor(neq_old)
203  CALL l3interp_pointer( &
204  equilibrium_iter%time_slice(1)%profiles_1d%gm3, &
205  equilibrium_iter%time_slice(1)%profiles_1d%rho_tor_norm, &
206  geometry%G1, geometry%RHO_NORM)
207  CALL l3interp_pointer( &
208  equilibrium_iter%time_slice(1)%profiles_1d%gm8, &
209  equilibrium_iter%time_slice(1)%profiles_1d%rho_tor_norm, &
210  geometry%G2, geometry%RHO_NORM)
211  CALL l3interp_pointer( &
212  equilibrium_old%time_slice(1)%profiles_1d%gm8, &
213  equilibrium_old%time_slice(1)%profiles_1d%rho_tor_norm, &
214  geometry%G2M, geometry%RHO_NORM)
215  CALL l3interp_pointer( &
216  equilibrium_iter%time_slice(1)%profiles_1d%gm2, &
217  equilibrium_iter%time_slice(1)%profiles_1d%rho_tor_norm, &
218  geometry%G3, geometry%RHO_NORM)
219  CALL l3interp_pointer( &
220  equilibrium_iter%time_slice(1)%profiles_1d%f, &
221  equilibrium_iter%time_slice(1)%profiles_1d%rho_tor_norm, &
222  geometry%FDIA,geometry%RHO_NORM)
223 
224  IF (geometry%VPR(1).LE.0._ids_real) geometry%VPR(1) = 0._ids_real
225  IF (geometry%VPRM(1).LE.0._ids_real) geometry%VPRM(1) = 0._ids_real
226 
227 
228 
229 
230 
231 
232 ! +++ Convert profiles:
233  ALLOCATE (fun1(nprof), dfun1(nprof), fun2(nprof_old), dfun2(nprof))
234 
235 ! PSI&CURRENT
236  profiles%PSI = coreprofiles_iter%profiles_1d(1)%grid%psi
237  CALL l3deriv_pointer( coreprofiles_iter%profiles_1d(1)%grid%psi, &
238  coreprofiles_iter%profiles_1d(1)%grid%rho_tor_norm, &
239  profiles%DPSI, geometry%RHO_NORM)
240  profiles%DPSI = profiles%DPSI/geometry%RHO(nprof)
241  CALL l3interp_pointer( coreprofiles_old%profiles_1d(1)%grid%psi, &
242  coreprofiles_old%profiles_1d(1)%grid%rho_tor_norm, &
243  profiles%PSIM, geometry%RHO_NORM)
244  CALL l3deriv_pointer( coreprofiles_old%profiles_1d(1)%grid%psi, &
245  coreprofiles_old%profiles_1d(1)%grid%rho_tor_norm, &
246  profiles%DPSIM, geometry%RHO_NORM)
247  profiles%DPSIM = profiles%DPSIM/geometry%RHO(nprof)
248  CALL l3interp_pointer( equilibrium_iter%time_slice(1)%profiles_1d%phi, &
249  equilibrium_iter%time_slice(1)%profiles_1d%rho_tor_norm, &
250  profiles%PHI, geometry%RHO_NORM)
251 
252  IF (ASSOCIATED(coreprofiles_iter%profiles_1d(1)%q)) THEN
253  profiles%QSF = coreprofiles_iter%profiles_1d(1)%q
254  END IF
255  IF (ASSOCIATED(coreprofiles_iter%profiles_1d(1)%magnetic_shear)) THEN
256  profiles%SHEAR = coreprofiles_iter%profiles_1d(1)%magnetic_shear
257  END IF
258  IF (ASSOCIATED(coreprofiles_iter%profiles_1d(1)%j_total)) THEN
259  profiles%CURR_PAR = coreprofiles_iter%profiles_1d(1)%j_total
260  END IF
261  IF (ASSOCIATED(coreprofiles_iter%profiles_1d(1)%j_tor)) THEN
262  profiles%CURR_TOR = coreprofiles_iter%profiles_1d(1)%j_tor
263  END IF
264  IF (ASSOCIATED(coreprofiles_iter%profiles_1d(1)%j_bootstrap)) THEN
265  profiles%JBOOT = coreprofiles_iter%profiles_1d(1)%j_bootstrap
266  END IF
267  IF (ASSOCIATED(coreprofiles_iter%profiles_1d(1)%j_ohmic)) THEN
268  profiles%JOH = coreprofiles_iter%profiles_1d(1)%j_ohmic
269  END IF
270  IF (ASSOCIATED(coreprofiles_iter%profiles_1d(1)%j_non_inductive)) THEN
271  profiles%JNI = coreprofiles_iter%profiles_1d(1)%j_non_inductive
272  END IF
273 
274 !NE
275 ! - THERMAL (NEW)
276  IF (ASSOCIATED(coreprofiles_iter%profiles_1d(1)%electrons%density_thermal)) THEN
277  fun1 = coreprofiles_iter%profiles_1d(1)%electrons%density_thermal
278  profiles%NE = fun1
279  IF (.NOT. acceptable_electron_density(profiles%NE)) then
280  diag%ERROR_INDEX = -1
281  diag%ERROR_MESSAGE = 'ERROR in CONVERT_IDS_TO_INTERNAL_TYPES, input PROFILES%NE includes unacceptably small value '//num2str(minval(profiles%NE))
282  WRITE(*,*)diag%ERROR_MESSAGE
283  RETURN
284  END IF
285  CALL l3deriv(fun1, &
286  coreprofiles_iter%profiles_1d(1)%grid%rho_tor_norm, nprof, &
287  profiles%DNE, geometry%RHO_NORM, geometry%NRHO)
288  profiles%DNE = profiles%DNE/geometry%RHO(nprof)
289  END IF
290 
291 ! - FAST (ITER)
292  IF(ASSOCIATED(coreprofiles_iter%profiles_1d(1)%electrons%density_fast)) THEN
293  profiles%NE_FAST = coreprofiles_iter%profiles_1d(1)%electrons%density_fast
294  END IF
295 
296 ! - FAST (OLD)
297  IF (ASSOCIATED(coreprofiles_old%profiles_1d(1)%electrons%density_thermal)) THEN
298  fun2 = coreprofiles_old%profiles_1d(1)%electrons%density_thermal
299  CALL l3interp(fun2, &
300  coreprofiles_old%profiles_1d(1)%grid%rho_tor_norm, nprof_old, &
301  profiles%NEM, geometry%RHO_NORM, geometry%NRHO)
302  CALL l3deriv(fun2, &
303  coreprofiles_old%profiles_1d(1)%grid%rho_tor_norm, nprof_old, &
304  profiles%DNEM, geometry%RHO_NORM, geometry%NRHO)
305  profiles%DNEM = profiles%DNEM/geometry%RHO(nprof)
306  END IF
307 
308 ! - ZEFF
309  IF (ASSOCIATED(coreprofiles_iter%profiles_1d(1)%zeff)) THEN
310  profiles%ZEFF = coreprofiles_iter%profiles_1d(1)%zeff
311  END IF
312 
313 !TE
314  IF (ASSOCIATED(coreprofiles_iter%profiles_1d(1)%electrons%temperature)) THEN
315  profiles%TE = coreprofiles_iter%profiles_1d(1)%electrons%temperature
316  IF (.NOT. acceptable_temperature(profiles%TE)) then
317  diag%ERROR_INDEX = -1
318  diag%ERROR_MESSAGE = 'ERROR in CONVERT_IDS_TO_INTERNAL_TYPES, input PROFILES%TE includes the unacceptably small value '//num2str(minval(profiles%TE))
319  WRITE(*,*)diag%ERROR_MESSAGE
320  RETURN
321  END IF
322  CALL l3deriv(coreprofiles_iter%profiles_1d(1)%electrons%temperature, &
323  coreprofiles_iter%profiles_1d(1)%grid%rho_tor_norm, nprof, &
324  profiles%DTE, geometry%RHO_NORM, geometry%NRHO)
325  profiles%DTE = profiles%DTE/geometry%RHO(nprof)
326  END IF
327  IF (ASSOCIATED(coreprofiles_old%profiles_1d(1)%electrons%temperature)) THEN
328  CALL l3interp(coreprofiles_old%profiles_1d(1)%electrons%temperature, &
329  coreprofiles_old%profiles_1d(1)%grid%rho_tor_norm, nprof_old, &
330  profiles%TEM, geometry%RHO_NORM, geometry%NRHO)
331  CALL l3deriv(coreprofiles_old%profiles_1d(1)%electrons%temperature, &
332  coreprofiles_old%profiles_1d(1)%grid%rho_tor_norm, nprof_old, &
333  profiles%DTEM, geometry%RHO_NORM, geometry%NRHO)
334  profiles%DTEM = profiles%DTEM/geometry%RHO(nprof)
335  END IF
336 
337 !PE
338  IF (ASSOCIATED(coreprofiles_iter%profiles_1d(1)%electrons%pressure_thermal)) THEN
339  profiles%PE = coreprofiles_iter%profiles_1d(1)%electrons%pressure_thermal
340  END IF
341  IF (ASSOCIATED(coreprofiles_iter%profiles_1d(1)%electrons%pressure_fast_perpendicular)) THEN
342  profiles%PE_FAST = coreprofiles_iter%profiles_1d(1)%electrons%pressure_fast_perpendicular
343  END IF
344  IF (ASSOCIATED(coreprofiles_iter%profiles_1d(1)%electrons%pressure_fast_parallel)) THEN
345  profiles%PE_FAST_PARALLEL = coreprofiles_iter%profiles_1d(1)%electrons%pressure_fast_parallel
346  END IF
347 
348 !NI
349  DO iion = 1, composition_iter%NION
350 ! IF (COMPOSITION_ITER%IONS(IION)%IDS_ion.NE.COMPOSITION_OLD%IONS(IION)%IDS_ion .OR. &
351 ! COMPOSITION_ITER%IONS(IION)%ELINDEX.NE.COMPOSITION_OLD%IONS(IION)%ELINDEX .OR. &
352 ! ABS(1.0_IDS_REAL - COMPOSITION_ITER%IONS(IION)%ZION &
353 ! / COMPOSITION_OLD%IONS(IION)%ZION) .GE. 0.1_IDS_REAL) &
354 ! CALL REPORT_ERROR (-1, "OLD and ITER compositions do not agree", DIAG)
355 ! IF(DIAG%ERROR_INDEX.LT.0) GOTO 666
356 
357  profiles%ZION(iion) = composition_iter%IONS(iion)%ZION
358  profiles%ZION2(iion) = composition_iter%IONS(iion)%ZION**2.0_ids_real
359  profiles%MION(iion) = 0.0_ids_real
360  ind = composition_iter%IONS(iion)%ELINDEX
361  DO i = 1, SIZE(composition_iter%ELEMENTS(ind)%NUCINDEX)
362  profiles%MION(iion) = profiles%MION(iion) + composition_iter%NUCLEI(composition_iter%ELEMENTS(ind)%NUCINDEX(i))%amn &
363  * composition_iter%ELEMENTS(ind)%NUCMULT(i)
364  END DO
365 
366  ind = composition_iter%IONS(iion)%IDS_ion
367 
368 ! NI
369  IF (ASSOCIATED(coreprofiles_iter%profiles_1d(1)%ion(ind)%density_fast)) THEN
370  profiles%NI_FAST(:,iion) = coreprofiles_iter%profiles_1d(1)%ion(ind)%density_fast(:)
371  END IF
372 
373  IF (ASSOCIATED(coreprofiles_iter%profiles_1d(1)%ion(ind)%density_thermal)) THEN
374  fun1(:) = coreprofiles_iter%profiles_1d(1)%ion(ind)%density_thermal(:)
375  profiles%NI(:,iion) = fun1(:)
376  IF (.NOT. acceptable_ion_density(profiles%NI(:,iion))) then
377  diag%ERROR_INDEX = -1
378  diag%ERROR_MESSAGE = 'ERROR in CONVERT_IDS_TO_INTERNAL_TYPES, input PROFILES%NI(:,'//num2str(iion)//') includes the unacceptably small value '//num2str(minval(profiles%NI(:,iion)))
379  WRITE(*,*)diag%ERROR_MESSAGE
380  RETURN
381  END IF
382  CALL l3deriv(fun1, &
383  coreprofiles_iter%profiles_1d(1)%grid%rho_tor_norm, nprof, &
384  dfun1, geometry%RHO_NORM, geometry%NRHO)
385  profiles%DNI(:,iion) = dfun1(:)/geometry%RHO(nprof)
386  END IF
387 
388  IF (ASSOCIATED(coreprofiles_old%profiles_1d(1)%ion(ind)%density_thermal)) THEN
389  fun2(:) = coreprofiles_old%profiles_1d(1)%ion(ind)%density_thermal(:)
390  CALL l3interp(fun2, &
391  coreprofiles_old%profiles_1d(1)%grid%rho_tor_norm, nprof_old, &
392  dfun2, geometry%RHO_NORM, geometry%NRHO)
393  profiles%NIM(:,iion) = dfun2(:)
394  CALL l3deriv(fun2, &
395  coreprofiles_old%profiles_1d(1)%grid%rho_tor_norm, nprof_old, &
396  dfun2, geometry%RHO_NORM, geometry%NRHO)
397  profiles%DNIM(:,iion) = dfun2(:)/geometry%RHO(nprof)
398  END IF
399 
400 !TI
401  IF (ASSOCIATED(coreprofiles_iter%profiles_1d(1)%ion(ind)%temperature)) THEN
402  profiles%TI(:,iion) = coreprofiles_iter%profiles_1d(1)%ion(ind)%temperature(:)
403  IF (.NOT. acceptable_temperature(profiles%TI(:,iion))) then
404  diag%ERROR_INDEX = -1
405  diag%ERROR_MESSAGE = 'ERROR in CONVERT_IDS_TO_INTERNAL_TYPES, input PROFILES%TI(:,'//num2str(iion)//') includes unacceptably small value '//num2str(minval(profiles%TI(:,iion)))
406  WRITE(*,*)diag%ERROR_MESSAGE
407  RETURN
408  END IF
409  fun1(:) = coreprofiles_iter%profiles_1d(1)%ion(ind)%temperature(:)
410  CALL l3deriv(fun1, &
411  coreprofiles_iter%profiles_1d(1)%grid%rho_tor_norm, nprof, &
412  dfun1, geometry%RHO_NORM, geometry%NRHO)
413  profiles%DTI(:,iion) = dfun1(:)/geometry%RHO(nprof)
414  END IF
415 
416  IF (ASSOCIATED(coreprofiles_old%profiles_1d(1)%ion(ind)%temperature)) THEN
417  fun2(:) = coreprofiles_old%profiles_1d(1)%ion(ind)%temperature(:)
418  CALL l3interp(fun2, &
419  coreprofiles_old%profiles_1d(1)%grid%rho_tor_norm, nprof_old, &
420  dfun2, geometry%RHO_NORM, geometry%NRHO)
421  profiles%TIM(:,iion) = dfun2(:)
422  CALL l3deriv(fun2, &
423  coreprofiles_old%profiles_1d(1)%grid%rho_tor_norm, nprof_old, &
424  dfun2, geometry%RHO_NORM, geometry%NRHO)
425  profiles%DTIM(:,iion) = dfun2(:)/geometry%RHO(nprof)
426  END IF
427 
428 !PI
429  IF(ASSOCIATED(coreprofiles_iter%profiles_1d(1)%ion(ind)%pressure_thermal)) THEN
430  profiles%PI(:,iion) = coreprofiles_iter%profiles_1d(1)%ion(ind)%pressure_thermal(:)
431  END IF
432  IF(ASSOCIATED(coreprofiles_iter%profiles_1d(1)%ion(ind)%pressure_fast_perpendicular)) THEN
433  profiles%PI_FAST(:,iion) = coreprofiles_iter%profiles_1d(1)%ion(ind)%pressure_fast_perpendicular(:)
434  END IF
435  IF(ASSOCIATED(coreprofiles_iter%profiles_1d(1)%ion(ind)%pressure_fast_parallel)) THEN
436  profiles%PI_FAST_PARALLEL(:,iion) = coreprofiles_iter%profiles_1d(1)%ion(ind)%pressure_fast_parallel(:)
437  END IF
438 
439 !VTOR
440  IF(ASSOCIATED(coreprofiles_iter%profiles_1d(1)%ion(ind)%velocity%toroidal)) THEN
441  profiles%VTOR(:,iion) = coreprofiles_iter%profiles_1d(1)%ion(ind)%velocity%toroidal(:)
442  fun1(:) = coreprofiles_iter%profiles_1d(1)%ion(ind)%velocity%toroidal(:)
443  ELSE
444  profiles%VTOR(:,iion) = 0.0_ids_real
445  fun1(:) = 0.0_ids_real
446  END IF
447 
448  CALL l3deriv(fun1, &
449  coreprofiles_iter%profiles_1d(1)%grid%rho_tor_norm, nprof, &
450  dfun1, geometry%RHO_NORM, geometry%NRHO)
451  profiles%DVTOR(:,iion)= dfun1(:)/geometry%RHO(nprof)
452 
453  IF(ASSOCIATED(coreprofiles_old%profiles_1d(1)%ion(ind)%velocity%toroidal)) THEN
454  fun2(:) = coreprofiles_old%profiles_1d(1)%ion(ind)%velocity%toroidal(:)
455  ELSE
456  fun2(:) = 0.0_ids_real
457  END IF
458 
459  CALL l3interp(fun2, &
460  coreprofiles_old%profiles_1d(1)%grid%rho_tor_norm, nprof_old, &
461  dfun2, geometry%RHO_NORM, geometry%NRHO)
462  profiles%VTORM(:,iion)= dfun2(:)
463  CALL l3deriv(fun2, &
464  coreprofiles_old%profiles_1d(1)%grid%rho_tor_norm, nprof_old, &
465  dfun2, geometry%RHO_NORM, geometry%NRHO)
466  profiles%DVTORM(:,iion) = dfun2(:)/geometry%RHO(nprof)
467  END DO
468 
469 ! +++ Convert impurities:
470  loop_over_input_impurities: DO iimp = 1, composition_iter%NIMP
471  ind = composition_iter%IMPURITY(iimp)%IDS_ion
472  impurity%NZIMP(iimp) = SIZE(composition_iter%IMPURITY(iimp)%IDS_state)
473  impurity%MAX_NZIMP = maxval(impurity%NZIMP)
474  loop_over_impurity_ionization_state: DO izimp = 1, impurity%NZIMP(iimp)
475 
476  IF (ids_is_valid(coreprofiles_iter%profiles_1d(1)%ion(ind)%state(izimp)%z_min)) then
477  impurity%ZMIN(iimp,izimp) = coreprofiles_iter%profiles_1d(1)%ion(ind)%state(izimp)%z_min
478  ELSE
479  diag%ERROR_INDEX = -1
480  diag%ERROR_MESSAGE = 'Minimum impurity charge not available in COREPROFILES_ITER%profiles_1d(1)%ion(IND)%state(IZIMP)%z_min'
481  return
482  END IF
483 
484  IF (ids_is_valid(coreprofiles_iter%profiles_1d(1)%ion(ind)%state(izimp)%z_max)) then
485  impurity%ZMAX(iimp,izimp) = coreprofiles_iter%profiles_1d(1)%ion(ind)%state(izimp)%z_max
486  ELSE
487  diag%ERROR_INDEX = -1
488  diag%ERROR_MESSAGE = 'Maximum impurity charge not available in COREPROFILES_ITER%profiles_1d(1)%ion(IND)%state(IZIMP)%z_max'
489  return
490  END IF
491 
492  IF (ASSOCIATED(coreprofiles_iter%profiles_1d(1)%ion(ind)%state(izimp)%z_average_1d)) THEN
493  impurity%Z_AVE(:,iimp,izimp) = coreprofiles_iter%profiles_1d(1)%ion(ind)%state(izimp)%z_average_1d(:)
494  ELSEIF (ids_is_valid(coreprofiles_iter%profiles_1d(1)%ion(ind)%state(izimp)%z_min) .AND. &
495  (coreprofiles_iter%profiles_1d(1)%ion(ind)%state(izimp)%z_min .EQ. &
496  coreprofiles_iter%profiles_1d(1)%ion(ind)%state(izimp)%z_max) ) THEN
497  impurity%Z_AVE(:,iimp,izimp) = coreprofiles_iter%profiles_1d(1)%ion(ind)%state(izimp)%z_min
498  ELSE
499  diag%ERROR_INDEX = -1
500  diag%ERROR_MESSAGE = 'Could not determine the average state charge for ion '//num2str(ind)//' and state '//num2str(izimp)
501  return
502  END IF
503 
504  IF (ASSOCIATED(coreprofiles_iter%profiles_1d(1)%ion(ind)%state(izimp)%z_average_square_1d)) THEN
505  impurity%Z2_AVE(:,iimp,izimp) = coreprofiles_iter%profiles_1d(1)%ion(ind)%state(izimp)%z_average_square_1d(:)
506  ELSEIF (ids_is_valid(coreprofiles_iter%profiles_1d(1)%ion(ind)%state(izimp)%z_min) .AND. &
507  (coreprofiles_iter%profiles_1d(1)%ion(ind)%state(izimp)%z_min .EQ. &
508  coreprofiles_iter%profiles_1d(1)%ion(ind)%state(izimp)%z_max) ) THEN
509  impurity%Z2_AVE(:,iimp,izimp) = coreprofiles_iter%profiles_1d(1)%ion(ind)%state(izimp)%z_min**2
510  ELSE
511  diag%ERROR_INDEX = -1
512  diag%ERROR_MESSAGE = 'Could not determine the average state charge square for ion '//num2str(ind)//' and state '//num2str(izimp)
513  return
514  END IF
515 
516  IF (ASSOCIATED(coreprofiles_iter%profiles_1d(1)%ion(ind)%state(izimp)%density)) THEN
517  fun1(:) = coreprofiles_iter%profiles_1d(1)%ion(ind)%state(izimp)%density_thermal(:)
518  impurity%NZ(:,iimp,izimp) = coreprofiles_iter%profiles_1d(1)%ion(ind)%state(izimp)%density(:)
519  CALL l3deriv(fun1, &
520  coreprofiles_iter%profiles_1d(1)%grid%rho_tor_norm, nprof, &
521  dfun1, geometry%RHO_NORM, geometry%NRHO)
522  impurity%DNZ(:,iimp,izimp) = dfun1(:)/geometry%RHO(nprof)
523  END IF
524 
525  IF (ASSOCIATED(coreprofiles_old%profiles_1d(1)%ion(ind)%state(izimp)%density)) THEN
526  fun2(:) = coreprofiles_old%profiles_1d(1)%ion(ind)%state(izimp)%density_thermal(:)
527  CALL l3interp(fun2, &
528  coreprofiles_old%profiles_1d(1)%grid%rho_tor_norm, nprof_old, &
529  dfun2, geometry%RHO_NORM, geometry%NRHO)
530  impurity%NZM(:,iimp,izimp) = dfun2(:)
531  CALL l3deriv(fun2, &
532  coreprofiles_old%profiles_1d(1)%grid%rho_tor_norm, nprof_old, &
533  dfun2, geometry%RHO_NORM, geometry%NRHO)
534  impurity%DNZM(:,iimp,izimp) = dfun2(:)/geometry%RHO(nprof)
535  END IF
536 
537  IF (.NOT. acceptable_impurity_density(impurity%NZ(:,iimp,izimp))) then
538  diag%ERROR_INDEX = -1
539  diag%ERROR_MESSAGE = 'ERROR in CONVERT_IDS_TO_INTERNAL_TYPES, input IMPURITY%NZ(:,'//num2str(iimp)//','//num2str(izimp)//') includes the unacceptably small value '//num2str(minval(impurity%NZ(:,iimp,izimp)))
540  WRITE(*,*)diag%ERROR_MESSAGE
541  RETURN
542  END IF
543 
544  END DO loop_over_impurity_ionization_state
545 
546  impurity%NZ(:,iimp,0) = 0.0_ids_real
547  IF (coreprofiles_iter%profiles_1d(1)%ion(ind)%neutral_index.GE.1_ids_int) THEN
548  IF(ASSOCIATED(coreprofiles_iter%profiles_1d(1)%neutral(coreprofiles_iter%profiles_1d(1)%ion(ind)%neutral_index)%density)) THEN
549  impurity%NZ(:,iimp,0) = coreprofiles_iter%profiles_1d(1)%neutral(coreprofiles_iter%profiles_1d(1)%ion(ind)%neutral_index)%density(:)
550  ELSE IF(ASSOCIATED(coreprofiles_iter%profiles_1d(1)%neutral(coreprofiles_iter%profiles_1d(1)%ion(ind)%neutral_index)%density_thermal)) THEN
551  impurity%NZ(:,iimp,0) = coreprofiles_iter%profiles_1d(1)%neutral(coreprofiles_iter%profiles_1d(1)%ion(ind)%neutral_index)%density_thermal(:)
552  END IF
553  END IF
554 
555  END DO loop_over_input_impurities
556 
557  DEALLOCATE (fun1, dfun1, fun2, dfun2)
558 
559 
560 ! +++ Boundary conditions
561 ! name & index convention:
562 ! psi : N.A.
563 ! density : 0-electrons; 1:-ions
564 ! temperature : 0-electrons; 1:-ions
565 ! rotation : 0-electrons; 1:-ions
566 
567  profiles%PSI_BND_TYPE = 0_ids_int
568  profiles%NE_BND_TYPE = 0_ids_int
569  profiles%NI_BND_TYPE = 0_ids_int
570  impurity%NZ_BND_TYPE = 0_ids_int
571  profiles%TE_BND_TYPE = 0_ids_int
572  profiles%TI_BND_TYPE = 0_ids_int
573  profiles%VTOR_BND_TYPE = 0_ids_int
574 
575  nsol = SIZE(coresolver_iter%solver_1d(1)%equation)
576 
577  loop_over_equations_to_be_solved: DO isol = 1, nsol
578  check_solver_has_identifier_name: IF ( ASSOCIATED(coresolver_iter%solver_1d(1)%equation(isol)%primary_quantity%identifier%name) ) THEN
579 
580  compare_sources_identifier: IF (coresolver_iter%solver_1d(1)%equation(isol)%primary_quantity%identifier%name(1) .EQ. "psi") THEN
581  ibnd = 1
582  DO i = 1, SIZE(coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition)
583  IF ( coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(i)%position.GE. &
584  coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%position ) THEN
585  ibnd = i
586  END IF
587  END DO
588  CALL ids_to_bc_type( &
589  coresolver_iter%solver_1d(1)%equation(isol)%computation_mode%index, &
590  coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%type%index, &
591  profiles%PSI_BND_TYPE )
592  profiles%PSI_BND = coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%value
593  profiles%PSI_BND_RHO_NORM = coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%position
594 
596 
597  ELSE IF (coresolver_iter%solver_1d(1)%equation(isol)%primary_quantity%identifier%name(1) .EQ. "density") THEN
598 
599  IF (coresolver_iter%solver_1d(1)%equation(isol)%primary_quantity%ion_index .EQ. 0) THEN
600 
601  !! This is an electron equation (coded by "ion_index .EQ. 0")
602 
603  ibnd = 1
604  DO i = 1, SIZE(coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition)
605  IF ( coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(i)%position.GE. &
606  coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%position )THEN
607  ibnd = i
608  END IF
609  END DO
610  CALL ids_to_bc_type( &
611  coresolver_iter%solver_1d(1)%equation(isol)%computation_mode%index, &
612  coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%type%index, &
613  profiles%NE_BND_TYPE )
614  profiles%NE_BND = coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%value
615  profiles%NE_BND_RHO_NORM = coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%position
616 
618 
619  ELSE IF (coresolver_iter%solver_1d(1)%equation(isol)%primary_quantity%ion_index .GT. 0) THEN
620 
621  !! This is an ion equation (coded by "ion_index .GT. 0")
622 
623  ifound = 0
624  DO iion = 1, composition_iter%NION
625  IF (composition_iter%IONS(iion)%IDS_ion .EQ. &
626  coresolver_iter%solver_1d(1)%equation(isol)%primary_quantity%ion_index) THEN
627  ifound = 1
628  ibnd = 1
629  DO ib0 = 1, SIZE(coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition)
630  IF ( coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ib0 )%position.GE. &
631  coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%position ) THEN
632  ibnd = ib0
633  END IF
634  END DO
635  CALL ids_to_bc_type( &
636  coresolver_iter%solver_1d(1)%equation(isol)%computation_mode%index, &
637  coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%type%index, &
638  profiles%NI_BND_TYPE(iion) )
639  profiles%NI_BND(:,iion) = coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%value(:)
640  profiles%NI_BND_RHO_NORM(iion) = coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%position
641 
643 
644  END IF
645  END DO
646 ! IF (IFOUND.EQ.0) CALL REPORT_ERROR (-1, "density BC for none existing ion" , DIAG)
647  IF (diag%ERROR_INDEX.LT.0) goto 666
648 !
649  DO iimp = 1, composition_iter%NIMP
650  IF (composition_iter%IMPURITY(iimp)%IDS_ion .EQ. &
651  coresolver_iter%solver_1d(1)%equation(isol)%primary_quantity%ion_index) THEN
652  i_index = coresolver_iter%solver_1d(1)%equation(isol)%primary_quantity%ion_index
653  ifound = 0
654  DO izimp = 1, SIZE(composition_iter%IMPURITY(iimp)%IDS_state)
655  IF (composition_iter%IMPURITY(iimp)%IDS_state(izimp) .EQ. &
656  coresolver_iter%solver_1d(1)%equation(isol)%primary_quantity%state_index) THEN
657 
658  ifound = 1
659  ibnd = 1
660  DO i = 1, SIZE(coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition)
661  IF (ASSOCIATED(coresolver_iter%solver_1d(1)%equation(isol)%computation_mode%name)) THEN
662  IF (coresolver_iter%solver_1d(1)%equation(isol)%computation_mode%name(1).EQ."predictive".AND.&
663  coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(i)%position.GE. &
664  coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%position) &
665  ibnd = i
666  END IF
667  END DO
668  impurity%NZ_BND(:,iimp,izimp) = coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%value(:)
669  impurity%NZ_BND_TYPE(iimp) = coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%type%index
670  impurity%NZ_BND_RHO_NORM(iimp) = coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%position
671  END IF
672  END DO
673  END IF
674  END DO
675 ! IF (IFOUND.EQ.0) CALL REPORT_ERROR (-1, "temperature BC for none existing ion" , DIAG)
676  IF (diag%ERROR_INDEX.LT.0) goto 666
677 
678  ELSE
679 ! CALL REPORT_ERROR (1, "not applicable primary_quantity equation index", DIAG)
680 
681  END IF
682 
683  ELSE IF (coresolver_iter%solver_1d(1)%equation(isol)%primary_quantity%identifier%name(1) .EQ. "temperature") THEN
684 
685  IF (coresolver_iter%solver_1d(1)%equation(isol)%primary_quantity%ion_index.EQ. 0) THEN
686 
687  !! This is an electron equation (coded by "ion_index .EQ. 0")
688 
689  ibnd = 1
690  DO i = 1, SIZE(coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition)
691  IF ( coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(i)%position.GE. &
692  coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%position ) THEN
693  ibnd = i
694  END IF
695  END DO
696  CALL ids_to_bc_type( &
697  coresolver_iter%solver_1d(1)%equation(isol)%computation_mode%index, &
698  coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%type%index, &
699  profiles%TE_BND_TYPE )
700  profiles%TE_BND = coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%value
701  profiles%TE_BND_RHO_NORM = coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%position
702 
704 
705  ELSE IF (coresolver_iter%solver_1d(1)%equation(isol)%primary_quantity%ion_index .GT. 0) THEN
706 
707  !! This is an ion equation (coded by "ion_index .GT. 0")
708 
709  ifound = 0
710  DO iion = 1, composition_iter%NION
711  IF (composition_iter%IONS(iion)%IDS_ion .EQ. &
712  coresolver_iter%solver_1d(1)%equation(isol)%primary_quantity%ion_index) THEN
713  ifound = 1
714  ibnd = 1
715  DO ib0 = 1, SIZE(coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition)
716  IF (coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ib0 )%position.GE. &
717  coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%position ) THEN
718  ibnd = ib0
719  END IF
720  END DO
721  CALL ids_to_bc_type( &
722  coresolver_iter%solver_1d(1)%equation(isol)%computation_mode%index, &
723  coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%type%index, &
724  profiles%TI_BND_TYPE(iion) )
725  profiles%TI_BND(:,iion) = coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%value(:)
726  profiles%TI_BND_RHO_NORM(iion) = coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%position
727 
729 
730  END IF
731  END DO
732 ! IF (IFOUND.EQ.0) CALL REPORT_ERROR (-1, "temperature BC for none existing ion" , DIAG)
733  IF (diag%ERROR_INDEX.LT.0) goto 666
734 
735  ELSE
736 ! CALL REPORT_ERROR (1, "not applicable primary_quantity equation index", DIAG)
737 
738  END IF
739 
740 
741  ELSE IF (coresolver_iter%solver_1d(1)%equation(isol)%primary_quantity%identifier%name(1) .EQ. "rotation") THEN
742 
743  IF (coresolver_iter%solver_1d(1)%equation(isol)%primary_quantity%ion_index .GT. 0) THEN
744 
745  !! This is an ion equation (coded by "ion_index .GT. 0")
746 
747  ifound = 0
748  DO iion = 1, composition_iter%NION
749  IF (composition_iter%IONS(iion)%IDS_ion .EQ. &
750  coresolver_iter%solver_1d(1)%equation(isol)%primary_quantity%ion_index) THEN
751  ifound = 1
752  ibnd = 1
753  DO i = 1, SIZE(coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition)
754  IF (coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(i)%position.GE. &
755  coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%position) &
756  ibnd = i
757  END DO
758  CALL ids_to_bc_type( &
759  coresolver_iter%solver_1d(1)%equation(isol)%computation_mode%index, &
760  coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%type%index, &
761  profiles%VTOR_BND_TYPE(iion) )
762  profiles%VTOR_BND(:,iion) = coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%value(:)
763  profiles%VTOR_BND_RHO_NORM(iion) = coresolver_iter%solver_1d(1)%equation(isol)%boundary_condition(ibnd)%position
764 
766 
767  END IF
768  END DO
769 ! IF (IFOUND.EQ.0) CALL REPORT_ERROR (-1, "rotation BC for none existing ion" , DIAG)
770  IF (diag%ERROR_INDEX.LT.0) goto 666
771  ELSE
772 ! CALL REPORT_ERROR (1, "not applicable primary_quantity equation index", DIAG)
773 
774  END IF
775 
776  ELSE
777 ! CALL REPORT_ERROR (1, "not applicable primary_quantity equation name", DIAG)
778  END IF compare_sources_identifier
779 
780  END IF check_solver_has_identifier_name
781 
782  END DO loop_over_equations_to_be_solved
783 
784 
785 ! +++ Convert control parameters:
786  numerics%TAU = coresolver_iter%time_step%data(SIZE(coresolver_iter%time_step%data))
787 
788  numerics%CONTROL%SOLVER_TYPE = coresolver_iter%solver%index
789 
790  associate(cp => coresolver_iter%solver_1d(1)%control_parameters)
791 
792  call get_control_parameter(cp, 'hyperdiffusion_active' , ival, diag%ERROR_INDEX)
793  profiles%HYPERDIFFUSION%ACTIVE = (ival==1)
794  if (diag%ERROR_INDEX < 0) then
795  diag%ERROR_MESSAGE = 'in convert, error from get_control_parameter, failed to set hyperdiffusion_active'
796  return
797  end if
798  call get_control_parameter(cp, 'hyperdiffusion_density_terms' , profiles%HYPERDIFFUSION%DENS%TERMS, diag%ERROR_INDEX)
799  if (diag%ERROR_INDEX < 0) then
800  diag%ERROR_MESSAGE = 'in convert, error from get_control_parameter, failed to set hyperdiffusion_density_terms'
801  return
802  end if
803  call get_control_parameter(cp, 'hyperdiffusion_temperature_terms', profiles%HYPERDIFFUSION%TEMP%TERMS, diag%ERROR_INDEX)
804  if (diag%ERROR_INDEX < 0) then
805  diag%ERROR_MESSAGE = 'in convert, error from get_control_parameter, failed to set hyperdiffusion_temperature_terms'
806  return
807  end if
808  call get_control_parameter(cp, 'hyperdiffusion_rotation_terms' , profiles%HYPERDIFFUSION%ROTATION%TERMS, diag%ERROR_INDEX)
809  if (diag%ERROR_INDEX < 0) then
810  diag%ERROR_MESSAGE = 'in convert, error from get_control_parameter, failed to set hyperdiffusion_rotation_terms'
811  return
812  end if
813 
814  call get_control_parameter(cp, 'hyperdiffusion_density_implicit' , profiles%HYPERDIFFUSION%DENS%IMPLICIT, diag%ERROR_INDEX)
815  if (diag%ERROR_INDEX < 0) then
816  diag%ERROR_MESSAGE = 'in convert, error from get_control_parameter, failed to set hyperdiffusion_density_implicit'
817  return
818  end if
819  call get_control_parameter(cp, 'hyperdiffusion_density_explicit' , profiles%HYPERDIFFUSION%DENS%EXPLICIT, diag%ERROR_INDEX)
820  if (diag%ERROR_INDEX < 0) then
821  diag%ERROR_MESSAGE = 'in convert, error from get_control_parameter, failed to set hyperdiffusion_density_explicit'
822  return
823  end if
824  call get_control_parameter(cp, 'hyperdiffusion_density_rho_tor_norm_cutoff' , profiles%HYPERDIFFUSION%DENS%RHO_NORM_CUTOFF, diag%ERROR_INDEX)
825  if (diag%ERROR_INDEX < 0) then
826  diag%ERROR_MESSAGE = 'in convert, error from get_control_parameter, failed to set hyperdiffusion_density_rho_tor_norm_cutoff'
827  return
828  end if
829  call get_control_parameter(cp, 'hyperdiffusion_temperature_implicit' , profiles%HYPERDIFFUSION%TEMP%IMPLICIT, diag%ERROR_INDEX)
830  if (diag%ERROR_INDEX < 0) then
831  diag%ERROR_MESSAGE = 'in convert, error from get_control_parameter, failed to set hyperdiffusion_temperature_implicit'
832  return
833  end if
834  call get_control_parameter(cp, 'hyperdiffusion_temperature_explicit' , profiles%HYPERDIFFUSION%TEMP%EXPLICIT, diag%ERROR_INDEX)
835  if (diag%ERROR_INDEX < 0) then
836  diag%ERROR_MESSAGE = 'in convert, error from get_control_parameter, failed to set hyperdiffusion_temperature_explicit'
837  return
838  end if
839  call get_control_parameter(cp, 'hyperdiffusion_temperature_rho_tor_norm_cutoff' , profiles%HYPERDIFFUSION%TEMP%RHO_NORM_CUTOFF, diag%ERROR_INDEX)
840  if (diag%ERROR_INDEX < 0) then
841  diag%ERROR_MESSAGE = 'in convert, error from get_control_parameter, failed to set hyperdiffusion_temperature_rho_tor_norm_cutoff'
842  return
843  end if
844  call get_control_parameter(cp, 'hyperdiffusion_rotation_implicit' , profiles%HYPERDIFFUSION%ROTATION%IMPLICIT, diag%ERROR_INDEX)
845  if (diag%ERROR_INDEX < 0) then
846  diag%ERROR_MESSAGE = 'in convert, error from get_control_parameter, failed to set hyperdiffusion_rotation_implicit'
847  return
848  end if
849  call get_control_parameter(cp, 'hyperdiffusion_rotation_explicit' , profiles%HYPERDIFFUSION%ROTATION%EXPLICIT, diag%ERROR_INDEX)
850  if (diag%ERROR_INDEX < 0) then
851  diag%ERROR_MESSAGE = 'in convert, error from get_control_parameter, failed to set hyperdiffusion_rotation_explicit'
852  return
853  end if
854  call get_control_parameter(cp, 'hyperdiffusion_rotation_rho_tor_norm_cutoff' , profiles%HYPERDIFFUSION%ROTATION%RHO_NORM_CUTOFF, diag%ERROR_INDEX)
855  if (diag%ERROR_INDEX < 0) then
856  diag%ERROR_MESSAGE = 'in convert, error from get_control_parameter, failed to set '
857  return
858  end if
859 
860 
861 
862  call get_control_parameter(cp, 'mixing_ratio_kinetic_profiles',numerics%CONTROL%AMIX, diag%ERROR_INDEX)
863  if (diag%ERROR_INDEX < 0) then
864  diag%ERROR_MESSAGE = 'in convert, error from get_control_parameter, failed to read mixing_ratio_kinetic_profiles'
865  return
866  end if
867 
868  call get_control_parameter(cp, 'mixing_ratio_transport',numerics%CONTROL%AMIX_TRANSP, diag%ERROR_INDEX)
869  if (diag%ERROR_INDEX < 0) then
870  diag%ERROR_MESSAGE = 'in convert, error from get_control_parameter_int, failed to read mixing_ratio_transport'
871  return
872  end if
873 
874  call get_control_parameter(cp, 'mixing_ratio_sources',numerics%CONTROL%AMIX_SRC, diag%ERROR_INDEX)
875  if (diag%ERROR_INDEX < 0) then
876  diag%ERROR_MESSAGE = 'in convert, error from get_control_parameter_int, failed to read mixing_ratio_sources'
877  return
878  end if
879 
880  call get_control_parameter(cp, 'ohmic_power_multiplier',numerics%CONTROL%MULTIPLIER_OH, diag%ERROR_INDEX)
881  if (diag%ERROR_INDEX < 0) then
882  diag%ERROR_MESSAGE = 'in convert, error from get_control_parameter_int, failed to read ohmic_power_multiplier'
883  return
884  end if
885  end associate
886 
887 ! +++ Convert transport:
888  transport%C1(1) = 0.0e0_ids_real
889  transport%C1(2) = 1.5e0_ids_real
890  transport%C1(3) = 2.5e0_ids_real
891 
892  transport%SIGMA = 0.0e0_ids_real
893  transport%DIFF_NE = 0.0e0_ids_real
894  transport%VCONV_NE = 0.0e0_ids_real
895  transport%DIFF_NI = 0.0e0_ids_real
896  transport%VCONV_NI = 0.0e0_ids_real
897  transport%DIFF_TE = 0.0e0_ids_real
898  transport%VCONV_TE = 0.0e0_ids_real
899  transport%DIFF_TI = 0.0e0_ids_real
900  transport%VCONV_TI = 0.0e0_ids_real
901  transport%DIFF_VTOR = 0.0e0_ids_real
902  transport%VCONV_VTOR = 0.0e0_ids_real
903 
904  CALL get_transp_composition(coretransport_iter, ref_combiner_transp, composition_transp)
905 
906  loop_over_models: DO imod = 1, SIZE(coretransport_iter%MODEL)
907 
908  ndiff = SIZE(coretransport_iter%MODEL(imod)%profiles_1d(1)%grid_d%rho_tor_norm)
909  nconv = SIZE(coretransport_iter%MODEL(imod)%profiles_1d(1)%grid_v%rho_tor_norm)
910  nflux = SIZE(coretransport_iter%MODEL(imod)%profiles_1d(1)%grid_flux%rho_tor_norm)
911 
912  ALLOCATE (fun2(ndiff),fun3(nconv),fun4(nflux), dfun2(nprof))
913 
914  check_combined_models: IF (trim(adjustl(coretransport_iter%MODEL(imod)%identifier%name(1))) .EQ. trim(adjustl(ref_combiner_transp))) THEN
915 
916  found_parallel_conductivity: IF (ASSOCIATED(coretransport_iter%MODEL(imod)%profiles_1d(1)%conductivity_parallel)) THEN
917  fun2(:) = coretransport_iter%MODEL(imod)%profiles_1d(1)%conductivity_parallel(:)
918  CALL l3interp(fun2, coretransport_iter%MODEL(imod)%profiles_1d(1)%grid_d%rho_tor_norm, ndiff, dfun2, geometry%RHO_NORM, nprof)
919  transport%SIGMA = transport%SIGMA + dfun2
920  END IF found_parallel_conductivity
921 
922 
923  found_electron_particle_duffusion: IF (ASSOCIATED(coretransport_iter%MODEL(imod)%profiles_1d(1)%electrons%particles%d)) THEN
924  fun2(:) = coretransport_iter%MODEL(imod)%profiles_1d(1)%electrons%particles%d(:)
925  CALL lininterp(fun2, coretransport_iter%MODEL(imod)%profiles_1d(1)%grid_d%rho_tor_norm, ndiff, dfun2, geometry%RHO_NORM, nprof)
926  IF(coretransport_iter%MODEL(imod)%flux_multiplier .EQ. 0.0_ids_real) transport%DIFF_NE_MODEL(:,1) = transport%DIFF_NE_MODEL(:,1) + dfun2
927  IF(coretransport_iter%MODEL(imod)%flux_multiplier .EQ. 1.5_ids_real) transport%DIFF_NE_MODEL(:,2) = transport%DIFF_NE_MODEL(:,2) + dfun2
928  IF(coretransport_iter%MODEL(imod)%flux_multiplier .EQ. 2.5_ids_real) transport%DIFF_NE_MODEL(:,3) = transport%DIFF_NE_MODEL(:,3) + dfun2
929  transport%DIFF_NE(:) = transport%DIFF_NE(:) + dfun2
930  END IF found_electron_particle_duffusion
931 
932 
933  found_electron_particle_pinch: IF (ASSOCIATED(coretransport_iter%MODEL(imod)%profiles_1d(1)%electrons%particles%v)) THEN
934  fun3(:) = coretransport_iter%MODEL(imod)%profiles_1d(1)%electrons%particles%v(:)
935  CALL lininterp(fun3, coretransport_iter%MODEL(imod)%profiles_1d(1)%grid_v%rho_tor_norm, nconv, dfun2, geometry%RHO_NORM, nprof)
936  IF(coretransport_iter%MODEL(imod)%flux_multiplier .EQ. 0.0_ids_real) transport%VCONV_NE_MODEL(:,1) = transport%VCONV_NE_MODEL(:,1) + dfun2
937  IF(coretransport_iter%MODEL(imod)%flux_multiplier .EQ. 1.5_ids_real) transport%VCONV_NE_MODEL(:,2) = transport%VCONV_NE_MODEL(:,2) + dfun2
938  IF(coretransport_iter%MODEL(imod)%flux_multiplier .EQ. 2.5_ids_real) transport%VCONV_NE_MODEL(:,3) = transport%VCONV_NE_MODEL(:,3) + dfun2
939  transport%VCONV_NE(:) = transport%VCONV_NE(:) + dfun2
940  END IF found_electron_particle_pinch
941 
942 
943  found_electron_energy_duffusion:IF (ASSOCIATED(coretransport_iter%MODEL(imod)%profiles_1d(1)%electrons%energy%d)) THEN
944  fun2(:) = coretransport_iter%MODEL(imod)%profiles_1d(1)%electrons%energy%d(:)
945  CALL lininterp(fun2, coretransport_iter%MODEL(imod)%profiles_1d(1)%grid_d%rho_tor_norm, ndiff, dfun2, geometry%RHO_NORM, nprof)
946  transport%DIFF_TE = transport%DIFF_TE + dfun2
947  END IF found_electron_energy_duffusion
948 
949 
950  found_electron_energy_pinch: IF (ASSOCIATED(coretransport_iter%MODEL(imod)%profiles_1d(1)%electrons%energy%v)) THEN
951  fun3(:) = coretransport_iter%MODEL(imod)%profiles_1d(1)%electrons%energy%v(:)
952  CALL lininterp(fun3, coretransport_iter%MODEL(imod)%profiles_1d(1)%grid_v%rho_tor_norm, nconv, dfun2, geometry%RHO_NORM, nprof)
953  transport%VCONV_TE = transport%VCONV_TE + dfun2
954  END IF found_electron_energy_pinch
955 
956 
957  loop_over_ets_ions: DO iion = 1, composition_iter%NION
958 
959  CALL find_ion(iion, composition_iter, composition_transp, ion_index)
960 
961  loop_over_found_ions: DO i = 1, SIZE(ion_index)
962  ions_match: IF (ion_index(i).GT.0 .AND. composition_transp%IONS(i)%IDS_model.EQ.imod) THEN
963  ids_iion = composition_transp%IONS(ion_index(i))%IDS_ion
964 
965  found_ion_particle_duffusion: IF (ASSOCIATED(coretransport_iter%MODEL(imod)%profiles_1d(1)%ion(ids_iion)%particles%d)) THEN
966  fun2(:) = coretransport_iter%MODEL(imod)%profiles_1d(1)%ion(ids_iion)%particles%d(:)
967  CALL lininterp(fun2, coretransport_iter%MODEL(imod)%profiles_1d(1)%grid_d%rho_tor_norm, ndiff, dfun2, geometry%RHO_NORM, nprof)
968  IF(coretransport_iter%MODEL(imod)%flux_multiplier .EQ. 0.0_ids_real) transport%DIFF_NI_MODEL(:,iion,1) = transport%DIFF_NI_MODEL(:,iion,1) + dfun2
969  IF(coretransport_iter%MODEL(imod)%flux_multiplier .EQ. 1.5_ids_real) transport%DIFF_NI_MODEL(:,iion,2) = transport%DIFF_NI_MODEL(:,iion,2) + dfun2
970  IF(coretransport_iter%MODEL(imod)%flux_multiplier .EQ. 2.5_ids_real) transport%DIFF_NI_MODEL(:,iion,3) = transport%DIFF_NI_MODEL(:,iion,3) + dfun2
971  transport%DIFF_NI(:,iion) = transport%DIFF_NI(:,iion) + dfun2
972  END IF found_ion_particle_duffusion
973 
974 
975  found_ion_particle_pinch: IF (ASSOCIATED(coretransport_iter%MODEL(imod)%profiles_1d(1)%ion(ids_iion)%particles%v)) THEN
976  fun3(:) = coretransport_iter%MODEL(imod)%profiles_1d(1)%ion(ids_iion)%particles%v(:)
977  CALL lininterp(fun3, coretransport_iter%MODEL(imod)%profiles_1d(1)%grid_v%rho_tor_norm, nconv, dfun2, geometry%RHO_NORM, nprof)
978  IF(coretransport_iter%MODEL(imod)%flux_multiplier .EQ. 0.0_ids_real) transport%VCONV_NI_MODEL(:,iion,1) = transport%VCONV_NI_MODEL(:,iion,1) + dfun2
979  IF(coretransport_iter%MODEL(imod)%flux_multiplier .EQ. 1.5_ids_real) transport%VCONV_NI_MODEL(:,iion,2) = transport%VCONV_NI_MODEL(:,iion,2) + dfun2
980  IF(coretransport_iter%MODEL(imod)%flux_multiplier .EQ. 2.5_ids_real) transport%VCONV_NI_MODEL(:,iion,3) = transport%VCONV_NI_MODEL(:,iion,3) + dfun2
981  transport%VCONV_NI(:,iion) = transport%VCONV_NI(:,iion) + dfun2
982  END IF found_ion_particle_pinch
983 
984 
985  found_ion_energy_duffusion: IF (ASSOCIATED(coretransport_iter%MODEL(imod)%profiles_1d(1)%ion(ids_iion)%energy%d)) THEN
986  fun2(:) = coretransport_iter%MODEL(imod)%profiles_1d(1)%ion(ids_iion)%energy%d(:)
987  CALL lininterp(fun2, coretransport_iter%MODEL(imod)%profiles_1d(1)%grid_d%rho_tor_norm, ndiff, dfun2, geometry%RHO_NORM, nprof)
988  transport%DIFF_TI(:,iion) = transport%DIFF_TI(:,iion) + dfun2
989  END IF found_ion_energy_duffusion
990 
991 
992  found_ion_energy_pinch: IF (ASSOCIATED(coretransport_iter%MODEL(imod)%profiles_1d(1)%ion(ids_iion)%energy%v)) THEN
993  fun3(:) = coretransport_iter%MODEL(imod)%profiles_1d(1)%ion(ids_iion)%energy%v(:)
994  CALL lininterp(fun3, coretransport_iter%MODEL(imod)%profiles_1d(1)%grid_v%rho_tor_norm, nconv, dfun2, geometry%RHO_NORM, nprof)
995  transport%VCONV_TI(:,iion) = transport%VCONV_TI(:,iion) + dfun2
996  END IF found_ion_energy_pinch
997 
998 
999 ! NO INDIVIDUAL TOROIDAL VELOCYTY TRANSPORT IMPLEMENTED IN IDS:
1000  transport%DIFF_VTOR = 0.0e0_ids_real
1001  transport%VCONV_VTOR = 0.0e0_ids_real
1002 
1003 
1004  END IF ions_match
1005 
1006  END DO loop_over_found_ions
1007 
1008  IF (ALLOCATED(ion_index)) DEALLOCATE(ion_index)
1009 
1010  END DO loop_over_ets_ions
1011 
1012 ! NO EXCHANGE TERMS IMPLEMENTED IN IDS:
1013  transport%QGI = 0.e0_ids_real
1014  transport%QGE = sum(transport%QGI, dim=2)
1015 
1016 
1017  END IF check_combined_models
1018 
1019 
1020 
1021 
1022 ! +++ Impurity fluxes/transport coefficients
1023  check_impurity_models: IF (trim(adjustl(coretransport_iter%MODEL(imod)%identifier%name(1))) .EQ. trim(adjustl(ref_combiner_transp))) THEN
1024 
1025  loop_over_ets_impurities: DO iimp = 1, composition_iter%NIMP
1026 
1027  CALL find_impurity(iimp, composition_iter, composition_transp, imp_index)
1028 
1029  loop_over_found_impurities: DO i = 1, SIZE(imp_index)
1030  impurities_match: IF (imp_index(i).GT.0 .AND. composition_transp%IMPURITY(imp_index(i))%IDS_model.EQ.imod) THEN
1031 
1032  ids_iion = composition_transp%IMPURITY(imp_index(i))%IDS_ion
1033 
1034  IF(ASSOCIATED(coretransport_iter%MODEL(imod)%profiles_1d(1)%ion(ids_iion)%state)) THEN
1035  loop_over_impurity_ionization_states: DO izimp = 1, SIZE(composition_transp%IMPURITY(iimp)%IDS_state)
1036  dfun2(:) = 0.0_ids_real
1037  IF(ASSOCIATED(coretransport_iter%MODEL(imod)%profiles_1d(1)%ion(ids_iion)%state(izimp)%particles%d)) THEN
1038  fun2(:) = coretransport_iter%MODEL(imod)%profiles_1d(1)%ion(ids_iion)%state(izimp)%particles%d(:)
1039  CALL l3interp(fun2, coretransport_iter%MODEL(imod)%profiles_1d(1)%grid_d%rho_tor_norm, ndiff, &
1040  dfun2, geometry%RHO_NORM, geometry%NRHO)
1041  END IF
1042  impurity%DIFF_NZ(:,iimp,izimp) = dfun2(:)
1043 
1044  dfun2(:) = 0.0_ids_real
1045  IF(ASSOCIATED(coretransport_iter%MODEL(imod)%profiles_1d(1)%ion(ids_iion)%state(izimp)%particles%v)) THEN
1046  fun3(:) = coretransport_iter%MODEL(imod)%profiles_1d(1)%ion(ids_iion)%state(izimp)%particles%v(:)
1047  CALL l3interp(fun3, coretransport_iter%MODEL(imod)%profiles_1d(1)%grid_v%rho_tor_norm, nconv, &
1048  dfun2, geometry%RHO_NORM, geometry%NRHO)
1049  END IF
1050  impurity%VCONV_NZ(:,iimp,izimp) = dfun2(:)
1051 
1052  dfun2(:) = 0.0_ids_real
1053  IF(ASSOCIATED(coretransport_iter%MODEL(imod)%profiles_1d(1)%ion(ids_iion)%state(izimp)%particles%flux)) THEN
1054  fun4(:) = coretransport_iter%MODEL(imod)%profiles_1d(1)%ion(ids_iion)%state(izimp)%particles%flux(:)
1055  CALL l3interp(fun4, coretransport_iter%MODEL(imod)%profiles_1d(1)%grid_flux%rho_tor_norm, nflux, &
1056  dfun2, geometry%RHO_NORM, geometry%NRHO)
1057  END IF
1058  impurity%FLUX_NZ(:,iimp,izimp) = dfun2(:)
1059 
1060  END DO loop_over_impurity_ionization_states
1061  END IF
1062  END IF impurities_match
1063 
1064  END DO loop_over_found_impurities
1065 
1066  IF (ALLOCATED(imp_index)) DEALLOCATE(imp_index)
1067 
1068  END DO loop_over_ets_impurities
1069 
1070 
1071  END IF check_impurity_models
1072 
1073  IF (ALLOCATED(fun2)) DEALLOCATE (fun2)
1074  IF (ALLOCATED(fun3)) DEALLOCATE (fun3)
1075  IF (ALLOCATED(fun4)) DEALLOCATE (fun4)
1076  IF (ALLOCATED(dfun2)) DEALLOCATE (dfun2)
1077 
1078  END DO loop_over_models
1079 
1080 
1081 
1082 
1083 
1084 ! +++ Convert sources:
1085  sources%SIGMA = 0.e0_ids_real
1086  sources%CURR_EXP = 0.e0_ids_real
1087  sources%CURR_IMP = 0.e0_ids_real
1088  sources%QE_EXP = 0.e0_ids_real
1089  sources%QE_IMP = 0.e0_ids_real
1090  sources%SE_EXP = 0.e0_ids_real
1091  sources%SE_IMP = 0.e0_ids_real
1092  sources%SI_EXP = 0.e0_ids_real
1093  sources%SI_IMP = 0.e0_ids_real
1094  sources%QI_EXP = 0.e0_ids_real
1095  sources%QI_IMP = 0.e0_ids_real
1096  sources%UI_EXP = 0.e0_ids_real
1097  sources%UI_IMP = 0.e0_ids_real
1098 
1099 
1100  CALL get_source_composition(coresource_iter, ref_combiner_src, composition_source)
1101 
1102  loop_over_sources: DO isrc = 1, SIZE(coresource_iter%SOURCE)
1103 
1104  nsrc = SIZE(coresource_iter%SOURCE(isrc)%profiles_1d(1)%grid%rho_tor_norm)
1105 
1106  ALLOCATE (fun2(nsrc), dfun2(nprof))
1107 
1108  check_combined_sources: IF (trim(adjustl(coresource_iter%SOURCE(isrc)%identifier%name(1))) .EQ. trim(adjustl(ref_combiner_src))) THEN
1109 
1110  found_parallel_conductivity2: IF (ASSOCIATED(coresource_iter%SOURCE(isrc)%profiles_1d(1)%conductivity_parallel)) THEN
1111  fun2(:) = coresource_iter%SOURCE(isrc)%profiles_1d(1)%conductivity_parallel(:)
1112  CALL l3interp(fun2, coresource_iter%SOURCE(isrc)%profiles_1d(1)%grid%rho_tor_norm, nsrc, dfun2, geometry%RHO_NORM, nprof)
1113  sources%SIGMA = sources%SIGMA + dfun2
1114  END IF found_parallel_conductivity2
1115 
1116 
1117  found_parallel_current: IF (ASSOCIATED(coresource_iter%SOURCE(isrc)%profiles_1d(1)%j_parallel)) THEN
1118  fun2(:) = coresource_iter%SOURCE(isrc)%profiles_1d(1)%j_parallel(:)
1119  CALL l3interp(fun2, coresource_iter%SOURCE(isrc)%profiles_1d(1)%grid%rho_tor_norm, nsrc, dfun2, geometry%RHO_NORM, nprof)
1120  sources%CURR_EXP = sources%CURR_EXP + dfun2
1121  END IF found_parallel_current
1122 
1123 
1124  found_electron_particle_exp_source: IF (ASSOCIATED(coresource_iter%SOURCE(isrc)%profiles_1d(1)%electrons%particles)) THEN
1125  fun2(:) = coresource_iter%SOURCE(isrc)%profiles_1d(1)%electrons%particles(:)
1126  CALL l3interp(fun2, coresource_iter%SOURCE(isrc)%profiles_1d(1)%grid%rho_tor_norm, nsrc, dfun2, geometry%RHO_NORM, nprof)
1127  sources%SE_EXP = sources%SE_EXP(:) + dfun2
1128  END IF found_electron_particle_exp_source
1129 
1130 
1131  found_electron_energy_exp_source: IF (ASSOCIATED(coresource_iter%SOURCE(isrc)%profiles_1d(1)%electrons%energy)) THEN
1132  fun2(:) = coresource_iter%SOURCE(isrc)%profiles_1d(1)%electrons%energy(:) / imas_constants%ev
1133  CALL l3interp(fun2, coresource_iter%SOURCE(isrc)%profiles_1d(1)%grid%rho_tor_norm, nsrc, dfun2, geometry%RHO_NORM, nprof)
1134  sources%QE_EXP = sources%QE_EXP(:) + dfun2
1135  END IF found_electron_energy_exp_source
1136 
1137 
1138  loop_over_ets_ions3: DO iion = 1, composition_iter%NION
1139 
1140  CALL find_ion(iion, composition_iter, composition_source, ion_index)
1141 
1142  loop_over_found_coresource_ions: DO i = 1, SIZE(ion_index)
1143  ions_match2: IF (ion_index(i).GT.0 .AND. composition_source%IONS(i)%IDS_model.EQ.isrc) THEN
1144  ids_iion = composition_source%IONS(ion_index(i))%IDS_ion
1145 
1146  found_ion_particle_source: IF (ASSOCIATED(coresource_iter%SOURCE(isrc)%profiles_1d(1)%ion(ids_iion)%particles)) THEN
1147  fun2(:) = coresource_iter%SOURCE(isrc)%profiles_1d(1)%ion(ids_iion)%particles(:)
1148  CALL l3interp(fun2, coresource_iter%SOURCE(isrc)%profiles_1d(1)%grid%rho_tor_norm, nsrc, dfun2, geometry%RHO_NORM, nprof)
1149  sources%SI_EXP(:,iion) = sources%SI_EXP(:,iion) + dfun2
1150  END IF found_ion_particle_source
1151 
1152 
1153  found_ion_energy_source: IF (ASSOCIATED(coresource_iter%SOURCE(isrc)%profiles_1d(1)%ion(ids_iion)%energy)) THEN
1154  fun2(:) = coresource_iter%SOURCE(isrc)%profiles_1d(1)%ion(ids_iion)%energy(:) / imas_constants%ev
1155  CALL l3interp(fun2, coresource_iter%SOURCE(isrc)%profiles_1d(1)%grid%rho_tor_norm, nsrc, dfun2, geometry%RHO_NORM, nprof)
1156  sources%QI_EXP(:,iion) = sources%QI_EXP(:,iion) + dfun2
1157  END IF found_ion_energy_source
1158 
1159 
1160  END IF ions_match2
1161 
1162  END DO loop_over_found_coresource_ions
1163 
1164  IF (ALLOCATED(ion_index)) DEALLOCATE(ion_index)
1165 
1166  END DO loop_over_ets_ions3
1167 
1168 
1169  loop_over_ets_impurities2: DO iimp = 1, composition_iter%NIMP
1170 
1171  CALL find_impurity(iimp, composition_iter, composition_source, imp_index)
1172 
1173  loop_over_found_impurities2: DO i = 1, SIZE(imp_index)
1174  impurities_match2: IF (imp_index(i).GT.0 .AND. composition_source%IMPURITY(imp_index(i))%IDS_model.EQ.imod) THEN
1175 
1176  ids_iion = composition_source%IMPURITY(imp_index(i))%IDS_ion
1177 
1178  IF(ASSOCIATED(coresource_iter%SOURCE(isrc)%profiles_1d(1)%ion(ids_iion)%state)) THEN
1179  loop_over_impurity_ionization_states2: DO izimp = 1, SIZE(composition_source%IMPURITY(iimp)%IDS_state)
1180  dfun2(:) = 0.0_ids_real
1181  IF(ASSOCIATED(coresource_iter%SOURCE(isrc)%profiles_1d(1)%ion(ids_iion)%state(izimp)%particles)) THEN
1182  fun2(:) = coresource_iter%SOURCE(isrc)%profiles_1d(1)%ion(ids_iion)%state(izimp)%particles(:)
1183  CALL l3interp(fun2, coresource_iter%SOURCE(isrc)%profiles_1d(1)%grid%rho_tor_norm, nsrc, &
1184  dfun2, geometry%RHO_NORM, geometry%NRHO)
1185  END IF
1186  impurity%SOURCE_NZ(:,iimp,izimp) = dfun2(:)
1187 
1188  dfun2(:) = 0.0_ids_real
1189  IF(ASSOCIATED(coresource_iter%SOURCE(isrc)%profiles_1d(1)%ion(ids_iion)%state(izimp)%particles_decomposed%explicit_part)) THEN
1190  fun2(:) = coresource_iter%SOURCE(isrc)%profiles_1d(1)%ion(ids_iion)%state(izimp)%particles_decomposed%explicit_part(:)
1191  CALL l3interp(fun2, coresource_iter%SOURCE(isrc)%profiles_1d(1)%grid%rho_tor_norm, nsrc, &
1192  dfun2, geometry%RHO_NORM, geometry%NRHO)
1193  END IF
1194  impurity%SOURCE_NZ_EXP(:,iimp,izimp) = dfun2(:)
1195 
1196  dfun2(:) = 0.0_ids_real
1197  IF(ASSOCIATED(coresource_iter%SOURCE(isrc)%profiles_1d(1)%ion(ids_iion)%state(izimp)%particles_decomposed%implicit_part)) THEN
1198  fun2(:) = coresource_iter%SOURCE(isrc)%profiles_1d(1)%ion(ids_iion)%state(izimp)%particles_decomposed%implicit_part(:)
1199  CALL l3interp(fun2, coresource_iter%SOURCE(isrc)%profiles_1d(1)%grid%rho_tor_norm, nsrc, &
1200  dfun2, geometry%RHO_NORM, geometry%NRHO)
1201  END IF
1202  impurity%SOURCE_NZ_IMP(:,iimp,izimp) = dfun2(:)
1203 
1204  IF(.NOT.ASSOCIATED(coresource_iter%SOURCE(isrc)%profiles_1d(1)%ion(ids_iion)%state(izimp)%particles_decomposed%explicit_part) .OR. &
1205  .NOT.ASSOCIATED(coresource_iter%SOURCE(isrc)%profiles_1d(1)%ion(ids_iion)%state(izimp)%particles_decomposed%implicit_part)) THEN
1206 
1207  impurity%SOURCE_NZ_EXP(:,iimp,izimp) = impurity%SOURCE_NZ(:,iimp,izimp)
1208  impurity%SOURCE_NZ_IMP(:,iimp,izimp) = 0.0_ids_real
1209 
1210  END IF
1211 
1212  END DO loop_over_impurity_ionization_states2
1213  END IF
1214  END IF impurities_match2
1215 
1216  END DO loop_over_found_impurities2
1217 
1218  IF (ALLOCATED(imp_index)) DEALLOCATE(imp_index)
1219 
1220  END DO loop_over_ets_impurities2
1221 
1222 
1223 ! irena
1224 
1225  loop_over_ets_impurities3: DO iimp = 1, composition_iter%NIMP
1226 
1227  CALL find_impurity(iimp, composition_iter, composition_source, imp_index)
1228 
1229  loop_over_found_impurity_sources: DO i = 1, SIZE(imp_index)
1230  impurities_match3: IF (imp_index(i).GT.0 .AND. composition_source%IMPURITY(imp_index(i))%IDS_model.EQ.isrc) THEN
1231 
1232  ids_iion = composition_source%IMPURITY(imp_index(i))%IDS_ion
1233 
1234  check_ion_sources: IF (trim(adjustl(coresource_iter%SOURCE(isrc)%identifier%name(1))) .EQ. 'ionisation') THEN
1235 
1236  IF(ASSOCIATED(coresource_iter%SOURCE(isrc)%profiles_1d(1)%ion(ids_iion)%state)) THEN
1237  loop_over_impurity_ionization_sources: DO izimp = 1, SIZE(composition_source%IMPURITY(iimp)%IDS_state)
1238  dfun2(:) = 0.0_ids_real
1239  IF(ASSOCIATED(coresource_iter%SOURCE(isrc)%profiles_1d(1)%ion(ids_iion)%state(izimp)%particles_decomposed%explicit_part)) THEN
1240  fun2(:) = coresource_iter%SOURCE(isrc)%profiles_1d(1)%ion(ids_iion)%state(izimp)%particles_decomposed%explicit_part
1241  CALL l3interp(fun2, coresource_iter%SOURCE(isrc)%profiles_1d(1)%grid%rho_tor_norm, nsrc, &
1242  dfun2, geometry%RHO_NORM, geometry%NRHO)
1243  END IF
1244  impurity%SOURCE_JON_NZ_IMP(:,iimp,izimp) = dfun2(:)
1245  END DO loop_over_impurity_ionization_sources
1246  END IF
1247  END IF check_ion_sources
1248 
1249 
1250  check_rec_sources: IF (trim(adjustl(coresource_iter%SOURCE(isrc)%identifier%name(1))) .EQ. 'recombination') THEN
1251 
1252  IF(ASSOCIATED(coresource_iter%SOURCE(isrc)%profiles_1d(1)%ion(ids_iion)%state)) THEN
1253  loop_over_impurity_recombination_sources: DO izimp = 1, SIZE(composition_source%IMPURITY(iimp)%IDS_state)
1254  dfun2(:) = 0.0_ids_real
1255  IF(ASSOCIATED(coresource_iter%SOURCE(isrc)%profiles_1d(1)%ion(ids_iion)%state(izimp)%particles_decomposed%explicit_part)) THEN
1256  fun2(:) = coresource_iter%SOURCE(isrc)%profiles_1d(1)%ion(ids_iion)%state(izimp)%particles_decomposed%explicit_part
1257  CALL l3interp(fun2, coresource_iter%SOURCE(isrc)%profiles_1d(1)%grid%rho_tor_norm, nsrc, &
1258  dfun2, geometry%RHO_NORM, geometry%NRHO)
1259  END IF
1260  impurity%SOURCE_REC_NZ_IMP(:,iimp,izimp) = dfun2(:)
1261  END DO loop_over_impurity_recombination_sources
1262 
1263  END IF
1264  END IF check_rec_sources
1265 
1266 
1267 !! Thomas: These should already be zero!
1268 !!$ IF(ASSOCIATED(CORESOURCE_ITER%SOURCE(ISRC)%profiles_1d(1)%ion(IDS_IION)%state)) THEN
1269 !!$ IF(.NOT.ASSOCIATED(CORESOURCE_ITER%SOURCE(ISRC)%profiles_1d(1)%ion(IDS_IION)%state(IZIMP)%particles_decomposed%explicit_part) .OR. &
1270 !!$ .NOT.ASSOCIATED(CORESOURCE_ITER%SOURCE(ISRC)%profiles_1d(1)%ion(IDS_IION)%state(IZIMP)%particles_decomposed%explicit_part)) THEN
1271 !!$
1272 !!$ IMPURITY%SOURCE_REC_NZ_IMP(:,IIMP,IZIMP) = 0.0_IDS_REAL
1273 !!$ IMPURITY%SOURCE_JON_NZ_IMP(:,IIMP,IZIMP) = 0.0_IDS_REAL
1274 !!$ END IF
1275 !!$ END IF
1276 
1277 
1278  END IF impurities_match3
1279 
1280  END DO loop_over_found_impurity_sources
1281  END DO loop_over_ets_impurities3
1282 
1283 
1284 
1285 ! irena
1286 
1287  END IF check_combined_sources
1288 
1289 
1290  check_source_has_identifier_neoclassical: IF ( ASSOCIATED(coresource_iter%SOURCE(isrc)%identifier%name) ) THEN
1291  check_neoclassical_sources: IF ( (trim(adjustl(coresource_iter%SOURCE(isrc)%identifier%name(1))) .EQ. 'bootstrap_current') .OR. &
1292  (trim(adjustl(coresource_iter%SOURCE(isrc)%identifier%name(1))) .EQ. 'neoclassical') ) THEN
1293 
1294  found_neoclassical_current: IF (ASSOCIATED(coresource_iter%SOURCE(isrc)%profiles_1d(1)%j_parallel)) THEN
1295  fun2(:) = coresource_iter%SOURCE(isrc)%profiles_1d(1)%j_parallel(:)
1296  CALL l3interp(fun2, coresource_iter%SOURCE(isrc)%profiles_1d(1)%grid%rho_tor_norm, nsrc, dfun2, geometry%RHO_NORM, nprof)
1297  profiles%JBOOT = dfun2
1298  END IF found_neoclassical_current
1299  END IF check_neoclassical_sources
1300 
1301  END IF check_source_has_identifier_neoclassical
1302 
1303 
1304  IF (ALLOCATED(fun2)) DEALLOCATE (fun2)
1305  IF (ALLOCATED(dfun2)) DEALLOCATE (dfun2)
1306 
1307  END DO loop_over_sources
1308 
1309 
1310  IF ( maxval(abs(sources%CURR_EXP + sources%CURR_IMP*profiles%PSI)) > 1.0e-20_ids_real ) THEN
1311  profiles%JNI = sources%CURR_EXP + sources%CURR_IMP*profiles%PSI
1312  END IF
1313 
1314 
1315 666 RETURN
1316 
1317 
1318  END SUBROUTINE convert_ids_to_internal_types
1319 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
1320 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
1321 
1322 
1323 
1324 
1325 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
1326 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
1327 !-------------------------------------------------------
1337 !-------------------------------------------------------
1339  (geometry, profiles, transport, sources, impurity, global, numerics, &
1340  composition, &
1341  coreprofiles_new, coretransport_new, coresource_new, coresolver_new, diag)
1342 
1343 
1344 
1345  USE imas_constants_module
1346  USE ids_schemas
1347  USE ids_routines
1348  USE imas_core_source_identifier, ONLY: core_source_identifier
1349  USE imas_core_transport_identifier, ONLY: core_transport_identifier
1350  USE l3interps
1351 
1352  USE ets_plasma
1354  USE get_composition
1357 
1358  IMPLICIT NONE
1359 
1360 
1361 ! +++ ETS derived types:
1362  TYPE (magnetic_geometry) :: geometry !contains all geometry quantities
1363  TYPE (plasma_profiles) :: profiles !contains profiles of plasma parameters
1364  TYPE (transport_coefficients) :: transport !contains profiles of trasport coefficients
1365  TYPE (sources_and_sinks) :: sources !contains profiles of sources
1366  TYPE (impurity_profiles) :: impurity !contains profiles of impurities
1367  TYPE (global_param) :: global !contains global plasma parameters
1368  TYPE (transport_solver_numerics) :: numerics !contains all parameters required by run
1369  TYPE (plasma_composition) :: composition
1370  TYPE (diagnostic) :: diag !contains warnings and error messages
1371 
1372 
1373 ! +++ IDS derived types:
1374  TYPE (ids_core_profiles) :: coreprofiles_new
1375  TYPE (ids_core_transport) :: coretransport_new
1376  TYPE (ids_core_sources) :: coresource_new
1377  TYPE (ids_transport_solver_numerics) :: coresolver_new
1378 
1379 
1380 ! +++ Dimensions & indexes:
1381  INTEGER (IDS_INT) :: nprof
1382  INTEGER (IDS_INT) :: ntra
1383  INTEGER (IDS_INT) :: nsrc
1384  INTEGER (IDS_INT) :: nsol
1385 
1386  INTEGER (IDS_INT) :: ntime, nmod, nsrs, neq, ncoef, nbc
1387  INTEGER (IDS_INT) :: iion,nion, iimp,nimp, inuc, nnuc, irho
1388  INTEGER (IDS_INT) :: izimp,nzimp, iel, nel, ind, isol, ieq, idim, ival
1389  INTEGER (IDS_INT) :: ibnd, ibnd0
1390 
1391  REAL (IDS_REAL) :: nztot
1392  REAL (IDS_REAL) :: arr(geometry%nrho)
1393  REAL (IDS_REAL) :: ni_sum(geometry%nrho)
1394  REAL (IDS_REAL) :: nith_sum(geometry%nrho)
1395  REAL (IDS_REAL) :: nith_ti_sum(geometry%nrho)
1396 
1397 
1398  ntime = 1 !Single slice output
1399 
1400 
1401 ! OUTPUT COREPROFILES:
1402  CALL allocate_coreprofiles_ids(ntime, profiles%NRHO, composition, coreprofiles_new)
1403 
1404  coreprofiles_new%ids_properties%comment(1) = "Output from a transport solver"
1405  coreprofiles_new%ids_properties%source(1) = "transport_solver"
1406  coreprofiles_new%ids_properties%homogeneous_time = 1
1407 
1408  coreprofiles_new%vacuum_toroidal_field%r0 = geometry%R0
1409  coreprofiles_new%vacuum_toroidal_field%b0(1) = geometry%B0
1410 
1411  coreprofiles_new%profiles_1d(1)%grid%rho_tor_norm = geometry%RHO_NORM
1412  coreprofiles_new%profiles_1d(1)%grid%rho_tor = geometry%RHO
1413  coreprofiles_new%profiles_1d(1)%grid%psi = profiles%PSI
1414  coreprofiles_new%profiles_1d(1)%grid%volume = geometry%VOL
1415  coreprofiles_new%profiles_1d(1)%grid%area = geometry%AREA
1416 
1417  ! Initialise, the value is accumulated below
1418  coreprofiles_new%profiles_1d(1)%pressure_thermal = 0.0_ids_real
1419  coreprofiles_new%profiles_1d(1)%pressure_ion_total = 0.0_ids_real
1420  coreprofiles_new%profiles_1d(1)%pressure_parallel = 0.0_ids_real
1421  coreprofiles_new%profiles_1d(1)%pressure_perpendicular = 0.0_ids_real
1422 
1423  coreprofiles_new%profiles_1d(1)%electrons%temperature = profiles%TE
1424  coreprofiles_new%profiles_1d(1)%electrons%density_thermal = profiles%NE
1425  coreprofiles_new%profiles_1d(1)%electrons%density_fast = profiles%NE_FAST
1426  coreprofiles_new%profiles_1d(1)%electrons%density = coreprofiles_new%profiles_1d(1)%electrons%density_thermal + &
1427  coreprofiles_new%profiles_1d(1)%electrons%density_fast
1428 
1429  coreprofiles_new%profiles_1d(1)%electrons%pressure_thermal = profiles%PE
1430  coreprofiles_new%profiles_1d(1)%electrons%pressure_fast_perpendicular = profiles%PE_FAST
1431  coreprofiles_new%profiles_1d(1)%electrons%pressure_fast_parallel = profiles%PE_FAST_PARALLEL
1432  coreprofiles_new%profiles_1d(1)%electrons%pressure = coreprofiles_new%profiles_1d(1)%electrons%pressure_thermal + &
1433  (1.0_ids_real/3.0_ids_real) * coreprofiles_new%profiles_1d(1)%electrons%pressure_fast_parallel + &
1434  (2.0_ids_real/3.0_ids_real) * coreprofiles_new%profiles_1d(1)%electrons%pressure_fast_perpendicular
1435 
1436  ! Add electron data to sums over electrons, ions, impurities
1437  coreprofiles_new%profiles_1d(1)%pressure_thermal = coreprofiles_new%profiles_1d(1)%pressure_thermal + &
1438  coreprofiles_new%profiles_1d(1)%electrons%pressure_thermal
1439  coreprofiles_new%profiles_1d(1)%pressure_parallel = coreprofiles_new%profiles_1d(1)%pressure_parallel + &
1440  coreprofiles_new%profiles_1d(1)%electrons%pressure_thermal + &
1441  coreprofiles_new%profiles_1d(1)%electrons%pressure_fast_parallel
1442  coreprofiles_new%profiles_1d(1)%pressure_perpendicular = coreprofiles_new%profiles_1d(1)%pressure_perpendicular + &
1443  coreprofiles_new%profiles_1d(1)%electrons%pressure_thermal + &
1444  coreprofiles_new%profiles_1d(1)%electrons%pressure_fast_perpendicular
1445 
1446  ions: DO iion = 1, composition%NION
1447  nel = composition%IONS(iion)%ELINDEX
1448  nnuc = SIZE(composition%ELEMENTS(nel)%NUCINDEX)
1449  DO inuc = 1, nnuc
1450  ind = composition%ELEMENTS(nel)%NUCINDEX(inuc)
1451  coreprofiles_new%profiles_1d(1)%ion(iion)%element(inuc)%a = composition%NUCLEI(ind)%AMN
1452  coreprofiles_new%profiles_1d(1)%ion(iion)%element(inuc)%z_n = composition%NUCLEI(ind)%ZN
1453  coreprofiles_new%profiles_1d(1)%ion(iion)%element(inuc)%atoms_n = composition%ELEMENTS(nel)%NUCMULT(inuc)
1454  END DO
1455  coreprofiles_new%profiles_1d(1)%ion(iion)%z_ion = composition%IONS(iion)%ZION
1456  coreprofiles_new%profiles_1d(1)%ion(iion)%z_ion_1d = composition%IONS(iion)%ZION
1457  coreprofiles_new%profiles_1d(1)%ion(iion)%z_ion_square_1d = composition%IONS(iion)%ZION**2
1458  coreprofiles_new%profiles_1d(1)%ion(iion)%neutral_index = composition%IONS(iion)%NEUTINDEX
1459 ! COREPROFILES_NEW%profiles_1d(1)%neutral(COMPOSITION%IONS(IION)%NEUTINDEX)%ion_index = IION
1460 
1461  coreprofiles_new%profiles_1d(1)%ion(iion)%temperature(:) = profiles%TI(:,iion)
1462  coreprofiles_new%profiles_1d(1)%ion(iion)%density_fast(:) = profiles%NI_FAST(:,iion)
1463  coreprofiles_new%profiles_1d(1)%ion(iion)%density_thermal(:) = profiles%NI(:,iion)
1464  coreprofiles_new%profiles_1d(1)%ion(iion)%density(:) = coreprofiles_new%profiles_1d(1)%ion(iion)%density_thermal(:) + &
1465  coreprofiles_new%profiles_1d(1)%ion(iion)%density_fast(:)
1466 
1467  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_thermal(:) = profiles%PI(:,iion)
1468  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_fast_perpendicular(:)= profiles%PI_FAST(:,iion)
1469  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_fast_parallel(:) = profiles%PI_FAST_PARALLEL(:,iion)
1470  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure(:) = coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_thermal + &
1471  (1.0_ids_real/3.0_ids_real) * coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_fast_parallel + &
1472  (2.0_ids_real/3.0_ids_real) * coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_fast_perpendicular
1473  !COREPROFILES_NEW%profiles_1d(1)%ion(IION)%velocity%poloidal = 0.0_IDS_REAL
1474  !COREPROFILES_NEW%profiles_1d(1)%ion(IION)%velocity%toroidal = 0.0_IDS_REAL
1475  deallocate(&
1476  coreprofiles_new%profiles_1d(1)%ion(iion)%velocity%poloidal, &
1477  coreprofiles_new%profiles_1d(1)%ion(iion)%velocity%toroidal)
1478 
1479  ! Add ion data to sums over electrons, ions and impurities, or sums over all ions and impurities
1480  coreprofiles_new%profiles_1d(1)%pressure_ion_total(:) = coreprofiles_new%profiles_1d(1)%pressure_ion_total(:) + &
1481  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure(:)
1482  coreprofiles_new%profiles_1d(1)%pressure_thermal(:) = coreprofiles_new%profiles_1d(1)%pressure_thermal(:) + &
1483  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_thermal(:)
1484  coreprofiles_new%profiles_1d(1)%pressure_parallel = coreprofiles_new%profiles_1d(1)%pressure_parallel + &
1485  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_thermal + &
1486  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_fast_parallel
1487  coreprofiles_new%profiles_1d(1)%pressure_perpendicular = coreprofiles_new%profiles_1d(1)%pressure_perpendicular + &
1488  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_thermal + &
1489  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_fast_perpendicular
1490  END DO ions
1491 
1492  impurities: DO iimp = 1, composition%NIMP
1493  nel = composition%IMPURITY(iimp)%ELINDEX
1494  nnuc = SIZE(composition%ELEMENTS(nel)%NUCINDEX)
1495 
1496  iion = composition%NION + iimp
1497 
1498  DO inuc = 1, nnuc
1499  ind = composition%ELEMENTS(nel)%NUCINDEX(inuc)
1500  coreprofiles_new%profiles_1d(1)%ion(iion)%element(inuc)%a = composition%NUCLEI(ind)%AMN
1501  coreprofiles_new%profiles_1d(1)%ion(iion)%element(inuc)%z_n = composition%NUCLEI(ind)%ZN
1502  coreprofiles_new%profiles_1d(1)%ion(iion)%element(inuc)%atoms_n = composition%ELEMENTS(nel)%NUCMULT(inuc)
1503  END DO
1504 
1505  coreprofiles_new%profiles_1d(1)%ion(iion)%z_ion_1d = 0.0_ids_real
1506  coreprofiles_new%profiles_1d(1)%ion(iion)%z_ion_square_1d = 0.0_ids_real
1507  coreprofiles_new%profiles_1d(1)%ion(iion)%neutral_index = composition%IMPURITY(iimp)%NEUTINDEX
1508 ! COREPROFILES_NEW%profiles_1d(1)%neutral(COMPOSITION%IMPURITY(IIMP)%NEUTINDEX)%ion_index = IION
1509  coreprofiles_new%profiles_1d(1)%ion(iion)%density = 0.0_ids_real
1510  coreprofiles_new%profiles_1d(1)%ion(iion)%density_fast = 0.0_ids_real
1511  coreprofiles_new%profiles_1d(1)%ion(iion)%density_thermal = 0.0_ids_real
1512  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure = 0.0_ids_real
1513  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_thermal = 0.0_ids_real
1514  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_fast_perpendicular = 0.0_ids_real
1515  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_fast_parallel = 0.0_ids_real
1516  coreprofiles_new%profiles_1d(1)%ion(iion)%velocity%poloidal = 0.0_ids_real
1517  coreprofiles_new%profiles_1d(1)%ion(iion)%velocity%toroidal = 0.0_ids_real
1518  nztot = 0.0_ids_real
1519 
1520  nzimp = SIZE(composition%IMPURITY(iimp)%IDS_state)
1521  states: DO izimp = 1, nzimp
1522  nztot = nztot + impurity%NZTOT(iimp,izimp)
1523  coreprofiles_new%profiles_1d(1)%ion(iion)%z_ion_1d(:) = coreprofiles_new%profiles_1d(1)%ion(iion)%z_ion_1d(:)&
1524  + impurity%Z_AVE(:,iimp,izimp) * impurity%NZ(:,iimp,izimp)
1525  coreprofiles_new%profiles_1d(1)%ion(iion)%z_ion_square_1d(:) = coreprofiles_new%profiles_1d(1)%ion(iion)%z_ion_square_1d(:)&
1526  + impurity%Z2_AVE(:,iimp,izimp) * impurity%NZ(:,iimp,izimp)
1527 
1528  coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%z_min = impurity%ZMIN(iimp,izimp)
1529  coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%z_max = impurity%ZMAX(iimp,izimp)
1530  coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%z_average_1d(:) = impurity%Z_AVE(:,iimp,izimp)
1531  coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%z_average_square_1d(:) = impurity%Z2_AVE(:,iimp,izimp)
1532  coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%temperature(:) = impurity%TZ(:,iimp,izimp)
1533  coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%density_thermal(:) = impurity%NZ(:,iimp,izimp)
1534  !COREPROFILES_NEW%profiles_1d(1)%ion(IION)%state(IZIMP)%density_fast = 0.0_IDS_REAL
1535  deallocate(coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%density_fast)
1536  coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%density(:) = coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%density_thermal(:) + &
1537  coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%density_fast(:)
1538 
1539  coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%pressure_thermal(:) = impurity%TZ(:,iimp,izimp)*impurity%NZ(:,iimp,izimp) * imas_constants%ev
1540  !COREPROFILES_NEW%profiles_1d(1)%ion(IION)%state(IZIMP)%pressure_fast_perpendicular(:) = 0.0_IDS_REAL ! Internal field is missing
1541  !COREPROFILES_NEW%profiles_1d(1)%ion(IION)%state(IZIMP)%pressure_fast_parallel(:) = 0.0_IDS_REAL ! Internal field is missing
1542  deallocate(&
1543  coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%pressure_fast_perpendicular, &
1544  coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%pressure_fast_parallel)
1545 
1546  coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%pressure(:) = coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%pressure_thermal(:) + &
1547  (1.0_ids_real/3.0_ids_real) * coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%pressure_fast_parallel(:) + &
1548  (2.0_ids_real/3.0_ids_real) * coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%pressure_fast_perpendicular(:)
1549  !COREPROFILES_NEW%profiles_1d(1)%ion(IION)%state(IZIMP)%velocity%poloidal = 0.0_IDS_REAL
1550  !COREPROFILES_NEW%profiles_1d(1)%ion(IION)%state(IZIMP)%velocity%toroidal = 0.0_IDS_REAL
1551  deallocate(&
1552  coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%velocity%poloidal, &
1553  coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%velocity%toroidal)
1554 
1555  ! Density, state -> ion
1556  coreprofiles_new%profiles_1d(1)%ion(iion)%density_thermal(:) = coreprofiles_new%profiles_1d(1)%ion(iion)%density_thermal(:) + &
1557  coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%density_thermal(:)
1558  coreprofiles_new%profiles_1d(1)%ion(iion)%density_fast(:) = coreprofiles_new%profiles_1d(1)%ion(iion)%density_fast(:) + &
1559  coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%density_fast(:)
1560  coreprofiles_new%profiles_1d(1)%ion(iion)%density(:) = coreprofiles_new%profiles_1d(1)%ion(iion)%density(:) + &
1561  coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%density(:)
1562 
1563  ! Pressure, state -> ion
1564  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure(:) = coreprofiles_new%profiles_1d(1)%ion(iion)%pressure(:) + &
1565  coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%pressure(:)
1566  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_thermal(:) = coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_thermal(:) + &
1567  coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%pressure_thermal(:)
1568  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_fast_perpendicular(:) = coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_fast_perpendicular(:) + &
1569  coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%pressure_fast_perpendicular(:)
1570  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_fast_parallel(:) = coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_fast_parallel(:) + &
1571  coreprofiles_new%profiles_1d(1)%ion(iion)%state(izimp)%pressure_fast_parallel(:)
1572 
1573  END DO states
1574 
1575  IF ( nztot > sqrt(epsilon(1.0e0_ids_real)) ) THEN
1576  coreprofiles_new%profiles_1d(1)%ion(iion)%z_ion = coreprofiles_new%profiles_1d(1)%ion(iion)%z_ion / nztot
1577  END IF
1578 
1579  DO irho = 1 , profiles%NRHO
1580  IF ( coreprofiles_new%profiles_1d(1)%ion(iion)%density(irho) > sqrt(epsilon(1.0e0_ids_real)) ) THEN
1581  coreprofiles_new%profiles_1d(1)%ion(iion)%z_ion_1d(irho) = coreprofiles_new%profiles_1d(1)%ion(iion)%z_ion_1d(irho) &
1582  / coreprofiles_new%profiles_1d(1)%ion(iion)%density(irho)
1583  coreprofiles_new%profiles_1d(1)%ion(iion)%z_ion_square_1d(irho) = coreprofiles_new%profiles_1d(1)%ion(iion)%z_ion_square_1d(irho) &
1584  / coreprofiles_new%profiles_1d(1)%ion(iion)%density(irho)
1585  coreprofiles_new%profiles_1d(1)%ion(iion)%temperature(irho) = coreprofiles_new%profiles_1d(1)%ion(iion)%pressure(irho) &
1586  / coreprofiles_new%profiles_1d(1)%ion(iion)%density(irho) &
1587  / imas_constants%ev
1588  END IF
1589  END DO
1590 
1591  ! Sum pressures over all impurities
1592  coreprofiles_new%profiles_1d(1)%pressure_ion_total(:) = coreprofiles_new%profiles_1d(1)%pressure_ion_total(:) + &
1593  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure(:)
1594  coreprofiles_new%profiles_1d(1)%pressure_thermal(:) = coreprofiles_new%profiles_1d(1)%pressure_thermal(:) + &
1595  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_thermal(:)
1596  coreprofiles_new%profiles_1d(1)%pressure_parallel(:) = coreprofiles_new%profiles_1d(1)%pressure_parallel(:) + &
1597  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_thermal(:) + &
1598  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_fast_parallel(:)
1599  coreprofiles_new%profiles_1d(1)%pressure_perpendicular(:) = coreprofiles_new%profiles_1d(1)%pressure_parallel(:) + &
1600  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_thermal(:) + &
1601  coreprofiles_new%profiles_1d(1)%ion(iion)%pressure_fast_perpendicular(:)
1602  END DO impurities
1603 
1604  ni_sum = 0.0_ids_real
1605  nith_sum = 0.0_ids_real
1606  nith_ti_sum = 0.0_ids_real
1607  DO iion=1,profiles%NION
1608  ni_sum = ni_sum + coreprofiles_new%profiles_1d(1)%ion(iion)%density
1609  nith_sum = nith_sum + coreprofiles_new%profiles_1d(1)%ion(iion)%density_thermal
1610  nith_ti_sum = nith_ti_sum + coreprofiles_new%profiles_1d(1)%ion(iion)%density_thermal * coreprofiles_new%profiles_1d(1)%ion(iion)%temperature
1611  END DO
1612 
1613  coreprofiles_new%profiles_1d(1)%t_i_average = nith_ti_sum / nith_sum
1614  coreprofiles_new%profiles_1d(1)%n_i_total_over_n_e = ni_sum / coreprofiles_new%profiles_1d(1)%electrons%density
1615 
1616  ! Not sure how to create momentum_tor
1617  ! Should be: COREPROFILES_NEW%profiles_1d(1)%ion(IION)%velocity%toroidal * <R^2> * mass
1618  ! But how to get <R^@>?.
1619  ! OLD: COREPROFILES_NEW%profiles_1d(1)%momentum_tor = 0.0_IDS_REAL
1620  ! Temporary solution is deallocate:
1621  deallocate(coreprofiles_new%profiles_1d(1)%momentum_tor)
1622 
1623  coreprofiles_new%profiles_1d(1)%zeff = profiles%ZEFF
1624 ! Filled already COREPROFILES_NEW%profiles_1d(1)%pressure_ion_total = 0.0_IDS_REAL
1625 ! Filled already COREPROFILES_NEW%profiles_1d(1)%pressure_thermal = 0.0_IDS_REAL
1626 ! Filled already COREPROFILES_NEW%profiles_1d(1)%pressure_perpendicular = 0.0_IDS_REAL
1627 ! Filled already COREPROFILES_NEW%profiles_1d(1)%pressure_parallel = 0.0_IDS_REAL
1628  coreprofiles_new%profiles_1d(1)%j_total = profiles%CURR_PAR
1629  coreprofiles_new%profiles_1d(1)%j_tor = profiles%CURR_TOR
1630  coreprofiles_new%profiles_1d(1)%j_ohmic = profiles%CURR_PAR - profiles%JNI
1631  coreprofiles_new%profiles_1d(1)%j_non_inductive = profiles%JNI
1632  coreprofiles_new%profiles_1d(1)%j_bootstrap = profiles%JBOOT
1633  coreprofiles_new%profiles_1d(1)%conductivity_parallel = transport%SIGMA
1634  coreprofiles_new%profiles_1d(1)%e_field_parallel = profiles%E_PAR
1635 ! COREPROFILES_NEW%profiles_1d(1)%e_field%radial = 0.0_IDS_REAL
1636 ! COREPROFILES_NEW%profiles_1d(1)%e_field%diamagnetic = 0.0_IDS_REAL
1637  coreprofiles_new%profiles_1d(1)%e_field%parallel = profiles%E_PAR
1638 ! COREPROFILES_NEW%profiles_1d(1)%e_field%poloidal = 0.0_IDS_REAL
1639 ! COREPROFILES_NEW%profiles_1d(1)%e_field%toroidal = 0.0_IDS_REAL
1640  deallocate(&
1641  coreprofiles_new%profiles_1d(1)%e_field%radial, &
1642  coreprofiles_new%profiles_1d(1)%e_field%diamagnetic, &
1643  coreprofiles_new%profiles_1d(1)%e_field%poloidal, &
1644  coreprofiles_new%profiles_1d(1)%e_field%toroidal)
1645  coreprofiles_new%profiles_1d(1)%q = profiles%QSF
1646  coreprofiles_new%profiles_1d(1)%magnetic_shear = profiles%SHEAR
1647 
1648 
1649  coreprofiles_new%global_quantities%ip(ntime) = global%CURRENT
1650  coreprofiles_new%global_quantities%current_non_inductive(ntime) = global%CURRENT_NI
1651  !COREPROFILES_NEW%global_quantities%current_bootstrap(NTIME) = 0.0_IDS_REAL
1652  deallocate(coreprofiles_new%global_quantities%current_bootstrap)
1653  !COREPROFILES_NEW%global_quantities%v_loop(NTIME) = GLOBAL%VLOOP
1654  deallocate(coreprofiles_new%global_quantities%v_loop)
1655  !COREPROFILES_NEW%global_quantities%li(NTIME) = GLOBAL%LI
1656  deallocate(coreprofiles_new%global_quantities%li)
1657  !COREPROFILES_NEW%global_quantities%li_3(NTIME) = GLOBAL%LI
1658  deallocate(coreprofiles_new%global_quantities%li_3)
1659  !COREPROFILES_NEW%global_quantities%beta_tor(NTIME) = GLOBAL%BETA_TOR
1660  deallocate(coreprofiles_new%global_quantities%beta_tor)
1661  !COREPROFILES_NEW%global_quantities%beta_tor_norm(NTIME) = GLOBAL%BETA_N
1662  deallocate(coreprofiles_new%global_quantities%beta_tor_norm)
1663  !COREPROFILES_NEW%global_quantities%beta_pol(NTIME) = GLOBAL%BETA_POL
1664  deallocate(coreprofiles_new%global_quantities%beta_pol)
1665  !COREPROFILES_NEW%global_quantities%energy_diamagnetic(NTIME) = GLOBAL%WDIA
1666  deallocate(coreprofiles_new%global_quantities%energy_diamagnetic)
1667 
1668  coreprofiles_new%time = 0.0_ids_real
1669 
1670 
1671 
1672 ! OUTPUT CORETRANSPORT:
1673  nmod = 1
1674 
1675  CALL allocate_coretransport_ids(ntime, nmod, profiles%NRHO, composition, coretransport_new)
1676 
1677  coretransport_new%ids_properties%comment(1) = "Output from a transport solver"
1678  coretransport_new%ids_properties%source(1) = "transport_solver"
1679  coretransport_new%ids_properties%homogeneous_time = 1
1680 
1681  coretransport_new%vacuum_toroidal_field%r0 = geometry%R0
1682  coretransport_new%vacuum_toroidal_field%b0(1) = geometry%B0
1683 
1684  coretransport_new%model(1)%identifier%index = core_transport_identifier%transport_solver
1685  coretransport_new%model(1)%identifier%name(1) = core_transport_identifier%name( core_transport_identifier%transport_solver )
1686  coretransport_new%model(1)%identifier%description(1) = core_transport_identifier%description( core_transport_identifier%transport_solver )
1687 
1688  coretransport_new%model(1)%flux_multiplier = 1.5_ids_real
1689 
1690 
1691  coretransport_new%model(1)%profiles_1d(1)%grid_d%rho_tor_norm = geometry%RHO_NORM
1692  coretransport_new%model(1)%profiles_1d(1)%grid_d%rho_tor = geometry%RHO
1693  coretransport_new%model(1)%profiles_1d(1)%grid_d%psi = profiles%PSI
1694  coretransport_new%model(1)%profiles_1d(1)%grid_d%volume = geometry%VOL
1695  coretransport_new%model(1)%profiles_1d(1)%grid_d%area = geometry%AREA
1696 
1697  coretransport_new%model(1)%profiles_1d(1)%grid_v%rho_tor_norm = geometry%RHO_NORM
1698  coretransport_new%model(1)%profiles_1d(1)%grid_v%rho_tor = geometry%RHO
1699  coretransport_new%model(1)%profiles_1d(1)%grid_v%psi = profiles%PSI
1700  coretransport_new%model(1)%profiles_1d(1)%grid_v%volume = geometry%VOL
1701  coretransport_new%model(1)%profiles_1d(1)%grid_v%area = geometry%AREA
1702 
1703  coretransport_new%model(1)%profiles_1d(1)%grid_flux%rho_tor_norm = geometry%RHO_NORM
1704  coretransport_new%model(1)%profiles_1d(1)%grid_flux%rho_tor = geometry%RHO
1705  coretransport_new%model(1)%profiles_1d(1)%grid_flux%psi = profiles%PSI
1706  coretransport_new%model(1)%profiles_1d(1)%grid_flux%volume = geometry%VOL
1707  coretransport_new%model(1)%profiles_1d(1)%grid_flux%area = geometry%AREA
1708 
1709  coretransport_new%model(1)%profiles_1d(1)%conductivity_parallel = transport%SIGMA
1710 
1711  coretransport_new%model(1)%profiles_1d(1)%electrons%particles%d(:) = transport%DIFF_NE(:)
1712  coretransport_new%model(1)%profiles_1d(1)%electrons%particles%v(:) = transport%VCONV_NE(:)
1713  coretransport_new%model(1)%profiles_1d(1)%electrons%particles%flux(:) = profiles%FLUX_NE(:) &
1714  / max(geometry%SURFACE,geometry%SURFACE(2)/3.0_ids_real) ! Flux per square meters
1715 
1716  coretransport_new%model(1)%profiles_1d(1)%electrons%energy%d(:) = transport%DIFF_TE(:)
1717  coretransport_new%model(1)%profiles_1d(1)%electrons%energy%v(:) = transport%VCONV_TE(:)
1718 
1719  ! ARR: The heat flux convected by particle transport from model 3
1720  arr(:) = 5.0_ids_real/6.0_ids_real * profiles%TE(:) & !Correction to 3/2
1721  * geometry%VPR(:)*geometry%G1(:) &
1722  *(profiles%NE(:) *transport%VCONV_NE_MODEL(:,3) &
1723  - profiles%DNE(:)*transport%DIFF_NE_MODEL(:,3))
1724 
1725  coretransport_new%model(1)%profiles_1d(1)%electrons%energy%flux(:) = (profiles%FLUX_TE(:)+arr(:)) * imas_constants%ev &
1726  / max(geometry%SURFACE,geometry%SURFACE(2)/3.0_ids_real) ! Flux per square meters
1727 
1728  ions_tr: DO iion = 1, composition%NION
1729  nel = composition%IONS(iion)%ELINDEX
1730  ions_tr_has_elements: IF (nel > 0) THEN
1731  nnuc = SIZE(composition%ELEMENTS(nel)%NUCINDEX)
1732  DO inuc = 1, nnuc
1733  ind = composition%ELEMENTS(nel)%NUCINDEX(inuc)
1734  coretransport_new%model(1)%profiles_1d(1)%ion(iion)%element(inuc)%a = composition%NUCLEI(ind)%AMN
1735  coretransport_new%model(1)%profiles_1d(1)%ion(iion)%element(inuc)%z_n = composition%NUCLEI(ind)%ZN
1736  coretransport_new%model(1)%profiles_1d(1)%ion(iion)%element(inuc)%atoms_n = composition%ELEMENTS(nel)%NUCMULT(inuc)
1737  END DO
1738  coretransport_new%model(1)%profiles_1d(1)%ion(iion)%z_ion = composition%IONS(iion)%ZION
1739 
1740  coretransport_new%model(1)%profiles_1d(1)%ion(iion)%particles%d(:) = transport%DIFF_NI(:,iion)
1741  coretransport_new%model(1)%profiles_1d(1)%ion(iion)%particles%v(:) = transport%VCONV_NI(:,iion)
1742  coretransport_new%model(1)%profiles_1d(1)%ion(iion)%particles%flux(:) = profiles%FLUX_NI(:,iion) &
1743  / max(geometry%SURFACE,geometry%SURFACE(2)/3.0_ids_real) ! Flux per square meters
1744 
1745  coretransport_new%model(1)%profiles_1d(1)%ion(iion)%energy%d(:) = transport%DIFF_TI(:,iion)
1746  coretransport_new%model(1)%profiles_1d(1)%ion(iion)%energy%v(:) = transport%VCONV_TI(:,iion)
1747 
1748  ! ARR: The heat flux convected by particle transport from model 3
1749  arr(:) = 5.0_ids_real/6.0_ids_real * profiles%TI(:,iion) & !Correction to 3/2
1750  * geometry%VPR(:)*geometry%G1(:) &
1751  *(profiles%NI(:,iion) *transport%VCONV_NI_MODEL(:,iion,3) &
1752  - profiles%DNI(:,iion)*transport%DIFF_NI_MODEL(:,iion,3))
1753 
1754  coretransport_new%model(1)%profiles_1d(1)%ion(iion)%energy%flux(:) = (profiles%FLUX_TI(:,iion) + arr(:)) * imas_constants%ev &
1755  / max(geometry%SURFACE,geometry%SURFACE(2)/3.0_ids_real)
1756 
1757  coretransport_new%model(1)%profiles_1d(1)%ion(iion)%momentum%toroidal%d(:) = transport%DIFF_VTOR(:,iion)
1758  coretransport_new%model(1)%profiles_1d(1)%ion(iion)%momentum%toroidal%v(:) = transport%VCONV_VTOR(:,iion)
1759  coretransport_new%model(1)%profiles_1d(1)%ion(iion)%momentum%toroidal%flux(:) = 0.0_ids_real ! PROFILES%FLUX_MTOR(:,IION) <- This has to be converted to the correct units
1760  END IF ions_tr_has_elements
1761  END DO ions_tr
1762 
1763 
1764 
1765 ! OUTPUT CORESOURCES:
1766 ! NSRC = 5 ! The ETS outputs: (1) transport_solver, (2) collisional_equipartition, (3) ohmic, (4) ionization, (5) srecombination
1767  nsrc = 3 ! The ETS outputs: (1) transport_solver, (2) collisional_equipartition, (3) ohmic
1768 
1769  CALL allocate_coresources_ids(ntime, nsrc, profiles%NRHO, composition, coresource_new)
1770 
1771  coresource_new%ids_properties%comment(1) = "Output from a transport solver"
1772  coresource_new%ids_properties%source(1) = "transport_solver"
1773  coresource_new%ids_properties%homogeneous_time = 1
1774 
1775  coresource_new%vacuum_toroidal_field%r0 = geometry%R0
1776  coresource_new%vacuum_toroidal_field%b0(1) = geometry%B0
1777 
1778  coresource_new%source(1)%profiles_1d(1)%grid%rho_tor_norm = geometry%RHO_NORM
1779  coresource_new%source(1)%profiles_1d(1)%grid%rho_tor = geometry%RHO
1780  coresource_new%source(1)%profiles_1d(1)%grid%psi = profiles%PSI
1781  coresource_new%source(1)%profiles_1d(1)%grid%volume = geometry%VOL
1782  coresource_new%source(1)%profiles_1d(1)%grid%area = geometry%AREA
1783 
1784 
1785  ! The "transport_solver" identifier exclusive to the TRANSPORT_SOLVER actor.
1786  ! It is not specified by the CORE_SOURCES_IDENTIFIERS and is therefore given the negative index -1.
1787  coresource_new%source(1)%identifier%name(1) = "transport_solver"
1788  coresource_new%source(1)%identifier%index = -1
1789  coresource_new%source(1)%identifier%description(1) = "Data used in the transport solver"
1790 
1791  coresource_new%source(1)%profiles_1d(1)%j_parallel = profiles%CURR_PAR
1792 
1793  coresource_new%source(1)%profiles_1d(1)%electrons%particles = profiles%SOURCE_NE
1794  coresource_new%source(1)%profiles_1d(1)%electrons%particles_decomposed%explicit_part = profiles%SOURCE_NE_EXP
1795  coresource_new%source(1)%profiles_1d(1)%electrons%particles_decomposed%implicit_part = profiles%SOURCE_NE_IMP
1796  coresource_new%source(1)%profiles_1d(1)%electrons%particles_inside = profiles%INT_SOURCE_NE
1797  coresource_new%source(1)%profiles_1d(1)%electrons%energy = profiles%SOURCE_TE * imas_constants%ev
1798  coresource_new%source(1)%profiles_1d(1)%electrons%energy_decomposed%explicit_part = profiles%SOURCE_TE_EXP * imas_constants%ev
1799  coresource_new%source(1)%profiles_1d(1)%electrons%energy_decomposed%implicit_part = profiles%SOURCE_TE_IMP
1800  coresource_new%source(1)%profiles_1d(1)%electrons%power_inside = profiles%INT_SOURCE_TE
1801 
1802 
1803  ions_sr1: DO iion = 1, composition%NION
1804  nel = composition%IONS(iion)%ELINDEX
1805  nnuc = SIZE(composition%ELEMENTS(nel)%NUCINDEX)
1806  DO inuc = 1, nnuc
1807  ind = composition%ELEMENTS(nel)%NUCINDEX(inuc)
1808  coresource_new%source(1)%profiles_1d(1)%ion(iion)%element(inuc)%a = composition%NUCLEI(ind)%AMN
1809  coresource_new%source(1)%profiles_1d(1)%ion(iion)%element(inuc)%z_n = composition%NUCLEI(ind)%ZN
1810  coresource_new%source(1)%profiles_1d(1)%ion(iion)%element(inuc)%atoms_n = composition%ELEMENTS(nel)%NUCMULT(inuc)
1811  END DO
1812  coresource_new%source(1)%profiles_1d(1)%ion(iion)%z_ion = composition%IONS(iion)%ZION
1813 
1814  coresource_new%source(1)%profiles_1d(1)%ion(iion)%particles(:) = profiles%SOURCE_NI(:,iion)
1815  coresource_new%source(1)%profiles_1d(1)%ion(iion)%particles_decomposed%explicit_part(:) = profiles%SOURCE_NI_EXP(:,iion)
1816  coresource_new%source(1)%profiles_1d(1)%ion(iion)%particles_decomposed%implicit_part(:) = profiles%SOURCE_NI_IMP(:,iion)
1817 ! CORESOURCE_NEW%source(1)%profiles_1d(1)%ion(IION)%particles_inside(:) = PROFILES%INT_SOURCE_NI(:,IION)
1818  coresource_new%source(1)%profiles_1d(1)%ion(iion)%energy(:) = profiles%SOURCE_TI(:,iion) * imas_constants%ev
1819  coresource_new%source(1)%profiles_1d(1)%ion(iion)%energy_decomposed%explicit_part(:) = profiles%SOURCE_TI_EXP(:,iion) * imas_constants%ev
1820  coresource_new%source(1)%profiles_1d(1)%ion(iion)%energy_decomposed%implicit_part(:) = profiles%SOURCE_TI_IMP(:,iion)
1821 ! CORESOURCE_NEW%source(1)%profiles_1d(1)%ion(IION)%power_inside(:) = PROFILES%INT_SOURCE_TI(:,IION)
1822  coresource_new%source(1)%profiles_1d(1)%ion(iion)%momentum%toroidal(:) = profiles%SOURCE_MTOR(:,iion)
1823  coresource_new%source(1)%profiles_1d(1)%ion(iion)%momentum%toroidal_decomposed%explicit_part(:) = profiles%SOURCE_MTOR_EXP(:,iion)
1824  coresource_new%source(1)%profiles_1d(1)%ion(iion)%momentum%toroidal_decomposed%implicit_part(:) = profiles%SOURCE_MTOR_IMP(:,iion)
1825 ! CORESOURCE_NEW%source(1)%profiles_1d(1)%ion(IION)%torque_inside(:) = PROFILES%INT_SOURCE_MTOR(:,IION)
1826 
1827  END DO ions_sr1
1828 
1829 
1830  coresource_new%source(2)%identifier%index = core_source_identifier%collisional_equipartition
1831  coresource_new%source(2)%identifier%name(1) = core_source_identifier%name( core_source_identifier%collisional_equipartition )
1832  coresource_new%source(2)%identifier%description(1) = core_source_identifier%description( core_source_identifier%collisional_equipartition )
1833 
1834  coresource_new%source(2)%profiles_1d(1)%grid%rho_tor_norm = geometry%RHO_NORM
1835  coresource_new%source(2)%profiles_1d(1)%grid%rho_tor = geometry%RHO
1836  coresource_new%source(2)%profiles_1d(1)%grid%psi = profiles%PSI
1837  coresource_new%source(2)%profiles_1d(1)%grid%volume = geometry%VOL
1838  coresource_new%source(2)%profiles_1d(1)%grid%area = geometry%AREA
1839 
1840 
1841 
1842  coresource_new%source(3)%identifier%index = core_source_identifier%ohmic
1843  coresource_new%source(3)%identifier%name(1) = core_source_identifier%name( core_source_identifier%ohmic )
1844  coresource_new%source(3)%identifier%description(1) = core_source_identifier%description( core_source_identifier%ohmic )
1845 
1846  coresource_new%source(3)%profiles_1d(1)%grid%rho_tor_norm = geometry%RHO_NORM
1847  coresource_new%source(3)%profiles_1d(1)%grid%rho_tor = geometry%RHO
1848  coresource_new%source(3)%profiles_1d(1)%grid%psi = profiles%PSI
1849  coresource_new%source(3)%profiles_1d(1)%grid%volume = geometry%VOL
1850  coresource_new%source(3)%profiles_1d(1)%grid%area = geometry%AREA
1851 
1852 ! From Irena
1853 ! Remove since these sources should be provided in the input IDS (Thomas, Sept 2020)
1854 !
1855 !!$ CORESOURCE_NEW%source(4)%identifier%name(1) = "ionisation"
1856 !!$ CORESOURCE_NEW%source(4)%identifier%index = 601
1857 !!$ CORESOURCE_NEW%source(4)%identifier%description(1) = "Source from ionisation processes"
1858 !!$
1859 !!$ CORESOURCE_NEW%source(4)%profiles_1d(1)%grid%rho_tor_norm = GEOMETRY%RHO_NORM
1860 !!$ CORESOURCE_NEW%source(4)%profiles_1d(1)%grid%rho_tor = GEOMETRY%RHO
1861 !!$ CORESOURCE_NEW%source(4)%profiles_1d(1)%grid%psi = PROFILES%PSI
1862 !!$ CORESOURCE_NEW%source(4)%profiles_1d(1)%grid%volume = GEOMETRY%VOL
1863 !!$ CORESOURCE_NEW%source(4)%profiles_1d(1)%grid%area = GEOMETRY%AREA
1864 !!$
1865 !!$ CORESOURCE_NEW%source(5)%identifier%name(1) = "recombination"
1866 !!$ CORESOURCE_NEW%source(5)%identifier%index = 602
1867 !!$ CORESOURCE_NEW%source(5)%identifier%description(1) = "Source from recombination processes "
1868 !!$
1869 !!$ CORESOURCE_NEW%source(5)%profiles_1d(1)%grid%rho_tor_norm = GEOMETRY%RHO_NORM
1870 !!$ CORESOURCE_NEW%source(5)%profiles_1d(1)%grid%rho_tor = GEOMETRY%RHO
1871 !!$ CORESOURCE_NEW%source(5)%profiles_1d(1)%grid%psi = PROFILES%PSI
1872 !!$ CORESOURCE_NEW%source(5)%profiles_1d(1)%grid%volume = GEOMETRY%VOL
1873 !!$ CORESOURCE_NEW%source(5)%profiles_1d(1)%grid%area = GEOMETRY%AREA
1874 
1875 
1876 !Irena
1877 
1878 
1879 
1880 ! OUTPUT CORESOLVER:
1881 
1882  ncoef = 8
1883  nbc = 2
1884 
1885  neq = 0
1886  DO isol = 1,SIZE(numerics%SOLVER)
1887  neq = neq + numerics%SOLVER(isol)%NDIM
1888  END DO
1889 
1890  CALL allocate_coresolver_ids(ntime, neq, ncoef, profiles%NRHO, nbc, coresolver_new)
1891 
1892 
1893  coresolver_new%ids_properties%comment(1) = "Output from a transport solver"
1894  coresolver_new%ids_properties%source(1) = "transport_solver"
1895  coresolver_new%ids_properties%homogeneous_time = 1
1896  ALLOCATE(coresolver_new%time_step%data(1))
1897  coresolver_new%time_step%data(1) = numerics%TAU
1898 
1899  coresolver_new%vacuum_toroidal_field%r0 = geometry%R0
1900  coresolver_new%vacuum_toroidal_field%b0(1) = geometry%B0
1901 
1902  coresolver_new%solver%index = numerics%CONTROL%SOLVER_TYPE
1903 
1904  coresolver_new%SOLVER_1D(1)%grid%rho_tor_norm = geometry%RHO_NORM
1905  coresolver_new%SOLVER_1D(1)%grid%rho_tor = geometry%RHO
1906  coresolver_new%SOLVER_1D(1)%grid%psi = profiles%PSI
1907  coresolver_new%SOLVER_1D(1)%grid%volume = geometry%VOL
1908  coresolver_new%SOLVER_1D(1)%grid%area = geometry%AREA
1909 
1910 
1911  associate(cp => coresolver_new%solver_1d(1)%control_parameters)
1912  call init_control_parameters(cp)
1913 
1914  if (profiles%HYPERDIFFUSION%ACTIVE) then
1915  ival=1
1916  else
1917  ival=0
1918  end if
1919  call set_control_parameter(cp, 'hyperdiffusion_active' , ival, diag%ERROR_INDEX)
1920  if (diag%ERROR_INDEX < 0) then
1921  diag%ERROR_MESSAGE = 'in convert, error from set_control_parameter, failed to set hyperdiffusion_active'
1922  return
1923  end if
1924  call set_control_parameter(cp, 'hyperdiffusion_density_terms' , profiles%HYPERDIFFUSION%DENS%TERMS, diag%ERROR_INDEX)
1925  if (diag%ERROR_INDEX < 0) then
1926  diag%ERROR_MESSAGE = 'in convert, error from set_control_parameter, failed to set hyperdiffusion_density_terms'
1927  return
1928  end if
1929  call set_control_parameter(cp, 'hyperdiffusion_temperature_terms', profiles%HYPERDIFFUSION%TEMP%TERMS, diag%ERROR_INDEX)
1930  if (diag%ERROR_INDEX < 0) then
1931  diag%ERROR_MESSAGE = 'in convert, error from set_control_parameter, failed to set hyperdiffusion_temperature_terms'
1932  return
1933  end if
1934  call set_control_parameter(cp, 'hyperdiffusion_rotation_terms' , profiles%HYPERDIFFUSION%ROTATION%TERMS, diag%ERROR_INDEX)
1935  if (diag%ERROR_INDEX < 0) then
1936  diag%ERROR_MESSAGE = 'in convert, error from set_control_parameter, failed to set hyperdiffusion_rotation_terms'
1937  return
1938  end if
1939 
1940  call set_control_parameter(cp, 'hyperdiffusion_density_implicit' , profiles%HYPERDIFFUSION%DENS%IMPLICIT, diag%ERROR_INDEX)
1941  if (diag%ERROR_INDEX < 0) then
1942  diag%ERROR_MESSAGE = 'in convert, error from set_control_parameter, failed to set hyperdiffusion_density_implicit'
1943  return
1944  end if
1945  call set_control_parameter(cp, 'hyperdiffusion_density_explicit' , profiles%HYPERDIFFUSION%DENS%EXPLICIT, diag%ERROR_INDEX)
1946  if (diag%ERROR_INDEX < 0) then
1947  diag%ERROR_MESSAGE = 'in convert, error from set_control_parameter, failed to set hyperdiffusion_density_explicit'
1948  return
1949  end if
1950  call set_control_parameter(cp, 'hyperdiffusion_density_rho_tor_norm_cutoff' , profiles%HYPERDIFFUSION%DENS%RHO_NORM_CUTOFF, diag%ERROR_INDEX)
1951  if (diag%ERROR_INDEX < 0) then
1952  diag%ERROR_MESSAGE = 'in convert, error from set_control_parameter, failed to set hyperdiffusion_density_rho_tor_norm_cutoff'
1953  return
1954  end if
1955  call set_control_parameter(cp, 'hyperdiffusion_temperature_implicit' , profiles%HYPERDIFFUSION%TEMP%IMPLICIT, diag%ERROR_INDEX)
1956  if (diag%ERROR_INDEX < 0) then
1957  diag%ERROR_MESSAGE = 'in convert, error from set_control_parameter, failed to set hyperdiffusion_temperature_implicit'
1958  return
1959  end if
1960  call set_control_parameter(cp, 'hyperdiffusion_temperature_explicit' , profiles%HYPERDIFFUSION%TEMP%EXPLICIT, diag%ERROR_INDEX)
1961  if (diag%ERROR_INDEX < 0) then
1962  diag%ERROR_MESSAGE = 'in convert, error from set_control_parameter, failed to set hyperdiffusion_temperature_explicit'
1963  return
1964  end if
1965  call set_control_parameter(cp, 'hyperdiffusion_temperature_rho_tor_norm_cutoff' , profiles%HYPERDIFFUSION%TEMP%RHO_NORM_CUTOFF, diag%ERROR_INDEX)
1966  if (diag%ERROR_INDEX < 0) then
1967  diag%ERROR_MESSAGE = 'in convert, error from set_control_parameter, failed to set hyperdiffusion_temperature_rho_tor_norm_cutoff'
1968  return
1969  end if
1970  call set_control_parameter(cp, 'hyperdiffusion_rotation_implicit' , profiles%HYPERDIFFUSION%ROTATION%IMPLICIT, diag%ERROR_INDEX)
1971  if (diag%ERROR_INDEX < 0) then
1972  diag%ERROR_MESSAGE = 'in convert, error from set_control_parameter, failed to set hyperdiffusion_rotation_implicit'
1973  return
1974  end if
1975  call set_control_parameter(cp, 'hyperdiffusion_rotation_explicit' , profiles%HYPERDIFFUSION%ROTATION%EXPLICIT, diag%ERROR_INDEX)
1976  if (diag%ERROR_INDEX < 0) then
1977  diag%ERROR_MESSAGE = 'in convert, error from set_control_parameter, failed to set hyperdiffusion_rotation_explicit'
1978  return
1979  end if
1980  call set_control_parameter(cp, 'hyperdiffusion_rotation_rho_tor_norm_cutoff' , profiles%HYPERDIFFUSION%ROTATION%RHO_NORM_CUTOFF, diag%ERROR_INDEX)
1981  if (diag%ERROR_INDEX < 0) then
1982  diag%ERROR_MESSAGE = 'in convert, error from set_control_parameter, failed to set '
1983  return
1984  end if
1985 
1986  call set_control_parameter(cp, 'mixing_ratio_kinetic_profiles', numerics%CONTROL%AMIX, diag%ERROR_INDEX)
1987  if (diag%ERROR_INDEX < 0) then
1988  diag%ERROR_MESSAGE = 'in convert, error from set_control_parameter, failed to set mixing_ratio_kinetic_profiles'
1989  return
1990  end if
1991  call set_control_parameter(cp, 'mixing_ratio_transport', numerics%CONTROL%AMIX_TRANSP, diag%ERROR_INDEX)
1992  if (diag%ERROR_INDEX < 0) then
1993  diag%ERROR_MESSAGE = 'in convert, error from set_control_parameter, failed to set mixing_ratio_transport'
1994  return
1995  end if
1996  call set_control_parameter(cp, 'mixing_ratio_sources', numerics%CONTROL%AMIX_SRC, diag%ERROR_INDEX)
1997  if (diag%ERROR_INDEX < 0) then
1998  diag%ERROR_MESSAGE = 'in convert, error from set_control_parameter, failed to set mixing_ratio_sources'
1999  return
2000  end if
2001  call set_control_parameter(cp, 'ohmic_power_multiplier', numerics%CONTROL%MULTIPLIER_OH, diag%ERROR_INDEX)
2002  if (diag%ERROR_INDEX < 0) then
2003  diag%ERROR_MESSAGE = 'in convert, error from set_control_parameter, failed to set ohmic_power_multiplier'
2004  return
2005  end if
2006  end associate
2007 
2008  ieq = 0
2009  DO isol = 1,SIZE(numerics%SOLVER)
2010  DO idim = 1, numerics%SOLVER(isol)%NDIM
2011  ieq = ieq + 1
2012 
2013  coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%primary_quantity%identifier%name = numerics%SOLVER(isol)%eq_name
2014  coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%primary_quantity%identifier%index = numerics%SOLVER(isol)%eq_ind
2015 
2016  coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%primary_quantity%ion_index = numerics%SOLVER(isol)%comp_ind(idim)
2017  coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%primary_quantity%state_index = numerics%SOLVER(isol)%state_ind(idim)
2018 
2019  coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%primary_quantity%profile(:) = numerics%SOLVER(isol)%y(idim,:)
2020  coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%primary_quantity%d_dr(:) = numerics%SOLVER(isol)%dy(idim,:)
2021  coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%primary_quantity%d2_dr2(:) = numerics%SOLVER(isol)%ddy(idim,:)
2022 
2023  coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%coefficient(1)%profile(:) = numerics%SOLVER(isol)%a(idim,:)
2024  coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%coefficient(2)%profile(:) = numerics%SOLVER(isol)%b(idim,:)
2025  coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%coefficient(3)%profile(:) = numerics%SOLVER(isol)%c(idim,:)
2026  coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%coefficient(4)%profile(:) = numerics%SOLVER(isol)%d(idim,:)
2027  coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%coefficient(5)%profile(:) = numerics%SOLVER(isol)%e(idim,:)
2028  coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%coefficient(6)%profile(:) = numerics%SOLVER(isol)%f(idim,:)
2029  coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%coefficient(7)%profile(:) = numerics%SOLVER(isol)%g(idim,:)
2030  coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%coefficient(8)%profile(:) = numerics%SOLVER(isol)%h
2031 
2032  DO ibnd=1,SIZE(coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%boundary_condition)
2033  IF (SIZE(coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%boundary_condition)==1) THEN
2034  ibnd0=2
2035  ELSE
2036  ibnd0=ibnd
2037  END IF
2038 
2039  CALL bc_type_to_ids( numerics%SOLVER(isol)%bc_type(idim,ibnd0), &
2040  coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%computation_mode%index, &
2041  coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%computation_mode%name(1), &
2042  coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%boundary_condition(ibnd)%type%index, &
2043  coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%boundary_condition(ibnd)%type%name(1) )
2044 
2045  coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%boundary_condition(ibnd)%position = numerics%SOLVER(isol)%BC_RHO_NORM(idim,ibnd0)
2046  coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%boundary_condition(ibnd)%value(:) = numerics%SOLVER(isol)%bc_value(idim,ibnd0,:)
2047  END DO
2048 
2049  coresolver_new%SOLVER_1D(1)%EQUATION(ieq)%convergence%delta_relative%value = numerics%SOLVER(isol)%ERROR(idim)
2050 
2051  END DO
2052  END DO
2053 
2054 
2055  RETURN
2056 
2057 
2058  END SUBROUTINE convert_internal_to_ids_types
2059 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
2060 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
2061 
2062 
2063  SUBROUTINE check_coreprof_grids (COREPROFILES)
2064 
2065  USE ids_schemas
2066 
2067  IMPLICIT NONE
2068 
2069  TYPE (ids_core_profiles) :: coreprofiles
2070  IF (.NOT.ASSOCIATED(coreprofiles%profiles_1d(1)%grid%rho_tor_norm)) THEN
2071  ALLOCATE (coreprofiles%profiles_1d(1)%grid%rho_tor_norm(SIZE(coreprofiles%profiles_1d(1)%grid%rho_tor)))
2072  coreprofiles%profiles_1d(1)%grid%rho_tor_norm = coreprofiles%profiles_1d(1)%grid%rho_tor &
2073  /coreprofiles%profiles_1d(1)%grid%rho_tor(SIZE(coreprofiles%profiles_1d(1)%grid%rho_tor))
2074  ELSE IF (maxval(coreprofiles%profiles_1d(1)%grid%rho_tor_norm).LE.0.0_ids_real) THEN
2075  coreprofiles%profiles_1d(1)%grid%rho_tor_norm = coreprofiles%profiles_1d(1)%grid%rho_tor &
2076  /coreprofiles%profiles_1d(1)%grid%rho_tor(SIZE(coreprofiles%profiles_1d(1)%grid%rho_tor))
2077  END IF
2078 
2079  END SUBROUTINE check_coreprof_grids
2080 
2081  SUBROUTINE check_coretransp_grids (CORETRANSPORT)
2082 
2083  USE ids_schemas
2084 
2085  IMPLICIT NONE
2086 
2087  TYPE (ids_core_transport) :: coretransport
2088  INTEGER :: i
2089  DO i = 1, SIZE(coretransport%MODEL)
2090  IF (.NOT.ASSOCIATED(coretransport%MODEL(i)%profiles_1d(1)%grid_d%rho_tor_norm)) THEN
2091  ALLOCATE (coretransport%MODEL(i)%profiles_1d(1)%grid_d%rho_tor_norm(SIZE(coretransport%MODEL(i)%profiles_1d(1)%grid_d%rho_tor)))
2092  coretransport%MODEL(i)%profiles_1d(1)%grid_d%rho_tor_norm = coretransport%MODEL(i)%profiles_1d(1)%grid_d%rho_tor &
2093  /coretransport%MODEL(i)%profiles_1d(1)%grid_d%rho_tor(SIZE(coretransport%MODEL(i)%profiles_1d(1)%grid_d%rho_tor))
2094  ELSE IF (maxval(coretransport%MODEL(i)%profiles_1d(1)%grid_d%rho_tor_norm).LE.0.0_ids_real) THEN
2095  coretransport%MODEL(i)%profiles_1d(1)%grid_d%rho_tor_norm = coretransport%MODEL(i)%profiles_1d(1)%grid_d%rho_tor &
2096  /coretransport%MODEL(i)%profiles_1d(1)%grid_d%rho_tor(SIZE(coretransport%MODEL(i)%profiles_1d(1)%grid_d%rho_tor))
2097  END IF
2098  IF (.NOT.ASSOCIATED(coretransport%MODEL(i)%profiles_1d(1)%grid_v%rho_tor_norm)) THEN
2099  ALLOCATE (coretransport%MODEL(i)%profiles_1d(1)%grid_v%rho_tor_norm(SIZE(coretransport%MODEL(i)%profiles_1d(1)%grid_v%rho_tor)))
2100  coretransport%MODEL(i)%profiles_1d(1)%grid_v%rho_tor_norm = coretransport%MODEL(i)%profiles_1d(1)%grid_v%rho_tor &
2101  /coretransport%MODEL(i)%profiles_1d(1)%grid_v%rho_tor(SIZE(coretransport%MODEL(i)%profiles_1d(1)%grid_v%rho_tor))
2102  ELSE IF (maxval(coretransport%MODEL(i)%profiles_1d(1)%grid_v%rho_tor_norm).LE.0.0_ids_real) THEN
2103  coretransport%MODEL(i)%profiles_1d(1)%grid_v%rho_tor_norm = coretransport%MODEL(i)%profiles_1d(1)%grid_v%rho_tor &
2104  /coretransport%MODEL(i)%profiles_1d(1)%grid_v%rho_tor(SIZE(coretransport%MODEL(i)%profiles_1d(1)%grid_v%rho_tor))
2105  END IF
2106  IF (.NOT.ASSOCIATED(coretransport%MODEL(i)%profiles_1d(1)%grid_flux%rho_tor_norm)) THEN
2107  ALLOCATE (coretransport%MODEL(i)%profiles_1d(1)%grid_flux%rho_tor_norm(SIZE(coretransport%MODEL(i)%profiles_1d(1)%grid_flux%rho_tor)))
2108  coretransport%MODEL(i)%profiles_1d(1)%grid_flux%rho_tor_norm = coretransport%MODEL(i)%profiles_1d(1)%grid_flux%rho_tor &
2109  /coretransport%MODEL(i)%profiles_1d(1)%grid_flux%rho_tor(SIZE(coretransport%MODEL(i)%profiles_1d(1)%grid_flux%rho_tor))
2110  ELSE IF (maxval(coretransport%MODEL(i)%profiles_1d(1)%grid_flux%rho_tor_norm).LE.0.0_ids_real) THEN
2111  coretransport%MODEL(i)%profiles_1d(1)%grid_flux%rho_tor_norm = coretransport%MODEL(i)%profiles_1d(1)%grid_flux%rho_tor &
2112  /coretransport%MODEL(i)%profiles_1d(1)%grid_flux%rho_tor(SIZE(coretransport%MODEL(i)%profiles_1d(1)%grid_flux%rho_tor))
2113  END IF
2114  END DO
2115 
2116  END SUBROUTINE check_coretransp_grids
2117 
2118  SUBROUTINE check_coresource_grids (CORESOURCES)
2119 
2120  USE ids_schemas
2121 
2122  IMPLICIT NONE
2123 
2124  TYPE (ids_core_sources) :: coresources
2125  INTEGER :: i
2126  DO i = 1, SIZE(coresources%SOURCE)
2127  IF (.NOT.ASSOCIATED(coresources%SOURCE(i)%profiles_1d(1)%grid%rho_tor_norm)) THEN
2128  ALLOCATE (coresources%SOURCE(i)%profiles_1d(1)%grid%rho_tor_norm(SIZE(coresources%SOURCE(i)%profiles_1d(1)%grid%rho_tor)))
2129  coresources%SOURCE(i)%profiles_1d(1)%grid%rho_tor_norm = coresources%SOURCE(i)%profiles_1d(1)%grid%rho_tor &
2130  /coresources%SOURCE(i)%profiles_1d(1)%grid%rho_tor(SIZE(coresources%SOURCE(i)%profiles_1d(1)%grid%rho_tor))
2131  ELSE IF (maxval(coresources%SOURCE(i)%profiles_1d(1)%grid%rho_tor_norm).LE.0.0_ids_real) THEN
2132  coresources%SOURCE(i)%profiles_1d(1)%grid%rho_tor_norm = coresources%SOURCE(i)%profiles_1d(1)%grid%rho_tor &
2133  /coresources%SOURCE(i)%profiles_1d(1)%grid%rho_tor(SIZE(coresources%SOURCE(i)%profiles_1d(1)%grid%rho_tor))
2134  END IF
2135  END DO
2136 
2137  END SUBROUTINE check_coresource_grids
2138 
2139 
2140 END MODULE convert
2141 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
2142 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
subroutine allocate_coreprofiles_ids(NTIME, NRHO, COMPOSITION, COREPROFILES)
This routine allocates COREPROFILES IDS.
subroutine allocate_magnetic_geometry(NRHO, GEOMETRY, ifail)
Definition: ets_plasma.f90:612
subroutine init_control_parameters(param, default_int, default_real)
IMPURITY.
Definition: impurity.f90:8
subroutine convert_ids_to_internal_types(EQUILIBRIUM_OLD, EQUILIBRIUM_ITER, COREPROFILES_OLD, COREPROFILES_ITER,CORETRANSPORT_ITER, CORESOURCE_ITER, CORESOLVER_ITER,
Definition: convert.f90:11
This module contains routines for allocation/deallocation if IDSs used in ETS.
subroutine allocate_plasma_profiles(NRHO, NION, PROFILES, ifail)
Definition: ets_plasma.f90:725
subroutine check_coretransp_grids(CORETRANSPORT)
Definition: convert.f90:2081
logical function, public acceptable_temperature(array)
subroutine get_source_composition(CORESOURCE, SOURCE, COMPOSITION)
subroutine allocate_sources_and_sinks(NRHO, NION, SOURCES, ifail)
logical function, public acceptable_impurity_density(array)
subroutine allocate_coretransport_ids(NTIME, NMOD, NRHO, COMPOSITION, CORETRANSPORT)
subroutine get_prof_composition(COREPROFILES, COMPOSITION)
This routine detects composition of CORE_PROFILES IDS.
subroutine find_ion(ION_INDEX_IN, COMPOSITION_IN, COMPOSITION_OUT, ION_INDEX_OUT)
subroutine convert_internal_to_ids_types(GEOMETRY, PROFILES, TRANSPORT, SOURCES, IMPURITY, GLOBAL, NUMERICS, COMPOSITION, COREPROFILES_NEW, CORETRANSPORT_NEW, CORESOURCE_NEW, CORESOLVER_NEW, DIAG)
This routine converts internal ETS into IDS derived types.
Definition: convert.f90:1338
subroutine ids_to_bc_type(IDS_COMPUTATION_MODE, IDS_BC, INTERNAL_BC)
subroutine allocate_impurity_profiles(NRHO, NIMP, MAX_NZIMP, IMPURITY, ifail)
This module contains routines for detecting plasma composition in IDSs.
subroutine allocate_coresources_ids(NTIME, NSRC, NRHO, COMPOSITION, CORESOURCES)
subroutine allocate_transport_coefficients(NRHO, NION, TRANSPORT, ifail)
logical function, public acceptable_ion_density(array)
logical function, public acceptable_electron_density(array)
subroutine l3deriv_pointer(YIN, XIN, YOUT, XOUT)
Calculate the derivative of YIN(XIN) using L3INTERP, but the input arrays are pointers.
subroutine find_impurity(IMP_INDEX_IN, COMPOSITION_IN, COMPOSITION_OUT, IMP_INDEX_OUT)
subroutine l3interp_pointer(YIN, XIN, YOUT, XOUT)
Translate number (reals or integers) into strings, i.e. character(*)
Definition: stringtools.f90:8
subroutine lininterp(YIN, XIN, NIN, YOUT, XOUT, NOUT)
Linear interpolatation and constant extrapolation. The input function is YIN given at XIN and with le...
subroutine check_coreprof_grids(COREPROFILES)
Definition: convert.f90:2063
Interpolations, derivatives and integrals of 1d functions. Includes wrappers to L3INTERP, L3DERIV and INTERPOS.
The module defines derived types used by ETS6-CoreActor and subroutines to allocate and deallocate in...
Definition: ets_plasma.f90:26
subroutine get_transp_composition(CORETRANSPORT, MODEL, COMPOSITION)
subroutine check_coresource_grids(CORESOURCES)
Definition: convert.f90:2118
subroutine bc_type_to_ids(INTERNAL_BC, IDS_COMPUTATION_MODE_INDEX, IDS_COMPUTATION_MODE_NAME, IDS_BC_INDEX, IDS_BC_NAME)
The module declares types of variables used by numerical solvers.
subroutine allocate_coresolver_ids(NTIME, NEQ, NCOEF, NRHO, NBC, CORESOLVER)