ETS  \$Id: Doxyfile 2162 2020-02-26 14:16:09Z g2dpc $
 All Classes Files Functions Variables Pages
neo_test.F90
Go to the documentation of this file.
1 program neo_test
2 
3  use itm_constants
4  USE euitm_schemas
5  use euitm_routines
6  use read_structures
7  use write_structures
8  use deallocate_structures
9 
10  implicit none
11 
12  type(type_coreprof), pointer :: coreprof(:)
13  type(type_coretransp), pointer :: coretransp(:)
14  type(type_coresource), pointer :: coresource(:)
15  type(type_coreimpur), pointer :: coreimpur(:)
16  type(type_equilibrium), pointer :: equilibrium(:)
17  type(type_toroidfield), pointer :: toroidfield(:)
18  type(type_neoclassic), pointer :: neoclassic_neo(:), neoclassic_neowes(:), neoclassic_itmneoart(:), neoclassic_neos(:)
19 
20  integer :: idx, shot, run, interpol, nrho, irho, nion, iion
21  REAL (R8) :: time
22  character(len=17) :: filename
23  logical :: exist
24 
25  INTERFACE
26  SUBROUTINE neo(EQUILIBRIUM,COREPROF,NEOCLASSIC)
27  USE euitm_schemas
28  USE itm_types
29  TYPE (type_equilibrium), POINTER :: equilibrium(:)
30  TYPE (type_coreprof), POINTER :: coreprof(:)
31  TYPE (type_neoclassic), POINTER :: neoclassic(:)
32  END SUBROUTINE neo
33  END INTERFACE
34 
35  INTERFACE
36  SUBROUTINE neowes(EQUILIBRIUM,COREPROF,NEOCLASSIC)
37  USE euitm_schemas
38  USE itm_types
39  TYPE (type_equilibrium), POINTER :: equilibrium(:)
40  TYPE (type_coreprof), POINTER :: coreprof(:)
41  TYPE (type_neoclassic), POINTER :: neoclassic(:)
42  END SUBROUTINE neowes
43  END INTERFACE
44 
45 #ifdef GOT_NEOART
46  INTERFACE
47  SUBROUTINE itmneoart(EQUILIBRIUM,COREPROF,NEOCLASSIC)
48  USE euitm_schemas
49  USE itm_types
50  TYPE (type_equilibrium), POINTER :: equilibrium(:)
51  TYPE (type_coreprof), POINTER :: coreprof(:)
52  TYPE (type_neoclassic), POINTER :: neoclassic(:)
53  END SUBROUTINE itmneoart
54  END INTERFACE
55 #endif
56 
57 #ifdef GOT_NEOS
58  INTERFACE
59  SUBROUTINE signeojbs(EQUILIBRIUM,COREPROF,NEOCLASSIC)
60  USE euitm_schemas
61  USE itm_types
62  TYPE (type_equilibrium), POINTER :: equilibrium(:)
63  TYPE (type_coreprof), POINTER :: coreprof(:)
64  TYPE (type_neoclassic), POINTER :: neoclassic(:)
65  END SUBROUTINE signeojbs
66  END INTERFACE
67 #endif
68 
69 ! shot = 5 ; run = 45 ! BDSEQ, nion=1, nimp=0
70 ! shot = 5 ; run = 47 ! HELENA, nion=1, nimp=0
71 ! shot = 5 ; run = 63 ! BDSEQ, nion=4, nimp=0
72 ! shot = 5 ; run = 64 ! HELENA, nion=4, nimp=0
73 ! shot = 5 ; run = 123 ! EMEQ, nion=1, nimp=1, W
74 ! shot = 5 ; run = 142 ! cyl, nion=1, nimp=0 (H)
75  shot = 5 ; run = 143 ! bdseq, nion=1, nimp=0 (H)
76 ! shot = 5 ; run = 144 ! emeq, nion=1, nimp=0 (H)
77 ! shot = 5 ; run = 145 ! helena, nion=1, nimp=0 (H)
78 
79  time = 10.0_r8
80 
81  interpol = 1
82 
83  allocate(coreprof(1),coretransp(1),coresource(1), &
84  coreimpur(1),equilibrium(1),toroidfield(1))
85 #ifdef UAL
86  CALL euitm_open('euitm', shot, run, idx)
87  CALL euitm_get_slice(idx, 'coreprof', coreprof(1), time, interpol)
88  CALL euitm_get_slice(idx, 'equilibrium', equilibrium(1), time, interpol)
89  CALL euitm_get_slice(idx, 'toroidfield', toroidfield(1), time, interpol)
90 ! the following are not needed for calling the neoclassical routines
91 ! but are needed if we want to write out the ASCII version
92  CALL euitm_get_slice(idx, 'coretransp', coretransp(1), time, interpol)
93  CALL euitm_get_slice(idx, 'coresource', coresource(1), time, interpol)
94  CALL euitm_get_slice(idx, 'coreimpur', coreimpur(1), time, interpol)
95 ! write out an ASCII version if one doesn't already exist
96  write(filename,'(''CPO_'',I6.6,''_'',I6.6)') shot, run
97  inquire(file=filename, exist=exist)
98  if(exist) then
99  write(*,*) filename,' already exists'
100  else
101  call open_write_file(1, filename)
102  call write_cpo(coreprof(1), 'coreprof')
103  call write_cpo(coretransp(1), 'coretransp')
104  call write_cpo(coresource(1), 'coresource')
105  call write_cpo(coreimpur(1), 'coreimpur')
106  call write_cpo(equilibrium(1), 'equilibrium')
107  call write_cpo(toroidfield(1), 'toroidfield')
108  call close_write_file
109  endif
110 #else
111  write(filename,'(''CPO_'',I6.6,''_'',I6.6)') shot, run
112  call open_read_file(1, filename)
113  call read_cpo(coreprof(1), 'coreprof')
114  call read_cpo(coretransp(1), 'coretransp')
115  call read_cpo(coresource(1), 'coresource')
116  call read_cpo(coreimpur(1), 'coreimpur')
117  call read_cpo(equilibrium(1), 'equilibrium')
118  call read_cpo(toroidfield(1), 'toroidfield')
119  call close_read_file
120 #endif
121 
122  nrho=size(coreprof(1)%rho_tor)
123  nion=size(coreprof(1)%ni%value,dim=2)
124  allocate(coreprof(1)%rho_tor_norm(nrho))
125  coreprof(1)%rho_tor_norm = coreprof(1)%rho_tor / coreprof(1)%rho_tor(nrho)
126  equilibrium(1)%eqgeometry%elongation=1.0_r8
127  if(.not.associated(coreprof(1)%profiles1d%pe%value)) then
128  write(*,*) 'Calculating coreprof(1)%profiles1d%pe%value'
129  allocate(coreprof(1)%profiles1d%pe%value(nrho))
130  coreprof(1)%profiles1d%pe%value = coreprof(1)%te%value * coreprof(1)%ne%value * itm_ev
131  endif
132  if(.not.associated(coreprof(1)%profiles1d%pi%value)) then
133  write(*,*) 'Calculating coreprof(1)%profiles1d%pi%value'
134  allocate(coreprof(1)%profiles1d%pi%value(nrho,nion))
135  coreprof(1)%profiles1d%pi%value = coreprof(1)%ti%value * coreprof(1)%ni%value * itm_ev
136  endif
137  if(.not.associated(coreprof(1)%profiles1d%pr_th%value)) then
138  write(*,*) 'Calculating coreprof(1)%profiles1d%pr_th%value'
139  allocate(coreprof(1)%profiles1d%pr_th%value(nrho))
140  coreprof(1)%profiles1d%pr_th%value = coreprof(1)%profiles1d%pe%value + sum(coreprof(1)%profiles1d%pi%value,dim=2)
141  endif
142 
143  write(*,*) (maxval(equilibrium(1)%coord_sys%position%z)-minval(equilibrium(1)%coord_sys%position%z))/ &
144  (maxval(equilibrium(1)%coord_sys%position%r)-minval(equilibrium(1)%coord_sys%position%r))
145 
146  call neo(equilibrium,coreprof,neoclassic_neo)
147  call neowes(equilibrium,coreprof,neoclassic_neowes)
148 #ifdef GOT_NEOS
149  call signeojbs(equilibrium,coreprof,neoclassic_neos)
150 #endif
151  if(coreprof(1)%rho_tor(1).eq.0.0_r8) then
152  coreprof(1)%rho_tor(1)=coreprof(1)%rho_tor(2)/1e10_r8
153  endif
154 #ifdef GOT_NEOART
155  call itmneoart(equilibrium,coreprof,neoclassic_itmneoart)
156 #endif
157 
158  open(10,file='neo.dat')
159  write(10,'(a,4(i4))') '# ',nrho, nion, shot, run
160  do irho=1,size(neoclassic_neo(1)%sigma)
161  write(10,1000) neoclassic_neo(1)%rho_tor(irho), &
162  neoclassic_neo(1)%sigma(irho), &
163  neoclassic_neo(1)%jboot(irho), &
164  neoclassic_neo(1)%ne_neo%flux(irho), &
165  neoclassic_neo(1)%ne_neo%diff_eff(irho), &
166  neoclassic_neo(1)%ne_neo%vconv_eff(irho), &
167  neoclassic_neo(1)%te_neo%flux(irho), &
168  neoclassic_neo(1)%te_neo%diff_eff(irho), &
169  neoclassic_neo(1)%te_neo%vconv_eff(irho), &
170  (neoclassic_neo(1)%ni_neo%flux(irho,iion),iion=1,nion), &
171  (neoclassic_neo(1)%ni_neo%diff_eff(irho,iion),iion=1,nion), &
172  (neoclassic_neo(1)%ni_neo%vconv_eff(irho,iion),iion=1,nion), &
173  (neoclassic_neo(1)%ti_neo%flux(irho,iion),iion=1,nion), &
174  (neoclassic_neo(1)%ti_neo%diff_eff(irho,iion),iion=1,nion), &
175  (neoclassic_neo(1)%ti_neo%vconv_eff(irho,iion),iion=1,nion)
176  enddo
177  close(10)
178 
179  open(10,file='neowes.dat')
180  write(10,'(a,4(i4))') '# ',nrho, nion, shot, run
181  do irho=1,size(neoclassic_neowes(1)%sigma)
182  write(10,1000) neoclassic_neowes(1)%rho_tor(irho), &
183  neoclassic_neowes(1)%sigma(irho), &
184  neoclassic_neowes(1)%jboot(irho), &
185  neoclassic_neowes(1)%ne_neo%flux(irho), &
186  neoclassic_neowes(1)%ne_neo%diff_eff(irho), &
187  neoclassic_neowes(1)%ne_neo%vconv_eff(irho), &
188  neoclassic_neowes(1)%te_neo%flux(irho), &
189  neoclassic_neowes(1)%te_neo%diff_eff(irho), &
190  neoclassic_neowes(1)%te_neo%vconv_eff(irho), &
191  (neoclassic_neowes(1)%ni_neo%flux(irho,iion),iion=1,nion), &
192  (neoclassic_neowes(1)%ni_neo%diff_eff(irho,iion),iion=1,nion), &
193  (neoclassic_neowes(1)%ni_neo%vconv_eff(irho,iion),iion=1,nion), &
194  (neoclassic_neowes(1)%ti_neo%flux(irho,iion),iion=1,nion), &
195  (neoclassic_neowes(1)%ti_neo%diff_eff(irho,iion),iion=1,nion), &
196  (neoclassic_neowes(1)%ti_neo%vconv_eff(irho,iion),iion=1,nion)
197  enddo
198  close (10)
199 
200 #ifdef GOT_NEOART
201  open(10,file='itmneoart.dat')
202  write(10,'(a,4(i4))') '# ',nrho, nion, shot, run
203  do irho=1,size(neoclassic_itmneoart(1)%sigma)
204  write(10,1000) neoclassic_itmneoart(1)%rho_tor(irho), &
205  neoclassic_itmneoart(1)%sigma(irho), &
206  neoclassic_itmneoart(1)%jboot(irho), &
207  neoclassic_itmneoart(1)%ne_neo%flux(irho), &
208  neoclassic_itmneoart(1)%ne_neo%diff_eff(irho), &
209  neoclassic_itmneoart(1)%ne_neo%vconv_eff(irho), &
210  neoclassic_itmneoart(1)%te_neo%flux(irho), &
211  neoclassic_itmneoart(1)%te_neo%diff_eff(irho), &
212  neoclassic_itmneoart(1)%te_neo%vconv_eff(irho), &
213  (neoclassic_itmneoart(1)%ni_neo%flux(irho,iion),iion=1,nion), &
214  (neoclassic_itmneoart(1)%ni_neo%diff_eff(irho,iion),iion=1,nion), &
215  (neoclassic_itmneoart(1)%ni_neo%vconv_eff(irho,iion),iion=1,nion), &
216  (neoclassic_itmneoart(1)%ti_neo%flux(irho,iion),iion=1,nion), &
217  (neoclassic_itmneoart(1)%ti_neo%diff_eff(irho,iion),iion=1,nion), &
218  (neoclassic_itmneoart(1)%ti_neo%vconv_eff(irho,iion),iion=1,nion)
219  enddo
220  close (10)
221 #ifdef GOT_NEOART
222 
223 #ifdef GOT_NEOS
224  open(10,file='neos.dat')
225  write(10,'(a,4(i4))') '# ',nrho, nion, shot, run
226  do irho=1,size(neoclassic_neos(1)%sigma)
227  write(10,1000) neoclassic_neos(1)%rho_tor(irho), &
228  neoclassic_neos(1)%sigma(irho), &
229  neoclassic_neos(1)%jboot(irho)
230 !, &
231 ! NEOCLASSIC_NEOS(1)%ne_neo%flux(irho), &
232 ! NEOCLASSIC_NEOS(1)%ne_neo%diff_eff(irho), &
233 ! NEOCLASSIC_NEOS(1)%ne_neo%vconv_eff(irho), &
234 ! NEOCLASSIC_NEOS(1)%te_neo%flux(irho), &
235 ! NEOCLASSIC_NEOS(1)%te_neo%diff_eff(irho), &
236 ! NEOCLASSIC_NEOS(1)%te_neo%vconv_eff(irho), &
237 ! (NEOCLASSIC_NEOS(1)%ni_neo%flux(irho,iion),iion=1,NION), &
238 ! (NEOCLASSIC_NEOS(1)%ni_neo%diff_eff(irho,iion),iion=1,NION), &
239 ! (NEOCLASSIC_NEOS(1)%ni_neo%vconv_eff(irho,iion),iion=1,NION), &
240 ! (NEOCLASSIC_NEOS(1)%ti_neo%flux(irho,iion),iion=1,NION), &
241 ! (NEOCLASSIC_NEOS(1)%ti_neo%diff_eff(irho,iion),iion=1,NION), &
242 ! (NEOCLASSIC_NEOS(1)%ti_neo%vconv_eff(irho,iion),iion=1,NION)
243  enddo
244  close(10)
245 #endif
246 
247  if(associated(coreprof(1)%psi%sigma_par%value)) then
248  open(10,file='sigma_par.dat')
249  do irho=1,size(coreprof(1)%psi%sigma_par%value)
250  write(10,1000) coreprof(1)%rho_tor(irho), coreprof(1)%psi%sigma_par%value(irho)
251  enddo
252  close (10)
253  endif
254 
255  open(10,file='rho_tor_te_ne.dat')
256  do irho=1,size(coreprof(1)%psi%sigma_par%value)
257  write(10,1000) coreprof(1)%rho_tor(irho), coreprof(1)%te%value(irho), coreprof(1)%ne%value(irho)
258  enddo
259  close (10)
260 
261  call deallocate_cpo(coreprof)
262  call deallocate_cpo(coretransp)
263  call deallocate_cpo(coresource)
264  call deallocate_cpo(coreimpur)
265  call deallocate_cpo(equilibrium)
266  call deallocate_cpo(toroidfield)
267  call deallocate_cpo(neoclassic_neo)
268 #ifdef GOT_NEOS
269  call deallocate_cpo(neoclassic_neos)
270 #endif
271  call deallocate_cpo(neoclassic_neowes)
272 #ifdef GOT_NEOART
273  call deallocate_cpo(neoclassic_itmneoart)
274 #endif
275 
276 1000 format(1p,100(100g15.6))
277 end program neo_test
subroutine neowes(eq, coreprof, neoclassic)
Definition: neowes.F90:4
real(r8) function p(a, x, xr, xs, yr, ys, psi, psir, F_dia)
subroutine neo(equilibrium, coreprof, neoclassic)
Definition: neo.f90:2817