ETS  \$Id: Doxyfile 2162 2020-02-26 14:16:09Z g2dpc $
 All Classes Files Functions Variables Pages
equilibrium_input.f90
Go to the documentation of this file.
1 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
12 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
14 
15 CONTAINS
16 
17 
18 
19 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
20 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
21 
22  SUBROUTINE equil_input &
23  (coreprof_in, toroidfield_in, equilibrium_in, &
24  equilibrium_out)
25 
26 !-------------------------------------------------------!
27 ! This routine provides consistent input for the !
28 ! equilibrium solver (EQUILIBRIUM_OUT) from the !
29 ! information saved in COREPROF_IN and !
30 ! EQUILIBRIUM_IN. !
31 !-------------------------------------------------------!
32 ! Source: --- !
33 ! Developers: D.Kalupin !
34 ! Kontacts: D.Kalupin@fz-juelich.de !
35 ! !
36 ! Comments: created for V&V between ETS and !
37 ! ASTRA !
38 ! !
39 !-------------------------------------------------------!
40 
41 
42  USE euitm_schemas
43  USE itm_constants
44  USE copy_structures
45  USE deallocate_structures
46  USE ets_plasma
47 
48  IMPLICIT NONE
49 
50 ! +++ CPO derived types:
51  TYPE (type_equilibrium), POINTER :: equilibrium_in(:) !input CPO with geometry quantities from previous time
52  TYPE (type_equilibrium), POINTER :: equilibrium_out(:) !output CPO with geometry quantities from previous iteration
53  TYPE (type_toroidfield), POINTER :: toroidfield_in(:) !toroidal field, major radius, total current
54  TYPE (type_coreprof), POINTER :: coreprof_in(:) !input CPO with plasma profiles
55 
56 ! +++ Local derived types:
57  TYPE (diagnostic) :: diag !diagnostic output
58 
59 
60 ! +++ Dimensions:
61  INTEGER :: nrho, irho !number and index of radial points
62  INTEGER :: nion, iion !number and index of ion components
63  INTEGER :: neq, ieq
64 
65  INTEGER :: ieq1, ieq2
66  REAL (R8) :: sign
67  REAL (R8), ALLOCATABLE :: pressure(:)
68  REAL (R8), DIMENSION (:), POINTER :: rho_tor, psi, rho_vol
69 
70  REAL (R8), ALLOCATABLE :: j_edge_mask(:)
71  REAL (R8) :: j_edge_pos, j_edge_pow
72  CHARACTER*3 :: char1, char2
73  CHARACTER*1000 :: output_diag
74  INTEGER :: nline, i
75  real (R8),allocatable,dimension(:):: psi_tmp,q_tmp,jtot_tmp,jphi_tmp
76  integer :: noll_cond
77 
78 
79 ! +++ Nullify diagnostic output:
80  diag%ERROR_MESSAGE = " "
81  diag%IERR = 0
82 
83 
84 ! +++ Set dimensions and allocate parameters:
85  nrho = SIZE (coreprof_in(1)%rho_tor, dim=1)
86  nion = SIZE (coreprof_in(1)%ni%value, dim=2)
87  neq = SIZE (equilibrium_in(1)%profiles_1d%psi, dim=1)
88 
89 
90 
91 ! +++ Fill the output EQUILIBRIUM CPO:
92  IF (ASSOCIATED(equilibrium_out)) &
93  CALL deallocate_cpo(equilibrium_out)
94  ALLOCATE (equilibrium_out(1))
95  CALL copy_cpo(equilibrium_in(1), equilibrium_out(1))
96 
97  IF (.NOT.ASSOCIATED(equilibrium_out(1)%profiles_1d%rho_tor)) ALLOCATE(equilibrium_out(1)%profiles_1d%rho_tor(neq))
98  IF (.NOT.ASSOCIATED(equilibrium_out(1)%profiles_1d%rho_vol)) ALLOCATE(equilibrium_out(1)%profiles_1d%rho_vol(neq))
99  IF (.NOT.ASSOCIATED(equilibrium_out(1)%profiles_1d%psi)) ALLOCATE(equilibrium_out(1)%profiles_1d%psi(neq))
100  IF (.NOT.ASSOCIATED(equilibrium_out(1)%profiles_1d%q)) ALLOCATE(equilibrium_out(1)%profiles_1d%q(neq))
101  IF (.NOT.ASSOCIATED(equilibrium_out(1)%profiles_1d%jphi)) ALLOCATE(equilibrium_out(1)%profiles_1d%jphi(neq))
102  IF (.NOT.ASSOCIATED(equilibrium_out(1)%profiles_1d%jparallel)) ALLOCATE(equilibrium_out(1)%profiles_1d%jparallel(neq))
103  IF (.NOT.ASSOCIATED(equilibrium_out(1)%profiles_1d%pressure)) ALLOCATE(equilibrium_out(1)%profiles_1d%pressure(neq))
104  IF (.NOT.ASSOCIATED(equilibrium_out(1)%profiles_1d%phi)) ALLOCATE(equilibrium_out(1)%profiles_1d%phi(neq))
105 
106  ALLOCATE (pressure(nrho))
107 
108 ! dpc
109  ALLOCATE(j_edge_mask(neq))
110  j_edge_pos = 0.02_r8
111  j_edge_pow = 8.0_r8
112 ! write(*,*) ' j_edge_pos = ', j_edge_pos
113 ! write(*,*) ' j_edge_pow = ', j_edge_pow
114  j_edge_mask = exp(-(j_edge_pos/(1.0_r8+1.0e-12_r8-equilibrium_in(1)%profiles_1d%rho_tor/equilibrium_in(1)%profiles_1d%rho_tor(neq)))**j_edge_pow)
115 ! cpd
116 
117 
118 ! +++ Fill the output EQUILIBRIUMvector%profiles_1d:
119  equilibrium_out(1)%profiles_1d%rho_tor = equilibrium_in(1)%profiles_1d%rho_tor
120 
121  IF (coreprof_in(1)%rho_tor(nrho).LT.equilibrium_in(1)%profiles_1d%rho_tor(neq)) &
122  coreprof_in(1)%rho_tor(nrho) = equilibrium_in(1)%profiles_1d%rho_tor(neq)
123 
124  IF ( ASSOCIATED(equilibrium_in(1)%profiles_1d%rho_vol )) THEN
125  equilibrium_out(1)%profiles_1d%rho_vol = equilibrium_in(1)%profiles_1d%rho_vol
126  ELSE IF (ASSOCIATED( equilibrium_in(1)%profiles_1d%volume)) THEN
127  equilibrium_out(1)%profiles_1d%rho_vol = sqrt(equilibrium_in(1)%profiles_1d%volume &
128  /equilibrium_in(1)%profiles_1d%volume(neq))
129  diag%ERROR_MESSAGE = trim(adjustl(diag%ERROR_MESSAGE))//"No rho_vol supplied; calculated from volume "
130  diag%IERR = diag%IERR + 1
131  ELSE
132  equilibrium_out(1)%profiles_1d%rho_vol = itm_invalid_float
133  diag%ERROR_MESSAGE = trim(adjustl(diag%ERROR_MESSAGE))//"No rho_vol supplied; set to itm_unvalid "
134  diag%IERR = diag%IERR + 1
135  END IF
136 
137  noll_cond=1
138  IF (ASSOCIATED (coreprof_in(1)%psi%value)) THEN
139  noll_cond=1
140  if (any(coreprof_in(1)%psi%value.ne. 0.0)) noll_cond=0
141  end if
142 
143  IF (noll_cond.eq.0) THEN
144  CALL l3interp(coreprof_in(1)%psi%value, coreprof_in(1)%rho_tor, nrho, &
145  equilibrium_out(1)%profiles_1d%psi, equilibrium_out(1)%profiles_1d%rho_tor, neq)
146  else if (associated(equilibrium_in(1)%profiles_1d%psi)) then
147  allocate(psi_tmp(nrho))
148  call l3interp(equilibrium_in(1)%profiles_1d%psi,equilibrium_in(1)%profiles_1d%rho_tor,neq, &
149  psi_tmp,coreprof_in(1)%rho_tor,nrho)
150  call l3interp(psi_tmp,coreprof_in(1)%rho_tor,nrho, &
151  equilibrium_out(1)%profiles_1d%psi,equilibrium_out(1)%profiles_1d%rho_tor,neq)
152  ELSE
153  equilibrium_out(1)%profiles_1d%psi = 0.0_r8
154  diag%ERROR_MESSAGE = trim(adjustl(diag%ERROR_MESSAGE))//"No Psi supplied; "
155  diag%IERR = diag%IERR + 1
156  END IF
157 
158  equilibrium_out(1)%global_param%psi_ax = equilibrium_out(1)%profiles_1d%psi(1)
159  equilibrium_out(1)%global_param%psi_bound = equilibrium_out(1)%profiles_1d%psi(neq)
160 
161  noll_cond=1
162  IF (ASSOCIATED (coreprof_in(1)%profiles1d%q%value)) THEN
163  noll_cond=1
164  if (any(coreprof_in(1)%profiles1d%q%value.ne. 0.0)) noll_cond=0
165  end if
166 
167  IF (noll_cond.eq.0) THEN
168  CALL l3interp(coreprof_in(1)%profiles1d%q%value, coreprof_in(1)%rho_tor, nrho, &
169  equilibrium_out(1)%profiles_1d%q, equilibrium_out(1)%profiles_1d%rho_tor, neq)
170  else if (associated(equilibrium_in(1)%profiles_1d%q)) then
171  allocate(q_tmp(nrho))
172  call l3interp(equilibrium_in(1)%profiles_1d%q,equilibrium_in(1)%profiles_1d%rho_tor,neq, &
173  q_tmp,coreprof_in(1)%rho_tor,nrho)
174  call l3interp(q_tmp,coreprof_in(1)%rho_tor,nrho, &
175  equilibrium_out(1)%profiles_1d%q,equilibrium_out(1)%profiles_1d%rho_tor,neq)
176  ELSE
177  equilibrium_out(1)%profiles_1d%q = 0.0_r8
178  diag%ERROR_MESSAGE = trim(adjustl(diag%ERROR_MESSAGE))//"No q supplied; "
179  diag%IERR = diag%IERR + 1
180  END IF
181 
182  noll_cond=1
183  IF (ASSOCIATED (coreprof_in(1)%profiles1d%jtot%value)) THEN
184  noll_cond=1
185  if (any(coreprof_in(1)%profiles1d%jtot%value.ne. 0.0)) noll_cond=0
186  end if
187 
188  IF (noll_cond.eq.0) THEN
189  CALL l3interp(coreprof_in(1)%profiles1d%jtot%value, coreprof_in(1)%rho_tor, nrho, &
190  equilibrium_out(1)%profiles_1d%jparallel,equilibrium_out(1)%profiles_1d%rho_tor, neq)
191  else if (associated(equilibrium_in(1)%profiles_1d%jparallel)) then
192  allocate(jtot_tmp(nrho))
193  call l3interp(equilibrium_in(1)%profiles_1d%jparallel,equilibrium_in(1)%profiles_1d%rho_tor,neq, &
194  jtot_tmp,coreprof_in(1)%rho_tor,nrho)
195  call l3interp(jtot_tmp,coreprof_in(1)%rho_tor,nrho, &
196  equilibrium_out(1)%profiles_1d%jparallel,equilibrium_out(1)%profiles_1d%rho_tor,neq)
197 ! dpc
198 ! IF(ALLOCATED(j_edge_mask)) THEN
199 ! write(*,*) 'Smashing edge jparallel'
200 ! EQUILIBRIUM_OUT(1)%profiles_1d%jparallel = EQUILIBRIUM_OUT(1)%profiles_1d%jparallel * j_edge_mask
201 ! ENDIF
202 ! cpd
203  ELSE
204  equilibrium_out(1)%profiles_1d%jparallel = 0.0_r8
205  diag%ERROR_MESSAGE = trim(adjustl(diag%ERROR_MESSAGE))//"No Jpar supplied; "
206  diag%IERR = diag%IERR + 1
207  END IF
208 
209 
210  noll_cond=1
211  IF (ASSOCIATED (coreprof_in(1)%profiles1d%jphi%value)) THEN
212  noll_cond=1
213  if (any(coreprof_in(1)%profiles1d%jphi%value.ne. 0.0)) noll_cond=0
214  end if
215 
216  IF (noll_cond.eq.0) THEN
217  CALL l3interp(coreprof_in(1)%profiles1d%jphi%value, coreprof_in(1)%rho_tor, nrho, &
218  equilibrium_out(1)%profiles_1d%jphi,equilibrium_out(1)%profiles_1d%rho_tor, neq)
219  else if (associated(equilibrium_in(1)%profiles_1d%jphi)) then
220  allocate(jphi_tmp(nrho))
221  call l3interp(equilibrium_in(1)%profiles_1d%jphi,equilibrium_in(1)%profiles_1d%rho_tor,neq, &
222  jphi_tmp,coreprof_in(1)%rho_tor,nrho)
223  call l3interp(jphi_tmp,coreprof_in(1)%rho_tor,nrho, &
224  equilibrium_out(1)%profiles_1d%jphi,equilibrium_out(1)%profiles_1d%rho_tor,neq)
225 ! dpc
226 ! IF(ALLOCATED(j_edge_mask)) THEN
227 ! write(*,*) 'Smashing edge jphi'
228 ! EQUILIBRIUM_OUT(1)%profiles_1d%jphi = EQUILIBRIUM_OUT(1)%profiles_1d%jphi * j_edge_mask
229 ! ENDIF
230 ! cpd
231  ELSE
232  equilibrium_out(1)%profiles_1d%jphi = 0.0_r8
233  diag%ERROR_MESSAGE = trim(adjustl(diag%ERROR_MESSAGE))//"No jphi supplied; "
234  diag%IERR = diag%IERR + 1
235  END IF
236 
237 
238 
239 
240 
241  IF(ASSOCIATED (coreprof_in(1)%profiles1d%pr_perp%value)) THEN
242  CALL l3interp(coreprof_in(1)%profiles1d%pr_perp%value, coreprof_in(1)%rho_tor, nrho, &
243  equilibrium_out(1)%profiles_1d%pressure, equilibrium_out(1)%profiles_1d%rho_tor, neq)
244  ELSE
245  pressure = coreprof_in(1)%ne%value * coreprof_in(1)%te%value * itm_ev
246 
247  rho_loop1: DO irho=1,nrho
248  ion_loop1: DO iion=1,nion
249  pressure(irho) = pressure(irho) + &
250  coreprof_in(1)%ni%value(irho,iion) * coreprof_in(1)%ti%value(irho,iion) * itm_ev
251  END DO ion_loop1
252  END DO rho_loop1
253 
254  CALL l3interp(pressure, coreprof_in(1)%rho_tor, nrho, &
255  equilibrium_out(1)%profiles_1d%pressure, equilibrium_out(1)%profiles_1d%rho_tor, neq)
256  diag%ERROR_MESSAGE = trim(adjustl(diag%ERROR_MESSAGE))//"Total perpendicular pressure not found in coreprof%profiles1d%pr_perp. Pressure calculated from thermal profiles (fast ion pressure is not included); "
257  diag%IERR = diag%IERR + 1
258  ENDIF
259 !JOFE>
260 
261 
262 
263 
264 ! dpc 2011-08-10 Force psi to be monotonic
265 
266  rho_tor => equilibrium_out(1)%profiles_1d%rho_tor
267  psi => equilibrium_out(1)%profiles_1d%psi
268 
269  IF(psi(1) .LT. psi(neq)) THEN
270  sign = +1.0_r8
271  ELSE IF(psi(1) .GT. psi(neq)) THEN
272  sign = -1.0_r8
273  ELSE
274  diag%ERROR_MESSAGE = trim(adjustl(diag%ERROR_MESSAGE))//"PSI end points equal; "
275  diag%IERR = - 1
276  goto 112
277  ENDIF
278  ieq = 1
279  DO WHILE(ieq .LT. neq)
280  IF(psi(ieq)*sign .LT. psi(ieq+1)*sign) THEN
281  ieq = ieq + 1
282  ELSE ! non-monotonic
283  ieq1 = ieq+1
284  IF(ieq1 .LT. neq) THEN
285  DO WHILE(ieq1 .LT. neq .AND. &
286  psi(ieq)*sign .GE. psi(ieq1)*sign)
287  ieq1 = ieq1 + 1
288  ENDDO
289  ENDIF
290  IF(ieq1 .LT. neq) THEN ! problem lies between IEQ and IEQ1
291  WRITE (char1,'(i3)') ieq
292  WRITE (char2,'(i3)') ieq1
293  diag%ERROR_MESSAGE = trim(adjustl(diag%ERROR_MESSAGE))//"fixed problem in PSI between "//trim(char1)//" and"//trim(char2)//"; "
294  diag%IERR = diag%IERR + 1
295  DO ieq2 = ieq+1, ieq1-1
296  psi(ieq2) = psi(ieq) + &
297  (rho_tor(ieq2) - rho_tor(ieq))/(rho_tor(ieq1) - rho_tor(ieq)) * &
298  (psi(ieq1) - psi(ieq))
299  ENDDO
300  ieq = ieq1
301  ELSE ! problem lies to the left of IEQ
302  ieq1 = ieq
303  DO WHILE(ieq1 .GT. 1 .AND. &
304  psi(ieq1)*sign .GE. psi(neq)*sign)
305  ieq1 = ieq1 - 1
306  ENDDO
307  IF(ieq1 .GT. 1) THEN ! problem lies between IEQ1 and NEQ
308  WRITE (char1,'(i3)') ieq
309  WRITE (char2,'(i3)') neq
310  diag%ERROR_MESSAGE = trim(adjustl(diag%ERROR_MESSAGE))//"fixed problem in PSI between "//trim(char1)//" and"//trim(char2)//"; "
311  diag%IERR = diag%IERR + 1
312  DO ieq2 = ieq1+1, neq-1
313  psi(ieq2) = psi(ieq1) + &
314  (rho_tor(ieq2) - rho_tor(ieq1))/(rho_tor(neq) - rho_tor(ieq1)) * &
315  (psi(neq) - psi(ieq1))
316  ENDDO
317  ieq = neq
318  ELSE
319  diag%ERROR_MESSAGE = trim(adjustl(diag%ERROR_MESSAGE))//"this should not happen! "
320  diag%IERR = - 1
321  goto 112
322  ENDIF
323  ENDIF
324  ENDIF
325  END DO
326 
327 ! Force rho_vol to be monotonic
328 
329  rho_tor => equilibrium_out(1)%profiles_1d%rho_tor
330  rho_vol => equilibrium_out(1)%profiles_1d%rho_vol
331 
332  IF(rho_vol(1) .LT. rho_vol(neq)) THEN
333  sign = +1.0_r8
334  ELSE IF(rho_vol(1) .GT. rho_vol(neq)) THEN
335  sign = -1.0_r8
336  ELSE
337  diag%ERROR_MESSAGE = trim(adjustl(diag%ERROR_MESSAGE))//"RHO_VOL end points equal; "
338  diag%IERR = diag%IERR + 1
339  ENDIF
340  ieq = 1
341  DO WHILE(ieq .LT. neq)
342  IF(rho_vol(ieq)*sign .LT. rho_vol(ieq+1)*sign) THEN
343  ieq = ieq + 1
344  ELSE ! non-monotonic
345  ieq1 = ieq+1
346  IF(ieq1 .LT. neq) THEN
347  DO WHILE(ieq1 .LT. neq .AND. &
348  rho_vol(ieq)*sign .GE. rho_vol(ieq1)*sign)
349  ieq1 = ieq1 + 1
350  ENDDO
351  ENDIF
352  IF(ieq1 .LT. neq) THEN ! problem lies between IEQ and IEQ1
353  WRITE (char1,'(i3)') ieq
354  WRITE (char2,'(i3)') ieq1
355  diag%ERROR_MESSAGE = trim(adjustl(diag%ERROR_MESSAGE))//"fixed problem in RHO_VOL between "//trim(char1)//" and"//trim(char2)//"; "
356  diag%IERR = diag%IERR + 1
357  DO ieq2 = ieq+1, ieq1-1
358  rho_vol(ieq2) = rho_vol(ieq) + &
359  (rho_tor(ieq2) - rho_tor(ieq))/(rho_tor(ieq1) - rho_tor(ieq)) * &
360  (rho_vol(ieq1) - rho_vol(ieq))
361  ENDDO
362  ieq = ieq1
363  ELSE ! problem lies to the left of IEQ
364  ieq1 = ieq
365  DO WHILE(ieq1 .GT. 1 .AND. &
366  rho_vol(ieq1)*sign .GE. rho_vol(neq)*sign)
367  ieq1 = ieq1 - 1
368  ENDDO
369  IF(ieq1 .GT. 1) THEN ! problem lies between IEQ1 and NEQ
370  WRITE (char1,'(i3)') ieq
371  WRITE (char2,'(i3)') neq
372  diag%ERROR_MESSAGE = trim(adjustl(diag%ERROR_MESSAGE))//"fixed problem in RHO_VOL between "//trim(char1)//" and"//trim(char2)//"; "
373  diag%IERR = diag%IERR + 1
374  DO ieq2 = ieq1+1, neq-1
375  rho_vol(ieq2) = rho_vol(ieq1) + &
376  (rho_tor(ieq2) - rho_tor(ieq1))/(rho_tor(neq) - rho_tor(ieq1)) * &
377  (rho_vol(neq) - rho_vol(ieq1))
378  ENDDO
379  ieq = neq
380  ELSE
381  diag%ERROR_MESSAGE = trim(adjustl(diag%ERROR_MESSAGE))//"this should not happen! "
382  diag%IERR = - 1
383  goto 112
384  ENDIF
385  ENDIF
386  ENDIF
387  END DO
388 
389 ! end of dpc addition
390 
391  equilibrium_out(1)%global_param%psi_ax = equilibrium_out(1)%profiles_1d%psi(1)
392  equilibrium_out(1)%global_param%psi_bound = equilibrium_out(1)%profiles_1d%psi(neq)
393  equilibrium_out(1)%global_param%volume = equilibrium_out(1)%profiles_1d%volume(neq)
394 
395  DEALLOCATE (pressure)
396 ! dpc
397  IF(ALLOCATED(j_edge_mask)) DEALLOCATE(j_edge_mask)
398 ! cpd
399 ! +++ CODEPARAM:
400 112 output_diag = "EQUILIBRIUM_INPUT: "//trim(adjustl(diag%ERROR_MESSAGE))
401  nline = floor(len_trim(adjustl(output_diag))/132.001)+1
402  IF(ASSOCIATED(equilibrium_out(1)%codeparam%codename)) &
403  DEALLOCATE(equilibrium_out(1)%codeparam%codename)
404  IF(ASSOCIATED(equilibrium_out(1)%codeparam%codeversion)) &
405  DEALLOCATE(equilibrium_out(1)%codeparam%codeversion)
406  IF(ASSOCIATED(equilibrium_out(1)%codeparam%output_diag)) &
407  DEALLOCATE(equilibrium_out(1)%codeparam%output_diag)
408  ALLOCATE (equilibrium_out(1)%codeparam%codename(1))
409  ALLOCATE (equilibrium_out(1)%codeparam%codeversion(1))
410  ALLOCATE (equilibrium_out(1)%codeparam%output_diag(nline))
411 
412  equilibrium_out(1)%codeparam%codename = 'EQILIBRIUM_INPUT'
413  equilibrium_out(1)%codeparam%codeversion = 'EQILIBRIUM_INPUT_4.10b.10'
414  equilibrium_out(1)%codeparam%output_flag = diag%IERR
415  DO i = 1,nline
416  equilibrium_out(1)%codeparam%output_diag(i) = output_diag(((i-1)*132+1):(i*132))
417  END DO
418 
419 
420  if(allocated(psi_tmp)) deallocate(psi_tmp)
421  if(allocated(q_tmp)) deallocate(q_tmp)
422  if(allocated(jtot_tmp)) deallocate(jtot_tmp)
423  if(allocated(jphi_tmp)) deallocate(jphi_tmp)
424 
425  RETURN
426 
427  END SUBROUTINE equil_input
428 
429 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
430 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
431 
432 
433 
434 
435 
436 
437 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
438 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
439 
440  SUBROUTINE ets_input &
441  (equilibrium_in, coreprof_out)
442 
443 !-------------------------------------------------------!
444 ! This routine updates the q-profile in the COREPROF! !
445 !-------------------------------------------------------!
446 ! Source: --- !
447 ! Developers: D.Kalupin !
448 ! Kontacts: D.Kalupin@fz-juelich.de !
449 ! !
450 ! Comments: created for V&V between ETS and !
451 ! ASTRA !
452 ! !
453 !-------------------------------------------------------!
454 
455 
456  USE euitm_schemas
457  USE itm_constants
458 
459  IMPLICIT NONE
460 
461 ! +++ CPO derived types:
462  TYPE (type_equilibrium), POINTER :: equilibrium_in(:) !input CPO with geometry quantities from previous time
463  TYPE (type_coreprof), POINTER :: coreprof_out(:) !output CPO with plasma profiles
464 
465 
466 ! +++ Dimensions:
467  INTEGER :: nrho !number of coreprof radial points
468  INTEGER :: neq !number and equilibrium radial points
469 
470 ! +++ Internal parameters(0-D parameters):
471 
472 ! Input
473 !! REAL (R8), ALLOCATABLE :: QSF(:)
474 
475 ! +++ Set dimensions and allocate parameters:
476  neq = SIZE (equilibrium_in(1)%profiles_1d%rho_tor, dim=1)
477  nrho = SIZE (coreprof_out(1)%rho_tor, dim=1)
478 
479 !! ALLOCATE (QSF(NRHO))
480 
481 !! QSF = EQUILIBRIUM_IN(1)%profiles_1d%q
482 
483 ! +++ Fill the output COREPROF_OUT CPO:
484 !! COREPROF_OUT(1)%profiles1d%q%value = QSF
485 
486  CALL l3interp(equilibrium_in(1)%profiles_1d%q, &
487  equilibrium_in(1)%profiles_1d%rho_tor / equilibrium_in(1)%profiles_1d%rho_tor(neq), &
488  neq, &
489  coreprof_out(1)%profiles1d%q%value, &
490  coreprof_out(1)%rho_tor / coreprof_out(1)%rho_tor(nrho), &
491  nrho)
492 
493 
494 !! DEALLOCATE (QSF)
495 
496 
497  RETURN
498 
499  END SUBROUTINE ets_input
500 
501 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
502 ! + + + + + + + + + + + + + + + + + + + + + + + + + + + +
503 
504 END MODULE equilibrium_input
subroutine equil_input(COREPROF_IN, TOROIDFIELD_IN, EQUILIBRIUM_IN, EQUILIBRIUM_OUT)
subroutine ets_input(EQUILIBRIUM_IN, COREPROF_OUT)
EQUILIBRIUM_INPUT.
subroutine l3interp(y_in, x_in, nr_in, y_out, x_out, nr_out)
Definition: l3interp.f90:1
The module declares types of variables used in ETS (transport code)
Definition: ets_plasma.f90:8
real(r8) function pressure(flux)