ETS  \$Id: Doxyfile 2162 2020-02-26 14:16:09Z g2dpc $
 All Classes Files Functions Variables Pages
analytics1.f90
Go to the documentation of this file.
1 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
7 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
8 MODULE analytics1
9 
10 CONTAINS
11 
12 
13 
14 
15 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
24 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
25  SUBROUTINE analytical_plasma ( &
26  time,coreprof_in, &
27  equilibrium,coreprof_analytic,coretransp,coresource,coreimpur, code_parameters)
28 
29 
30 !-------------------------------------------------------!
31 ! This routine manufactures the solution for the !
32 ! set of transport equations describing the main !
33 ! plasma. It provides sources and transport !
34 ! coefficients for all transport equations. !
35 !-------------------------------------------------------!
36 ! Source: --- !
37 ! Developers: R.Stankiewicz, D.Kalupin !
38 ! Kontacts: Roman.Stankiewich@gmail.com !
39 ! D.Kalupin@fz-juelich.de !
40 ! !
41 ! Comments: --- !
42 ! !
43 !-------------------------------------------------------!
44 
45 
46  USE itm_types
47  USE euitm_schemas
48  USE ets_plasma
51  USE copy_cpo_ets
52  USE is_set_structures
53  USE copy_structures
54  use deallocate_structures
55 !! use euitm_printcpo
56 
57  IMPLICIT NONE
58 
59 
60  INTEGER, SAVE :: shot_no !shot mumber
61  INTEGER, SAVE :: run_no !run mumber
62 
63  INTEGER, SAVE :: nrho !number of radial points (input)
64  INTEGER, SAVE :: nnucl !number of ion species (input)
65  INTEGER, SAVE :: nneut !number of ion species (input)
66  INTEGER, SAVE :: nion !number of ion species (input)
67  INTEGER, SAVE :: nimp !number of impurity species (input)
68  INTEGER, ALLOCATABLE, SAVE :: nzimp(:) !number of charge states for each impurity (input)
69  INTEGER, ALLOCATABLE, SAVE :: ncomp(:) !number of charge states for each impurity (input)
70  INTEGER, ALLOCATABLE, SAVE :: ntype(:) !number of charge states for each impurity (input)
71  INTEGER, SAVE :: npoints !number of charge states for each impurity (input)
72  INTEGER, SAVE :: ntime !number of time points (input)
73  INTEGER, PARAMETER :: nocur = 1 !number of CPO ocurancies in the work flow
74 
75  INTEGER :: irho !current radial knot
76  INTEGER :: iion !current ion type
77  INTEGER :: itime !current time step
78 
79  INTEGER, SAVE :: nsol !Number of analytical example
80  INTEGER, SAVE :: solver_type !representation of transport equations
81  !1-"standard"; 2-"integral"(default)
82  INTEGER, SAVE :: sigma_source !origin of Plasma electrical conductivity
83  !0: plasma collisions;1: transport module;2: source module
84 
85  REAL (R8), SAVE :: convrec !required convergency
86  REAL (R8) :: time !Time
87  REAL (R8), SAVE :: tau !time step, and mixing coefficient
88  REAL (R8), SAVE :: amix !mixing factor
89  INTEGER :: iter !iteration index
90  INTEGER, PARAMETER :: maxiter=1000 !maximum number of convergence iterations
91 
92  INTEGER, SAVE :: psi_bnd_type !Type of boundary conditions current
93  INTEGER, SAVE :: ni_bnd_type !Type of boundary conditions ion density
94  INTEGER, SAVE :: ti_bnd_type !Type of boundary conditions ion temperature
95  INTEGER, SAVE :: te_bnd_type !Type of boundary conditions electron temperature
96  INTEGER, SAVE :: vtor_bnd_type !Type of boundary conditions toroidal rotation
97 
98  TYPE (run_control) :: control !contains all parameters required by run
99 
100  INTEGER :: ifail
101 
102 
103 ! +++ CPO derived types:
104  TYPE (type_equilibrium),POINTER :: equilibrium(:) !input CPO with geometry
105  TYPE (type_coreprof), POINTER :: coreprof_in(:) !CPO with internal ETS parameters profiles from previous time step
106  TYPE (type_coreprof), POINTER :: coreprof_analytic(:) !CPO with profiles
107  TYPE (type_coretransp), POINTER :: coretransp(:) !CPO with transport coefficients
108  TYPE (type_coresource), POINTER :: coresource(:) !CPO with sources
109  TYPE (type_coreimpur), POINTER :: coreimpur(:) !CPO with impurities
110  TYPE (type_param) :: code_parameters
111 
112 
113 ! +++ Local derived types:
114  TYPE (magnetic_geometry) :: geometry1 !contains all geometry quantities
115  TYPE (plasma_profiles) :: profiles1 !contains profiles of plasma parameters
116  TYPE (transport_coefficients) :: transport1 !contains profiles of trasport coefficients
117  TYPE (sources_and_sinks) :: sources1 !contains profiles of sources
118  TYPE (global_param) :: global !contains global
119 
120  CHARACTER (len=32) :: database_format
121 
122  LOGICAL, SAVE :: first=.true.
123 
124 ! +++ Set up dimensions:
125 
126  WRITE(6,*) time
127 !! call is_set_cpo(coreprof_in(1),'coreprof_in')
128 !!! flush(6)
129 
130  IF(first) THEN
131  IF (.NOT. ASSOCIATED(code_parameters%parameters)) THEN
132  WRITE(6, *) 'ERROR: code parameters not associated!'
133  stop
134  END IF
135  CALL process_xml(solver_type,sigma_source,tau,amix,convrec,nrho,nion,nimp,nzimp,ntime,nsol, &
136  psi_bnd_type,ni_bnd_type,ti_bnd_type,te_bnd_type,vtor_bnd_type, &
137  shot_no, run_no, code_parameters, database_format)
138  first=.false.
139  ENDIF
140 
141  npoints = 101
142  nnucl = nion + nimp
143  nneut = nion + nimp
144  if(.not. allocated(ntype)) ALLOCATE (ntype(nneut))
145  if(.not. allocated(ncomp)) ALLOCATE (ncomp(nneut))
146  ntype = 1
147  ncomp = 1
148 
149  CALL allocate_equilibrium_cpo(nocur, nrho, nrho, 101, npoints, equilibrium )
150 
151  CALL allocate_coreprof_cpo(nocur, nrho, nnucl, nion, nimp, nzimp, nneut, ntype, ncomp, coreprof_analytic )
152  CALL allocate_coretransp_cpo(nocur, nrho, nnucl, nion, nimp, nzimp, nneut, ntype, ncomp, coretransp )
153  CALL allocate_coresource_cpo(nocur, nrho, nnucl, nion, nimp, nzimp, nneut, ntype, ncomp, coresource )
154 
155  if(nimp .GT. 0) &
156  CALL allocate_coreimpur_cpo(nocur, nrho, nnucl, nion, nimp, nzimp, nneut, ntype, ncomp, coreimpur )
157 
158 
159  CALL copy_codeparam(coreprof_in(1)%codeparam, coreprof_analytic(1)%codeparam)
160  call deallocate_cpo(coreprof_analytic(1)%COMPOSITIONS)
161  CALL copy_cpo(coreprof_in(1)%COMPOSITIONS, coreprof_analytic(1)%COMPOSITIONS)
162  call deallocate_cpo(coretransp(1)%COMPOSITIONS)
163  CALL copy_cpo(coreprof_in(1)%COMPOSITIONS, coretransp(1)%COMPOSITIONS)
164  call deallocate_cpo(coresource(1)%COMPOSITIONS)
165  CALL copy_cpo(coreprof_in(1)%COMPOSITIONS, coresource(1)%COMPOSITIONS)
166 
167 ! +++ Set up parameters for the test work flow:
168 ! CALL perfon ('setcpo')
169  CALL set_cpo(nrho,nion,nimp,nzimp,ntime,nsol, &
170  psi_bnd_type,ni_bnd_type,ti_bnd_type,te_bnd_type,vtor_bnd_type, &
171  equilibrium, coreprof_in, coretransp)
172 ! CALL perfoff
173 
174 
175 ! +++ Allocation of local derived types:
176  CALL allocate_magnetic_geometry(nrho, geometry1, ifail)
177  CALL allocate_plasma_profiles(nrho,nion, profiles1, ifail)
178  CALL allocate_transport_coefficients(nrho,nion, transport1, ifail)
179  CALL allocate_sources_and_sinks(nrho,nion, sources1, ifail)
180  CALL allocate_global_param( global, ifail)
181 
182 
183 
184 ! +++ Copy CPOs in local derived types:
185  CALL convert_cpo_to_local_types(coreprof_in, profiles1, transport1)
186 
187 
188 
189 ! +++ Call MAIN_PLASMA with internal ETS derived types:
190  CALL analytics_full(time, geometry1, profiles1, transport1, sources1)
191 
192 
193 
194 ! +++ Copy local derived types in CPOs:
195  CALL convert_local_to_cpo_types(geometry1, profiles1, transport1, sources1, &
196  equilibrium, coreprof_analytic, coretransp, coresource)
197 
198 
199 ! +++ Deallocation of ETS derived types:
200  CALL deallocate_magnetic_geometry(geometry1, ifail)
201  CALL deallocate_plasma_profiles(profiles1, ifail)
202  CALL deallocate_transport_coefficients(transport1, ifail)
203  CALL deallocate_sources_and_sinks(sources1, ifail)
204 
205 
206  WRITE(6,*) 'EQUILIBRIUM'
207  IF(.NOT.ASSOCIATED(equilibrium(1)%codeparam%codename)) THEN
208  ALLOCATE(equilibrium(1)%codeparam%codename(1)) ! For a string of 132 characters max.
209  equilibrium(1)%codeparam%codename(1) = 'ANALYTIC_SOLVER_TEST'
210  ENDIF
211  IF(.NOT.ASSOCIATED(equilibrium(1)%codeparam%codeversion)) THEN
212  ALLOCATE(equilibrium(1)%codeparam%codeversion(1)) ! For a string of 132 characters max.
213  equilibrium(1)%codeparam%codeversion(1) = '?'
214  ENDIF
215  IF(.NOT.ASSOCIATED(equilibrium(1)%codeparam%parameters)) THEN
216  ALLOCATE(equilibrium(1)%codeparam%parameters(1)) ! For a string of 132 characters max.
217  equilibrium(1)%codeparam%parameters(1) = 'my_code_specific_parameters'
218  ENDIF
219  IF(.NOT.ASSOCIATED(equilibrium(1)%codeparam%output_diag)) THEN
220  ALLOCATE(equilibrium(1)%codeparam%output_diag(1)) ! For a string of 132 characters max.
221  equilibrium(1)%codeparam%output_diag(1) = 'my_output_diag'
222  ENDIF
223  equilibrium(1)%codeparam%output_flag = 0 ! Integer output flag, 0 means the run was successful and can be used in the rest of the workflow, <0 means failure
224  equilibrium(1)%time=time
225 !! call printcpo(EQUILIBRIUM(1))
226  WRITE(6,*) 'COREPROF_ANALYTIC'
227  IF(.NOT.ASSOCIATED(coreprof_analytic(1)%codeparam%codename)) THEN
228  ALLOCATE(coreprof_analytic(1)%codeparam%codename(1)) ! For a string of 132 characters max.
229  coreprof_analytic(1)%codeparam%codename(1) = 'ANALYTIC_SOLVER_TEST'
230  ENDIF
231  IF(.NOT.ASSOCIATED(coreprof_analytic(1)%codeparam%codeversion)) THEN
232  ALLOCATE(coreprof_analytic(1)%codeparam%codeversion(1)) ! For a string of 132 characters max.
233  coreprof_analytic(1)%codeparam%codeversion(1) = '?'
234  ENDIF
235  IF(.NOT.ASSOCIATED(coreprof_analytic(1)%codeparam%parameters)) THEN
236  ALLOCATE(coreprof_analytic(1)%codeparam%parameters(1)) ! For a string of 132 characters max.
237  coreprof_analytic(1)%codeparam%parameters(1) = 'my_code_specific_parameters'
238  ENDIF
239  IF(.NOT.ASSOCIATED(coreprof_analytic(1)%codeparam%output_diag)) THEN
240  ALLOCATE(coreprof_analytic(1)%codeparam%output_diag(1)) ! For a string of 132 characters max.
241  coreprof_analytic(1)%codeparam%output_diag(1) = 'my_output_diag'
242  ENDIF
243  coreprof_analytic(1)%codeparam%output_flag = 0 ! Integer output flag, 0 means the run was successful and can be used in the rest of the workflow, <0 means failure
244  coreprof_analytic(1)%time=time
245 !! call printcpo(COREPROF_ANALYTIC(1))
246  WRITE(6,*) 'CORETRANSP'
247  IF(.NOT.ASSOCIATED(coretransp(1)%codeparam%codename)) THEN
248  ALLOCATE(coretransp(1)%codeparam%codename(1)) ! For a string of 132 characters max.
249  coretransp(1)%codeparam%codename(1) = 'ANALYTIC_SOLVER_TEST'
250  ENDIF
251  IF(.NOT.ASSOCIATED(coretransp(1)%codeparam%codeversion)) THEN
252  ALLOCATE(coretransp(1)%codeparam%codeversion(1)) ! For a string of 132 characters max.
253  coretransp(1)%codeparam%codeversion(1) = '?'
254  ENDIF
255  IF(.NOT.ASSOCIATED(coretransp(1)%codeparam%parameters)) THEN
256  ALLOCATE(coretransp(1)%codeparam%parameters(1)) ! For a string of 132 characters max.
257  coretransp(1)%codeparam%parameters(1) = 'my_code_specific_parameters'
258  ENDIF
259  IF(.NOT.ASSOCIATED(coretransp(1)%codeparam%output_diag)) THEN
260  ALLOCATE(coretransp(1)%codeparam%output_diag(1)) ! For a string of 132 characters max.
261  coretransp(1)%codeparam%output_diag(1) = 'my_output_diag'
262  ENDIF
263  coretransp(1)%codeparam%output_flag = 0 ! Integer output flag, 0 means the run was successful and can be used in the rest of the workflow, <0 means failure
264  coretransp(1)%time=time
265 !! call printcpo(CORETRANSP(1))
266  WRITE(6,*) 'CORESOURCE'
267  IF(.NOT.ASSOCIATED(coresource(1)%codeparam%codename)) THEN
268  ALLOCATE(coresource(1)%codeparam%codename(1)) ! For a string of 132 characters max.
269  coresource(1)%codeparam%codename(1) = 'ANALYTIC_SOLVER_TEST'
270  ENDIF
271  IF(.NOT.ASSOCIATED(coresource(1)%codeparam%codeversion)) THEN
272  ALLOCATE(coresource(1)%codeparam%codeversion(1)) ! For a string of 132 characters max.
273  coresource(1)%codeparam%codeversion(1) = '?'
274  ENDIF
275  IF(.NOT.ASSOCIATED(coresource(1)%codeparam%parameters)) THEN
276  ALLOCATE(coresource(1)%codeparam%parameters(1)) ! For a string of 132 characters max.
277  coresource(1)%codeparam%parameters(1) = 'my_code_specific_parameters'
278  ENDIF
279  IF(.NOT.ASSOCIATED(coresource(1)%codeparam%output_diag)) THEN
280  ALLOCATE(coresource(1)%codeparam%output_diag(1)) ! For a string of 132 characters max.
281  coresource(1)%codeparam%output_diag(1) = 'my_output_diag'
282  ENDIF
283  coresource(1)%codeparam%output_flag = 0 ! Integer output flag, 0 means the run was successful and can be used in the rest of the workflow, <0 means failure
284  coresource(1)%time=time
285 !! call printcpo(CORESOURCE(1))
286  if(nimp.gt.0) then
287  WRITE(6,*) 'COREIMPUR'
288  IF(.NOT.ASSOCIATED(coreimpur(1)%codeparam%codename)) THEN
289  ALLOCATE(coreimpur(1)%codeparam%codename(1)) ! For a string of 132 characters max.
290  coreimpur(1)%codeparam%codename(1) = 'ANALYTIC_SOLVER_TEST'
291  ENDIF
292  IF(.NOT.ASSOCIATED(coreimpur(1)%codeparam%codeversion)) THEN
293  ALLOCATE(coreimpur(1)%codeparam%codeversion(1)) ! For a string of 132 characters max.
294  coreimpur(1)%codeparam%codeversion(1) = '?'
295  ENDIF
296  IF(.NOT.ASSOCIATED(coreimpur(1)%codeparam%parameters)) THEN
297  ALLOCATE(coreimpur(1)%codeparam%parameters(1)) ! For a string of 132 characters max.
298  coreimpur(1)%codeparam%parameters(1) = 'my_code_specific_parameters'
299  ENDIF
300  IF(.NOT.ASSOCIATED(coreimpur(1)%codeparam%output_diag)) THEN
301  ALLOCATE(coreimpur(1)%codeparam%output_diag(1)) ! For a string of 132 characters max.
302  coreimpur(1)%codeparam%output_diag(1) = 'my_output_diag'
303  ENDIF
304  coreimpur(1)%codeparam%output_flag = 0 ! Integer output flag, 0 means the run was successful and can be used in the rest of the workflow, <0 means failure
305  coreimpur(1)%time=time
306 !! call printcpo(COREIMPUR(1))
307  endif
308  WRITE(6,*) 'Finished in ANALYTICAL_PLASMA'
309 !!! flush(6)
310 
311  RETURN
312 
313 
314 
315  END SUBROUTINE analytical_plasma
316 
317 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
318 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
319 
320 
321 
322 
323 
324 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
330 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
331  SUBROUTINE convert_cpo_to_local_types (COREPROF, PROFILES1, TRANSPORT1)
332 
333 !-------------------------------------------------------!
334 ! This routine converts CPOs into the local !
335 ! derived types. !
336 !-------------------------------------------------------!
337 ! Source: --- !
338 ! Developers: D.Kalupin !
339 ! Kontacts: D.Kalupin@fz-juelich.de !
340 ! !
341 ! Comments: --- !
342 ! !
343 !-------------------------------------------------------!
344 
345 
346  USE euitm_schemas
347 
348  USE ets_plasma
349 
350  IMPLICIT NONE
351 
352  INTEGER :: iion, nion
353 
354 ! +++ CPO derived types:
355  TYPE (type_coreprof), POINTER :: coreprof(:) !input CPO with plasma profiles
356 ! +++ Local derived types:
357  TYPE (plasma_profiles), INTENT(INOUT) :: profiles1 !contains profiles of plasma parameters
358  TYPE (transport_coefficients), INTENT(INOUT) :: transport1 !contains profiles of trasport coefficients
359 
360 
361 
362 ! +++ Convert profiles:
363  nion = size(coreprof(1)%compositions%ions)
364  DO iion =1, nion
365  profiles1%ZION(iion) = coreprof(1)%compositions%ions(iion)%zion
366  profiles1%ZION2(iion) = coreprof(1)%compositions%ions(iion)%zion**2
367  profiles1%MION(iion) = coreprof(1)%compositions%nuclei(coreprof(1)%compositions%ions(iion)%nucindex)%amn
368  END DO
369 
370  profiles1%PSI_BND_TYPE = coreprof(1)%psi%boundary%type
371 
372  profiles1%NI_BND_TYPE = coreprof(1)%ni%boundary%type
373 
374  profiles1%TI_BND_TYPE = coreprof(1)%ti%boundary%type
375 
376  profiles1%TE_BND_TYPE = coreprof(1)%te%boundary%type
377 
378  profiles1%VTOR_BND_TYPE = coreprof(1)%vtor%boundary%type
379 
380 
381 
382 ! +++ Convert transport:
383  transport1%C1(1) = 0.0e0_r8
384  transport1%C1(2) = 1.5e0_r8
385  transport1%C1(3) = 2.5e0_r8
386 
387 
388 
389 
390  RETURN
391 
392 
393  END SUBROUTINE convert_cpo_to_local_types
394 
395 
396 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
397 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
398 
399 
400 
401 
402 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
408 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
409  SUBROUTINE convert_local_to_cpo_types (GEOMETRY1, PROFILES1, TRANSPORT1, SOURCES1, &
410  equilibrium, coreprof, coretransp, coresource)
411 
412 !-------------------------------------------------------!
413 ! This routine converts local into the CPOs !
414 ! derived types. !
415 !-------------------------------------------------------!
416 ! Source: --- !
417 ! Developers: D.Kalupin !
418 ! Kontacts: D.Kalupin@fz-juelich.de !
419 ! !
420 ! Comments: --- !
421 ! !
422 !-------------------------------------------------------!
423 
424 
425  USE euitm_schemas
426 
427  USE ets_plasma
428 
429  USE itm_constants
430 
431  IMPLICIT NONE
432 
433 ! +++ Local derived types:
434  TYPE (magnetic_geometry), INTENT(IN) :: geometry1 !contains all geometry quantities
435  TYPE (plasma_profiles), INTENT(IN) :: profiles1 !contains profiles of plasma parameters
436  TYPE (transport_coefficients), INTENT(IN) :: transport1 !contains profiles of trasport coefficients
437  TYPE (sources_and_sinks), INTENT(IN) :: sources1 !contains profiles of sources
438 
439 ! +++ CPO derived types:
440  TYPE (type_equilibrium), POINTER :: equilibrium(:) !output CPO with geometry
441  TYPE (type_coreprof), POINTER :: coreprof(:) !output CPO with plasma profiles
442  TYPE (type_coretransp), POINTER :: coretransp(:) !output CPO with transport coefficients
443  TYPE (type_coresource), POINTER :: coresource(:) !output CPO with sources
444 
445 
446 
447 ! +++ Convert geometry:
448  equilibrium(1)%profiles_1d%rho_tor = geometry1%RHO
449 !!!DPC-EQ-4.09a
450  equilibrium(1)%profiles_1d%vprime = geometry1%VPR
451  equilibrium(1)%profiles_1d%gm3 = geometry1%G1
452  equilibrium(1)%profiles_1d%gm8 = geometry1%G2 ! added DPC 2012-04-27
453  equilibrium(1)%profiles_1d%gm2 = geometry1%G3
454  equilibrium(1)%profiles_1d%F_dia = geometry1%FDIA
455 
456  coreprof(1)%rho_tor = geometry1%RHO
457  coreprof(1)%toroid_field%r0 = geometry1%RGEO
458  coreprof(1)%toroid_field%b0 = geometry1%BGEO
459 
460  coretransp(1)%values(1)%rho_tor = geometry1%RHO
461  coresource(1)%values(1)%rho_tor = geometry1%RHO
462 
463 
464 
465 ! +++ Convert profiles:
466 ! allocate(COREPROF(1)%compositions%ions(NION))
467 ! DO IION =1, NION
468 ! COREPROF(1)%compositions%ions(iion)%zion = PROFILES1%ZION(IION)
469 ! COREPROF(1)%composition%amn = PROFILES1%MION
470 ! END DO
471 
472  coreprof(1)%psi%value = profiles1%PSI
473  coreprof(1)%profiles1d%q%value = profiles1%QSF
474  coreprof(1)%psi%boundary%value = profiles1%PSI_BND
475  coreprof(1)%psi%boundary%type = profiles1%PSI_BND_TYPE
476 
477  coreprof(1)%ni%value = profiles1%NI
478  coreprof(1)%ni%boundary%value = profiles1%NI_BND
479  coreprof(1)%ni%boundary%type = profiles1%NI_BND_TYPE
480 !DPC
481 ! COREPROF(1)%ni%transp_coef%diff = sum(TRANSPORT1%DIFF_NI,dim=3)
482 ! COREPROF(1)%ni%transp_coef%vconv = sum(TRANSPORT1%VCONV_NI,dim=3)
483 !DPC.end
484 
485  coreprof(1)%ne%value = profiles1%NE
486 !DPC
487 ! COREPROF(1)%ne%transp_coef%diff = TRANSPORT1%DIFF_NI
488 ! COREPROF(1)%ne%transp_coef%vconv = TRANSPORT1%VCONV_NI
489 !DPC.end
490 
491  coreprof(1)%ti%value = profiles1%TI
492  coreprof(1)%ti%boundary%value = profiles1%TI_BND
493  coreprof(1)%ti%boundary%type = profiles1%TI_BND_TYPE
494 !DPC
495  coreprof(1)%ti%transp_coef%diff = transport1%DIFF_TI
496  coreprof(1)%ti%transp_coef%vconv = transport1%VCONV_TI
497 !DPC.end
498 
499  coreprof(1)%te%value = profiles1%TE
500  coreprof(1)%te%boundary%value = profiles1%TE_BND
501  coreprof(1)%te%boundary%type = profiles1%TE_BND_TYPE
502 !DPC
503  coreprof(1)%te%transp_coef%diff = transport1%DIFF_TE
504  coreprof(1)%te%transp_coef%vconv = transport1%VCONV_TE
505 !DPC.end
506 
507  coreprof(1)%vtor%value = profiles1%VTOR
508  coreprof(1)%vtor%boundary%value = profiles1%VTOR_BND
509  coreprof(1)%vtor%boundary%type = profiles1%VTOR_BND_TYPE
510 !DPC
511  coreprof(1)%vtor%transp_coef%diff = transport1%DIFF_VTOR
512  coreprof(1)%vtor%transp_coef%vconv = transport1%VCONV_VTOR
513 !DPC.end
514 
515 !DPC set fluxes to 0 (otherwise we write garbage)
516  coreprof(1)%ni%flux%flux_dv = 0.0_r8
517  coreprof(1)%ne%flux%flux_dv = 0.0_r8
518  coreprof(1)%ti%flux%flux_dv = 0.0_r8
519  coreprof(1)%te%flux%flux_dv = 0.0_r8
520  coreprof(1)%vtor%flux%flux_dv = 0.0_r8
521 !DPC.end
522 
523 ! +++ Convert transport:
524  coretransp(1)%values(1)%sigma = transport1%SIGMA
525 
526  coretransp(1)%values(1)%te_transp%diff_eff = transport1%DIFF_TE
527  coretransp(1)%values(1)%te_transp%vconv_eff = transport1%VCONV_TE
528 
529  coretransp(1)%values(1)%ni_transp%diff_eff = transport1%DIFF_NI
530  coretransp(1)%values(1)%ni_transp%vconv_eff = transport1%VCONV_NI
531 
532  coretransp(1)%values(1)%ti_transp%diff_eff = transport1%DIFF_TI
533  coretransp(1)%values(1)%ti_transp%vconv_eff = transport1%VCONV_TI
534 
535  coretransp(1)%values(1)%vtor_transp%diff_eff = transport1%DIFF_VTOR
536  coretransp(1)%values(1)%vtor_transp%vconv_eff = transport1%VCONV_VTOR
537 
538 
539 
540 ! +++ Convert sources:
541  coresource(1)%VALUES(1)%j = sources1%CURR_EXP
542 
543  coresource(1)%VALUES(1)%qe%exp = sources1%QE_EXP * itm_ev !! DPC 2009-07-02
544  coresource(1)%VALUES(1)%qe%imp = sources1%QE_IMP
545 
546  coresource(1)%VALUES(1)%si%exp = sources1%SI_EXP
547  coresource(1)%VALUES(1)%si%imp = sources1%SI_IMP
548 
549  coresource(1)%VALUES(1)%qi%exp = sources1%QI_EXP * itm_ev !! DPC 2009-07-02
550  coresource(1)%VALUES(1)%qi%imp = sources1%QI_IMP
551 
552  coresource(1)%VALUES(1)%ui%exp = sources1%UI_EXP
553  coresource(1)%VALUES(1)%ui%imp = sources1%UI_IMP
554 
555 
556 
557  RETURN
558 
559 
560  END SUBROUTINE convert_local_to_cpo_types
561 
562 
563 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
564 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
565 
566 
567 
568 
569 
570 
571 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
580 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
581  SUBROUTINE analytics_full (TIME, GEOMETRY1, PROFILES1, TRANSPORT1, SOURCES1)
582 
583 
584 !-------------------------------------------------------!
585 ! This routine manufactures the solution for the !
586 ! set of transport equations describing the main !
587 ! plasma. It provides sources and transport !
588 ! coefficients for all transport equations. !
589 !-------------------------------------------------------!
590 ! Source: --- !
591 ! Developers: R.Stankiewicz, D.Kalupin !
592 ! Kontacts: Roman.Stankiewich@gmail.com !
593 ! D.Kalupin@fz-juelich.de !
594 ! !
595 ! Comments: --- !
596 ! !
597 !-------------------------------------------------------!
598 
599 
600  USE itm_types
601  USE ets_plasma
604 
605  IMPLICIT NONE
606 
607  REAL (R8) :: time !Time
608 
609 
610  INTEGER :: nrho, irho !number of radial points (input) and current radial knot
611  INTEGER :: nion, iion !number of ion species (input) and current ion type
612  INTEGER :: sigma_source
613  INTEGER :: ifail
614 
615 
616 ! +++ Local derived types:
617  TYPE (magnetic_geometry), INTENT(INOUT) :: geometry1 !contains all geometry quantities
618  TYPE (plasma_profiles), INTENT(INOUT) :: profiles1 !contains profiles of plasma parameters
619  TYPE (transport_coefficients), INTENT(INOUT) :: transport1 !contains profiles of trasport coefficients
620  TYPE (sources_and_sinks), INTENT(INOUT) :: sources1 !contains profiles of sources
621  TYPE (collisionality) :: collisions1 !contains all terms determined by plasma collisions
622 
623 
624 ! +++ Local variables:
625 ! Geometry:
626  REAL (R8) :: rho(profiles1%nrho), hrho(profiles1%nrho)
627  REAL (R8) :: fdia(profiles1%nrho), vpr(profiles1%nrho)
628  REAL (R8) :: g1(profiles1%nrho), g2(profiles1%nrho)
629  REAL (R8) :: g3(profiles1%nrho)
630 
631 ! Plasma composition:
632  REAL (R8) :: mion(profiles1%nion), zion(profiles1%nion)
633 
634 ! Plasma profiles:
635  REAL (R8) :: psi(profiles1%nrho)
636  REAL (R8) :: ne(profiles1%nrho)
637  REAL (R8) :: ni(profiles1%nrho,profiles1%nion)
638  REAL (R8) :: te(profiles1%nrho)
639  REAL (R8) :: ti(profiles1%nrho,profiles1%nion)
640  REAL (R8) :: vtor(profiles1%nrho,profiles1%nion)
641 
642 ! Transport coefficients:
643  REAL (R8) :: sigma(profiles1%nrho)
644  REAL (R8) :: diff_ni(profiles1%nrho,profiles1%nion)
645  REAL (R8) :: vconv_ni(profiles1%nrho,profiles1%nion)
646  REAL (R8) :: diff_ti(profiles1%nrho,profiles1%nion)
647  REAL (R8) :: vconv_ti(profiles1%nrho,profiles1%nion)
648  REAL (R8) :: diff_te(profiles1%nrho)
649  REAL (R8) :: vconv_te(profiles1%nrho)
650  REAL (R8) :: diff_vtor(profiles1%nrho,profiles1%nion)
651  REAL (R8) :: vconv_vtor(profiles1%nrho,profiles1%nion)
652 
653 ! Sources:
654  REAL (R8) :: vie(profiles1%nrho), vei(profiles1%nrho)
655  REAL (R8) :: qie(profiles1%nrho), qei(profiles1%nrho)
656  REAL (R8) :: vzi(profiles1%nrho), uzi(profiles1%nrho)
657  REAL (R8) :: qzi(profiles1%nrho), wzi(profiles1%nrho)
658 
659 ! Other quantities:
660  REAL (R8) :: dne(profiles1%nrho), dtne(profiles1%nrho)
661  REAL (R8) :: c1
662  REAL (R8) :: bt, r
663 
664 
665 ! +++ Functions:
666  REAL (R8) :: x, t
667  REAL (R8) :: aval1, aval2, aval3, aval4, aval5, aval6, aval7
668  REAL (R8) :: gamma(profiles1%nion,profiles1%nrho)
669  REAL (R8) :: dgamma(profiles1%nion,profiles1%nrho)
670  REAL (R8) :: gammae(profiles1%nrho)
671  REAL (R8) :: dgammae(profiles1%nrho)
672 
673 
674 ! +++ Set up dimensions:
675  nrho = SIZE(profiles1%NI)
676  nion = SIZE(profiles1%MION)
677  c1 = transport1%C1(1)
678  sigma_source = 0
679 
680 
681  CALL allocate_collisionality(nrho, nion, collisions1, ifail)
682 
683 
684 
685 
686 !=================================================================================
687 ! +++ MAGNETIC GEOMETRY:
688 
689 !--------------------------------------------------------------------------------!
690 ! This part sets up the geometry used by analytical solution. !
691 !--------------------------------------------------------------------------------!
692 
693  t = time
694  r = 1.0e0_r8
695  bt = abt(t)
696 
697 
698  DO irho=1,nrho
699  x = arho(irho,nrho)
700  rho(irho) = x
701  vpr(irho) = avpr(x,t)
702  g1(irho) = ag1(x,t)
703  g2(irho) = ag2(x,t)
704  g3(irho) = ag3(x,t)
705  fdia(irho) = afdia(x,t)
706 
707  geometry1%RHO(irho) = x
708  geometry1%VPR(irho) = vpr(irho)
709  geometry1%G1(irho) = g1(irho)
710  geometry1%G2(irho) = g2(irho)
711  geometry1%G3(irho) = g3(irho)
712  geometry1%FDIA(irho) = fdia(irho)
713  ENDDO
714 
715  geometry1%RGEO = r
716  geometry1%BGEO = bt
717 
718 
719 
720 
721 
722 
723 !=================================================================================
724 ! +++ PROFILES OF PLASMA PARAMETERS:
725 
726 !--------------------------------------------------------------------------------!
727 ! This part describes analytical profiles for all plasma parameters. !
728 !--------------------------------------------------------------------------------!
729 
730  t = time
731 
732  DO iion=1,nion
733  zion(iion) = profiles1%ZION(iion)
734  mion(iion) = profiles1%MION(iion)!*MP
735  ENDDO
736 
737 
738  DO irho=1,nrho
739  x = rho(irho)
740  te(irho) = ate(x,t)
741  psi(irho) = apsi(x,t)
742 
743  profiles1%TE(irho) = te(irho)
744  profiles1%PSI(irho) = psi(irho)
745  profiles1%NE(irho) = 0.e0_r8
746 
747  DO iion=1,nion
748  ti(irho,iion) = ati(iion,x,t)
749  ni(irho,iion) = ani(iion,x,t)
750  vtor(irho,iion) = avtor(iion,x,t)
751 
752  profiles1%NE(irho) = profiles1%NE(irho)+profiles1%ZION(iion)*ni(irho,iion)
753  profiles1%NI(irho,iion) = ni(irho,iion)
754  profiles1%TI(irho,iion) = ti(irho,iion)
755  profiles1%VTOR(irho,iion) = vtor(irho,iion)
756  ENDDO
757 
758  ENDDO
759 
760 
761 
762 
763 
764 
765 !=================================================================================
766 ! +++ CALCULATING COLLISIONAL TERMS:
767  CALL plasma_collisions(geometry1,profiles1,collisions1,ifail)
768 
769 
770 
771 
772 
773 
774 !=================================================================================
775 ! +++ TRANSPORT_COEFFICIENT
776 
777 !--------------------------------------------------------------------------------!
778 ! This part describes analytical profiles of transport coefficients. !
779 !--------------------------------------------------------------------------------!
780 
781  DO irho=1,nrho
782  x = rho(irho)
783 
784 
785  IF(sigma_source.EQ.1) THEN
786  sigma(irho) = asigma(x,t)
787 
788  transport1%SIGMA(irho) = sigma(irho)
789  ENDIF
790 
791 
792  IF(sigma_source.EQ.0) THEN
793  sigma(irho) = collisions1%SIGMA(irho)
794  transport1%SIGMA(irho) = sigma(irho)
795  ENDIF
796 
797  ENDDO
798 
799 
800 
801 
802  DO irho=1,nrho
803  x = rho(irho)
804  diff_te(irho) = akate(x,t)
805  vconv_te(irho) = avte(x,t)
806 
807  transport1%DIFF_TE(irho) = diff_te(irho)
808  transport1%VCONV_TE(irho) = vconv_te(irho)
809 
810 
811  DO iion=1,nion
812  diff_ni(irho,iion) = ad(iion,x,t)
813  vconv_ni(irho,iion) = av(iion,x,t)
814  diff_ti(irho,iion) = akappae(iion,x,t)
815  vconv_ti(irho,iion) = avti(iion,x,t)
816  diff_vtor(irho,iion) = advtor(iion,x,t)
817  vconv_vtor(irho,iion) = aconvtor(iion,x,t)
818 
819 
820  IF(c1.EQ.0.0e0_r8) THEN
821  transport1%DIFF_NI(irho,iion,1) = diff_ni(irho,iion)
822  transport1%VCONV_NI(irho,iion,1)= vconv_ni(irho,iion)
823  ELSE IF(c1.EQ.1.5e0_r8) THEN
824  transport1%DIFF_NI(irho,iion,2) = diff_ni(irho,iion)
825  transport1%VCONV_NI(irho,iion,2)= vconv_ni(irho,iion)
826  ELSE IF(c1.EQ.2.5e0_r8) THEN
827  transport1%DIFF_NI(irho,iion,3) = diff_ni(irho,iion)
828  transport1%VCONV_NI(irho,iion,3)= vconv_ni(irho,iion)
829  END IF
830  transport1%DIFF_TI(irho,iion) = diff_ti(irho,iion)
831  transport1%VCONV_TI(irho,iion) = vconv_ti(irho,iion)
832  transport1%DIFF_VTOR(irho,iion) = diff_vtor(irho,iion)
833  transport1%VCONV_VTOR(irho,iion) = vconv_vtor(irho,iion)
834  ENDDO
835 
836  ENDDO
837 
838 
839 
840 !=================================================================================
841 ! +++ SOURCES_AND_SINKS
842 !--------------------------------------------------------------------------------!
843 ! CALCULATION OF EXPLICIT SOURCES FROM TRANSPORT EQUATION USING THE !
844 ! ANALATICAL FORMULAE FOR SOLUTIONS AND TRANSPORT COEFFICIENTS !
845 !--------------------------------------------------------------------------------!
846 
847 
848 
849 
850 
851 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
852 ! +++ CURRENT TRANSPORT:
853  t = time
854 
855  DO irho=2,nrho
856  x = rho(irho)
857 
858  aval1 = sigma(irho)*dtapsi(x,t)
859 
860  aval2 = -x*dtabt(t)/2.0e0_r8/abt(t)*dapsi(x,t)*sigma(irho)
861  aval3 = fdia(irho)**2/itm_mu0/abt(t)/x
862 
863  aval4 = vpr(irho)/(4.0e0_r8*itm_pi**2)*g3(irho)/fdia(irho)
864  aval5 = davpr(x,t)/(4.0e0_r8*itm_pi**2)*g3(irho)/fdia(irho)
865  aval5 = aval5+avpr(x,t)/(4.0e0_r8*itm_pi**2)*dag3(x,t)/fdia(irho)
866  aval5 = aval5-avpr(x,t)/(4.0e0_r8*itm_pi**2)*g3(irho)/fdia(irho)**2*dafdia(x,t)
867  aval1 = aval1+aval2-aval3*aval4*ddapsi(x,t)-aval3*aval5*dapsi(x,t)
868 
869 
870  aval1 = aval1*(2.e0_r8*itm_pi*x)/vpr(irho)+acurr_imp(x,t)*psi(irho)
871 
872 
873  sources1%CURR_EXP(irho)= -aval1
874  sources1%CURR_IMP(irho)= acurr_imp(x,t)
875 
876  ENDDO
877 
878 
879 
880  sources1%CURR_EXP(1) = sources1%CURR_EXP(2)
881  sources1%CURR_IMP(1) = acurr_imp(0.e0_r8,t)
882 
883 
884 
885 ! +++ BOUNDARY CONDITION
886  IF(profiles1%PSI_BND_TYPE.EQ.0) THEN
887  DO irho=2,nrho
888  profiles1%QSF(irho) = 2.e0_r8*itm_pi*rho(irho)*abt(t)/dapsi(rho(irho),t)
889  ENDDO
890 
891  profiles1%QSF(1) = 2.e0_r8*itm_pi*abt(t)/ddapsi(rho(1),t)
892  ENDIF
893 
894 
895  IF(profiles1%PSI_BND_TYPE.EQ.1) THEN
896  profiles1%PSI_BND(1) = psi(nrho)
897  ENDIF
898 
899  IF(profiles1%PSI_BND_TYPE.EQ.2) THEN
900  profiles1%PSI_BND(1) = vpr(nrho)*g3(nrho)*dapsi(rho(nrho),t)/4.e0_r8/itm_pi**2/itm_mu0
901  ENDIF
902 
903  IF(profiles1%PSI_BND_TYPE.EQ.3) THEN
904  profiles1%PSI_BND(1) = dtapsi(rho(nrho),t)
905  ENDIF
906 
907 
908 
909 
910 
911 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
912 ! +++ PARTICLE TRANSPORT
913 
914  DO iion=1,nion
915  t = time
916  gamma(iion,1) = 0.0e0_r8 !DPC 2009-01-19
917  dgamma(iion,1) = 0.0e0_r8 !DPC 2009-01-19
918  DO irho=2,nrho
919  x = rho(irho)
920  aval1 = dtavpr(x,t)*ani(iion,x,t)+avpr(x,t)*dtani(iion,x,t)
921  aval2 = -dtabt(t)/2.e0_r8/abt(t)*(avpr(x,t)* ani(iion,x,t)+ &
922  x *(davpr(x,t)* ani(iion,x,t)+avpr(x,t)* dani(iion,x,t)))
923  gamma(iion,irho) = avpr(x,t)* ag1(x,t)*(-ad(iion,x,t)* dani(iion,x,t)+ani(iion,x,t)*av(iion,x,t))
924  dgamma(iion,irho) = avpr(x,t)* ag1(x,t)*(-ad(iion,x,t)* ddani(iion,x,t)+dani(iion,x,t)*av(iion,x,t))
925 
926  dgamma(iion,irho) = dgamma(iion,irho)+ &
927  davpr(x,t)* ag1(x,t)*(-ad(iion,x,t)* dani(iion,x,t)+av(iion,x,t)*ani(iion,x,t))
928  dgamma(iion,irho) = dgamma(iion,irho)+ &
929  avpr(x,t)* dag1(x,t)*(-ad(iion,x,t)* dani(iion,x,t)+av(iion,x,t)*ani(iion,x,t))
930  dgamma(iion,irho) = dgamma(iion,irho)+ &
931  avpr(x,t)* ag1(x,t)*(-dad(iion,x,t)* dani(iion,x,t)+dav(iion,x,t)*ani(iion,x,t))
932 
933  aval1 = aval1+aval2+dgamma(iion,irho)+avpr(x,t)*ani(iion,x,t)*asi_imp(iion,x,t)
934  sources1%SI_EXP(irho,iion) = aval1/avpr(x,t)
935  sources1%SI_IMP(irho,iion) = asi_imp(iion,x,t)
936  ENDDO
937 
938 
939  sources1%SI_EXP(1,iion) = sources1%SI_EXP(2,iion)
940  sources1%SI_IMP(1,iion) = asi_imp(iion,0.e0_r8,t)
941 
942 
943 ! +++ BOUNDARY CONDITION
944  IF (profiles1%NI_BND_TYPE(iion).EQ.1) THEN
945  profiles1%NI_BND(1,iion) = ani(iion,rho(nrho),t)
946  ENDIF
947 
948  IF (profiles1%NI_BND_TYPE(iion).EQ.2) THEN
949  profiles1%NI_BND(1,iion) = -dani(iion,rho(nrho),t)
950  ENDIF
951 
952  IF (profiles1%NI_BND_TYPE(iion).EQ.3) THEN
953  profiles1%NI_BND(1,iion) = -ani(iion,rho(nrho),t)/dani(iion,rho(nrho),t)
954  ENDIF
955 
956  IF (profiles1%NI_BND_TYPE(iion).EQ.4) THEN
957  profiles1%NI_BND(1,iion) = gamma(iion,nrho)
958  ENDIF
959 
960  ENDDO
961 
962 
963 
964 
965 
966 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
967 ! +++ ELECTRON DENSISTY AND FLUXES
968  DO irho =1,nrho
969  x = rho(irho)
970  ne(irho) = 0.e0_r8
971  dne(irho) = 0.e0_r8
972  dtne(irho) = 0.e0_r8
973  gammae(irho) = 0.e0_r8
974  dgammae(irho) = 0.e0_r8
975 
976  DO iion=1,nion
977  ne(irho) = ne(irho) + zion(iion)*ni(irho,iion)
978  dtne(irho) = dtne(irho) + zion(iion)*dtani(iion,x,t)
979  dne(irho) = dne(irho) + zion(iion)*dani(iion,x,t)
980  gammae(irho) = gammae(irho) + gamma(iion,irho)*zion(iion)
981  dgammae(irho) = dgammae(irho) + dgamma(iion,irho)*zion(iion)
982  ENDDO
983 
984  ENDDO
985 
986 
987 
988 
989 
990 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
991 ! +++ ION ENERGY TRANSPORT
992  DO iion=1,nion
993  t = time
994 
995  DO irho=2,nrho
996  x = rho(irho)
997  aval1 = dtavpr(x,t)*avpr(x,t)**(2.e0_r8/3.e0_r8)*ani(iion,x,t)*ati(iion,x,t)*5.e0_r8/3.e0_r8
998  aval1 = aval1+avpr(x,t)**(5.e0_r8/3.e0_r8)*dtani(iion,x,t)*ati(iion,x,t)
999  aval1 = aval1+avpr(x,t)**(5.e0_r8/3.e0_r8)*ani(iion,x,t)*dtati(iion,x,t)
1000 
1001  aval2 = -dtabt(t)/2.e0_r8/abt(t)
1002 
1003  aval3 = ani(iion,x,t)*ati(iion,x,t)*avpr(x,t)**(5.e0_r8/3.e0_r8)
1004  aval3 = aval3+ x*dani(iion,x,t)*ati(iion,x,t)*avpr(x,t)**(5.e0_r8/3.e0_r8)
1005  aval3 = aval3+ x*ani(iion,x,t)*dati(iion,x,t)*avpr(x,t)**(5.e0_r8/3.e0_r8)
1006  aval3 = aval3+ x*ani(iion,x,t)*ati(iion,x,t)*avpr(x,t)**(2.e0_r8/3.e0_r8)*davpr(x,t)*5.e0_r8/3.e0_r8
1007 
1008  aval1 = aval1+aval2*aval3
1009  aval1 = aval1*3.e0_r8/2.e0_r8
1010 
1011  aval2 = avpr(x,t)*ag1(x,t)*ani(iion,x,t)
1012 
1013  aval3 = -akappae(iion,x,t)*dati(iion,x,t)+ati(iion,x,t)*avti(iion,x,t)
1014 
1015  aval4 = davpr(x,t)*ag1(x,t)*ani(iion,x,t)+ &
1016  avpr(x,t)*dag1(x,t)*ani(iion,x,t)+avpr(x,t)*ag1(x,t)*dani(iion,x,t)
1017 
1018  aval5 = -dakappae(iion,x,t)*dati(iion,x,t)+ati(iion,x,t)*davti(iion,x,t)
1019  aval5 = aval5-akappae(iion,x,t)*ddati(iion,x,t)+dati(iion,x,t)*avti(iion,x,t)
1020 
1021  aval6 = aval2*aval5+aval4*aval3
1022 
1023  aval1 = aval1+aval6*avpr(x,t)**(2.e0_r8/3.e0_r8)
1024 
1025 
1026  aval2 = dati(iion,x,t)*gamma(iion,irho)+ati(iion,x,t)*dgamma(iion,irho)
1027 
1028  aval1 = aval1+c1*aval2*avpr(x,t)**(2.e0_r8/3.e0_r8)
1029  aval1 = aval1/avpr(x,t)**(5.e0_r8/3.e0_r8)
1030 
1031  vei(irho) = collisions1%VEI(irho,iion)
1032  qei(irho) = collisions1%QEI(irho,iion)
1033  vzi(irho) = collisions1%VZI(irho,iion)
1034  qzi(irho) = collisions1%QZI(irho,iion)
1035 
1036  aval2 = (aqi_imp(iion,x,t) + vei(irho) + vzi(irho) )
1037 
1038  sources1%QI_EXP(irho,iion) = aval1- qei(irho) - qzi(irho) +aval2*ati(iion,x,t)
1039  sources1%QI_IMP(irho,iion) = aqi_imp(iion,x,t)
1040  ENDDO
1041 
1042  sources1%QI_EXP(1,iion) = sources1%QI_EXP(2,iion)
1043  sources1%QI_IMP(1,iion) = aqi_imp(iion,0.e0_r8,t)
1044  ENDDO
1045 
1046 
1047 
1048 
1049 ! +++ BOUNDARY CONDITION
1050  DO iion=1,nion
1051  IF (profiles1%TI_BND_TYPE(iion).EQ.1)THEN
1052  profiles1%TI_BND(1,iion) = ati(iion,rho(nrho),t)
1053  ENDIF
1054 
1055  IF (profiles1%TI_BND_TYPE(iion).EQ.2)THEN
1056  profiles1%TI_BND(1,iion) = -dati(iion,rho(nrho),t)
1057  ENDIF
1058 
1059  IF (profiles1%TI_BND_TYPE(iion).EQ.3)THEN
1060  profiles1%TI_BND(1,iion) = -ati(iion,rho(nrho),t)/dati(iion,rho(nrho),t)
1061  ENDIF
1062 
1063  IF (profiles1%TI_BND_TYPE(iion).EQ.4)THEN
1064  x = rho(nrho)
1065  aval2 = avpr(x,t)*ag1(x,t)*ani(iion,x,t)
1066  aval3 = (-akappae(iion,x,t)*dati(iion,x,t)+ati(iion,x,t)*avti(iion,x,t))* &
1067  aval2+c1*ati(iion,x,t)*gamma(iion,nrho)
1068  profiles1%TI_BND(1,iion) = aval3
1069  ENDIF
1070 
1071  ENDDO
1072 
1073 
1074 
1075 
1076 
1077 
1078 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
1079 ! +++ ELECTRON ENERGY
1080 
1081  t = time
1082 
1083  DO irho=2,nrho
1084  x = rho(irho)
1085  aval1 = dtavpr(x,t)*avpr(x,t)**(2.e0_r8/3.e0_r8)*ne(irho)*ate(x,t)*5.e0_r8/3.e0_r8
1086  aval1 = aval1+avpr(x,t)**(5.e0_r8/3.e0_r8)*dtne(irho)*ate(x,t)
1087  aval1 = aval1+avpr(x,t)**(5.e0_r8/3.e0_r8)*ne(irho)*dtate(x,t)
1088  aval2 = -dtabt(t)/2.e0_r8/abt(t)
1089  aval3 = ne(irho)*ate(x,t)*avpr(x,t)**(5.e0_r8/3.e0_r8)
1090  aval3 = aval3+ x*dne(irho)*ate(x,t)*avpr(x,t)**(5.e0_r8/3.e0_r8)
1091  aval3 = aval3+ x*ne(irho)*date(x,t)*avpr(x,t)**(5.e0_r8/3.e0_r8)
1092  aval3 = aval3+ x*ne(irho)*ate(x,t)*avpr(x,t)**(2.e0_r8/3.e0_r8)*davpr(x,t)*5.e0_r8/3.e0_r8
1093  aval1 = aval1+aval2*aval3
1094  aval1 = aval1*3.e0_r8/2.e0_r8
1095  aval2 = avpr(x,t)*ag1(x,t)*ne(irho)
1096  aval3 = -akate(x,t)*date(x,t)+ate(x,t)*avte(x,t)
1097  aval4 = davpr(x,t)*ag1(x,t)*ne(irho)+avpr(x,t)*dag1(x,t)*ne(irho)+avpr(x,t)*ag1(x,t)*dne(irho)
1098  aval5 = -dakate(x,t)*date(x,t)+ate(x,t)*davte(x,t)
1099  aval5 = aval5-akate(x,t)*ddate(x,t)+date(x,t)*avte(x,t)
1100  aval6 = aval2*aval5+aval4*aval3
1101  aval1 = aval1+aval6*avpr(x,t)**(2.e0_r8/3.e0_r8)
1102  c1 = 0.e0_r8
1103  aval2 = c1*(avpr(x,t)**(2.e0_r8/3.e0_r8)*(date(x,t)*gammae(irho)+ ate(x,t)*dgammae(irho)))
1104  aval1 = aval1+aval2
1105  aval1 = aval1/avpr(x,t)**(5.e0_r8/3.e0_r8)
1106  vie(irho) = collisions1%VIE(irho)
1107  qie(irho) = collisions1%QIE(irho)
1108  aval2 = (aqe_imp(x,t) + vie(irho) )
1109  sources1%QE_EXP(irho) = aval1 +aval2*ate(x,t)-qie(irho)
1110  sources1%QE_IMP(irho) = aqe_imp(x,t)
1111  ENDDO
1112 
1113  sources1%QE_EXP(1) = sources1%QE_EXP(2)
1114  sources1%QE_IMP(1) = sources1%QE_IMP(2)
1115 
1116 
1117 
1118 ! +++ BOUNDARY CONDITION
1119  IF (profiles1%TE_BND_TYPE.EQ.1)THEN
1120  profiles1%TE_BND(1) = te(nrho)
1121  ENDIF
1122 
1123  IF (profiles1%TE_BND_TYPE.EQ.2)THEN
1124  profiles1%TE_BND(1) = -date(rho(nrho),t)
1125  ENDIF
1126 
1127  IF (profiles1%TE_BND_TYPE.EQ.3)THEN
1128  profiles1%TE_BND(1) = -te(nrho)/date(rho(nrho),t)
1129  ENDIF
1130 
1131  IF (profiles1%TE_BND_TYPE.EQ.4)THEN
1132  x = rho(nrho)
1133  aval2 = avpr(x,t)*ag1(x,t)*ne(nrho)
1134  aval3 = -akate(x,t)*date(x,t)+ate(rho(nrho),t)*avte(x,t)
1135  profiles1%TE_BND(1) = aval2*aval3+c1*ate(rho(nrho),t)*gammae(nrho)
1136  ENDIF
1137 
1138 
1139 
1140 
1141 
1142 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
1143 ! +++ ROTATION TRANSPORT
1144  DO iion=1,nion
1145  t = time
1146 
1147  DO irho=2,nrho
1148  x = rho(irho)
1149  aval1 = dtavpr(x,t)*ani(iion,x,t)*avtor(iion,x,t)*ag2(x,t)
1150  aval1 = aval1+avpr(x,t)*dtani(iion,x,t)*avtor(iion,x,t)*ag2(x,t)
1151  aval1 = aval1+avpr(x,t)*ani(iion,x,t)*dtavtor(iion,x,t)*ag2(x,t)
1152  aval1 = aval1+avpr(x,t)*ani(iion,x,t)*avtor(iion,x,t)*dtag2(x,t)
1153 
1154  aval2 = -dtabt(t)/2.e0_r8/abt(t)
1155 
1156  aval3 = ani(iion,x,t)*avtor(iion,x,t)*avpr(x,t)*ag2(x,t)
1157  aval3 = aval3+ x*dani(iion,x,t)*avtor(iion,x,t)*avpr(x,t)*ag2(x,t)
1158  aval3 = aval3+ x*ani(iion,x,t)*davtor(iion,x,t)*avpr(x,t)*ag2(x,t)
1159  aval3 = aval3+ x*ani(iion,x,t)*avtor(iion,x,t)*davpr(x,t)*ag2(x,t)
1160  aval3 = aval3+ x*ani(iion,x,t)*avtor(iion,x,t)*davpr(x,t)*dag2(x,t)
1161 
1162  aval1 = aval1+aval2*aval3
1163 
1164  aval2 = avpr(x,t)*ag2(x,t)*ani(iion,x,t)*ag1(x,t)
1165 
1166  aval3 = -advtor(iion,x,t)*davtor(iion,x,t)+avtor(iion,x,t)*aconvtor(iion,x,t)
1167 
1168  aval4 = davpr(x,t)*ag2(x,t)*ani(iion,x,t)+ &
1169  avpr(x,t)*dag2(x,t)*ani(iion,x,t)+avpr(x,t)*ag2(x,t)*dani(iion,x,t)
1170  aval4 = aval4*ag1(x,t)+avpr(x,t)*ag2(x,t)*ani(iion,x,t)*dag1(x,t)
1171 
1172  aval5 = -dadvtor(iion,x,t)*davtor(iion,x,t)+avtor(iion,x,t)*daconvtor(iion,x,t)
1173  aval5 = aval5-advtor(iion,x,t)*ddavtor(iion,x,t)+davtor(iion,x,t)*aconvtor(iion,x,t)
1174 
1175  aval6 = aval2*aval5+aval4*aval3 +dag2(x,t)*gamma(iion,irho)*avtor(iion,x,t)!*AVAL4
1176  aval6 = aval6+(ag2(x,t)*dgamma(iion,irho)*avtor(iion,x,t)+ &
1177  ag2(x,t)*gamma(iion,irho)*davtor(iion,x,t))!*AVAL4
1178 
1179  aval1 = aval1+aval6
1180  aval1 = aval1*profiles1%MION(iion)/avpr(x,t)
1181 
1182 
1183 
1184  uzi(irho) = collisions1%UZI(irho,iion)
1185  wzi(irho) = collisions1%WZI(irho,iion)
1186  aval2 = aui_imp(iion,x,t) + wzi(irho)
1187 
1188  sources1%UI_EXP(irho,iion) = aval1 +aval2*avtor(iion,x,t) - uzi(irho)
1189  sources1%UI_IMP(irho,iion) = aui_imp(iion,x,t)
1190  ENDDO
1191 
1192  sources1%UI_EXP(1,iion) = sources1%UI_EXP(2,iion)
1193  sources1%UI_EXP(1,iion) = aui_imp(iion,0.e0_r8,t)
1194  sources1%UI_IMP(1,iion) = aui_imp(iion,0.e0_r8,t)
1195  ENDDO
1196 
1197 
1198 
1199 
1200 ! +++ BOUNDARY CONDITION
1201  DO iion=1,nion
1202  IF (profiles1%VTOR_BND_TYPE(iion).EQ.1)THEN
1203  profiles1%VTOR_BND(1,iion) = avtor(iion,rho(nrho),t)
1204  ENDIF
1205 
1206  IF (profiles1%VTOR_BND_TYPE(iion).EQ.2)THEN
1207  profiles1%VTOR_BND(1,iion) = -davtor(iion,rho(nrho),t)
1208  ENDIF
1209 
1210  IF (profiles1%VTOR_BND_TYPE(iion).EQ.3)THEN
1211  profiles1%VTOR_BND(1,iion) = -avtor(iion,rho(nrho),t)/davtor(iion,rho(nrho),t)
1212  ENDIF
1213 
1214  IF (profiles1%VTOR_BND_TYPE(iion).EQ.4)THEN
1215  x = rho(nrho)
1216  aval2 = avpr(x,t)*ag2(x,t)*ani(iion,x,t)*ag1(x,t)
1217  aval3 = -advtor(iion,x,t)*davtor(iion,x,t)+avtor(iion,rho(nrho),t)*aconvtor(iion,x,t)
1218  profiles1%VTOR_BND(1,iion) = aval2*aval3 +avtor(iion,rho(nrho),t)*ag2(x,t)*gamma(iion,nrho)
1219  profiles1%VTOR_BND(1,iion) = profiles1%VTOR_BND(1,iion)*mion(iion)
1220  ENDIF
1221 
1222  ENDDO
1223 
1224 
1225 
1226  CALL deallocate_collisionality(collisions1, ifail)
1227 
1228 
1229  RETURN
1230 
1231 
1232 
1233  END SUBROUTINE analytics_full
1234 
1235 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
1236 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
1237 
1238 
1239 
1240 
1241 
1242 END MODULE analytics1
subroutine copy_codeparam(CODEPARAM_IN, CODEPARAM_OUT)
COPY CODEPARAM.
REAL(R8) function avpr(X, T)
Module provides routines for testing.
subroutine process_xml(SOLVER_TYPE, SIGMA_SOURCE, TAU, AMIX, CONVREC, NRHO, NION, NIMP, NZIMP, NTIME, NSOL, PSI_BND_TYPE, NI_BND_TYPE, TI_BND_TYPE, TE_BND_TYPE, VTOR_BND_TYPE, shot_no, run_no, codeparam, database_format)
process the xml version of the input file from codeparam
REAL(R8) function dakate(X, T)
Module provides routines for copying parts of CPOs (COREPROF and EQUILIBRIUM)
Definition: copy_cpo_ets.f90:8
REAL(R8) function dadvtor(IION, X, T)
subroutine allocate_magnetic_geometry(NRHO, GEOMETRY, ifail)
Definition: ets_plasma.f90:502
REAL(R8) function date(X, T)
subroutine allocate_coreimpur_cpo(NSLICE, NRHO, NNUCL, NION, NIMP, NZIMP, NNEUT, NTYPE, NCOMP, COREIMPUR)
This routine allocates COREIMPUR CPO.
REAL(R8) function dtati(IION, X, T)
REAL(R8) function ddati(IION, X, T)
REAL(R8) function davte(X, T)
subroutine convert_cpo_to_local_types(COREPROF, PROFILES1, TRANSPORT1)
This routine converts CPOs into the local derived types.
Definition: analytics1.f90:331
REAL(R8) function aui_imp(IION, X, T)
REAL(R8) function ag3(X, T)
REAL(R8) function dtag2(X, T)
REAL(R8) function dakappae(IION, X, T)
REAL(R8) function dad(IION, X, T)
REAL(R8) function advtor(IION, X, T)
REAL(R8) function ag1(X, T)
REAL(R8) function abt(T)
subroutine allocate_coreprof_cpo(NSLICE, NRHO, NNUCL, NION, NIMP, NZIMP, NNEUT, NTYPE, NCOMP, COREPROF)
This routine allocates COREPROF CPO.
REAL(R8) function apsi(X, T)
subroutine allocate_global_param(GLOBAL, ifail)
Definition: ets_plasma.f90:630
subroutine allocate_plasma_profiles(NRHO, NION, PROFILES, ifail)
Definition: ets_plasma.f90:717
subroutine set_cpo(NRHO, NION, NIMP, NZIMP, NTIME, NSOL, PSI_BND_TYPE, NI_BND_TYPE, TI_BND_TYPE, TE_BND_TYPE, VTOR_BND_TYPE, EQUILIBRIUM, COREPROF, CORETRANSP)
REAL(R8) function aconvtor(IION, X, T)
subroutine convert_local_to_cpo_types(GEOMETRY1, PROFILES1, TRANSPORT1, SOURCES1, EQUILIBRIUM, COREPROF, CORETRANSP, CORESOURCE)
This routine converts local into the CPOs derived types.
Definition: analytics1.f90:409
REAL(R8) function asigma(X, T)
REAL(R8) function av(IION, X, T)
REAL(R8) function davti(IION, X, T)
REAL(R8) function dtate(X, T)
subroutine allocate_sources_and_sinks(NRHO, NION, SOURCES, ifail)
Allocate profiles of sources needed by the transport solver.
REAL(R8) function ate(X, T)
REAL(R8) function arho(I, NRHO)
REAL(R8) function dav(IION, X, T)
REAL(R8) function dtapsi(X, T)
REAL(R8) function akate(X, T)
REAL(R8) function ddani(IION, X, T)
subroutine analytics_full(TIME, GEOMETRY1, PROFILES1, TRANSPORT1, SOURCES1)
This routine manufactures the solution for the set of transport equations describing the main plasma...
Definition: analytics1.f90:581
REAL(R8) function ad(IION, X, T)
This routine calculates the collision frquencies and various exchange terms determined by collisions...
REAL(R8) function dapsi(X, T)
REAL(R8) function daconvtor(IION, X, T)
REAL(R8) function ati(IION, X, T)
This module contains routines for allocation/deallocation if CPOs used in ETS.
REAL(R8) function dtavtor(IION, X, T)
subroutine allocate_coresource_cpo(NSLICE, NRHO, NNUCL, NION, NIMP, NZIMP, NNEUT, NTYPE, NCOMP, CORESOURCE)
This routine allocates CORESOURCE CPO.
REAL(R8) function avtor(IION, X, T)
subroutine allocate_equilibrium_cpo(NSLICE, NPSI, NDIM1, NDIM2, NPOINTS, EQUILIBRIUM)
This routine allocates EQUILIBRIUM CPO.
real(r8) function fdia(psi_n)
REAL(R8) function dtavpr(X, T)
subroutine allocate_transport_coefficients(NRHO, NION, TRANSPORT, ifail)
Allocate profiles of transport coefficients needed by the transport solver.
subroutine deallocate_magnetic_geometry(GEOMETRY, ifail)
Definition: ets_plasma.f90:575
subroutine deallocate_transport_coefficients(TRANSPORT, ifail)
Deallocate plasma profiles needed by the transport solver.
REAL(R8) function dafdia(X, T)
REAL(R8) function dag3(X, T)
REAL(R8) function ddate(X, T)
REAL(R8) function avti(IION, X, T)
REAL(R8) function ag2(X, T)
REAL(R8) function davtor(IION, X, T)
REAL(R8) function akappae(IION, X, T)
subroutine plasma_collisions(GEOMETRY, PROFILES, COLLISIONS, ifail)
REAL(R8) function ddavtor(IION, X, T)
REAL(R8) function acurr_imp(X, T)
REAL(R8) function avte(X, T)
REAL(R8) function afdia(X, T)
subroutine deallocate_sources_and_sinks(SOURCES, ifail)
Deallocate plasma profiles needed by the transport solver.
subroutine deallocate_collisionality(COLLISIONS, ifail)
Deallocate plasma profiles needed by the transport solver.
REAL(R8) function dag2(X, T)
subroutine allocate_collisionality(NRHO, NION, COLLISIONS, ifail)
Allocate profiles of sources needed by the transport solver???
REAL(R8) function asi_imp(IION, X, T)
The module declares types of variables used in ETS (transport code)
Definition: ets_plasma.f90:8
subroutine dati(D, T)
Definition: ppplib.f:4383
REAL(R8) function dtabt(T)
Analytical functions for the calculation of the &quot;analytic&quot; solution.
REAL(R8) function ddapsi(X, T)
subroutine analytical_plasma(TIME, COREPROF_in, EQUILIBRIUM, COREPROF_ANALYTIC, CORETRANSP, CORESOURCE, COREIMPUR, code_parameters)
This routine manufactures the solution for the set of transport equations describing the main plasma...
Definition: analytics1.f90:25
subroutine deallocate_plasma_profiles(PROFILES, ifail)
Deallocate plasma profiles needed by the transport solver.
Definition: ets_plasma.f90:992
REAL(R8) function dag1(X, T)
REAL(R8) function davpr(X, T)
REAL(R8) function aqe_imp(X, T)
REAL(R8) function dani(IION, X, T)
REAL(R8) function dtani(IION, X, T)
subroutine allocate_coretransp_cpo(NSLICE, NRHO, NNUCL, NION, NIMP, NZIMP, NNEUT, NTYPE, NCOMP, CORETRANSP)
This routine allocates CORETRANSP CPO.
REAL(R8) function aqi_imp(IION, X, T)
REAL(R8) function ani(IION, X, T)
Module for manufacture of a test case for the ETS.
Definition: analytics1.f90:8