This repository has been archived by the owner on Jul 29, 2022. It is now read-only.
-
Notifications
You must be signed in to change notification settings - Fork 0
/
stella_diagnostics.f90
1439 lines (1169 loc) · 52.6 KB
/
stella_diagnostics.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
module stella_diagnostics
implicit none
public :: init_stella_diagnostics, finish_stella_diagnostics
public :: diagnose_stella
public :: nsave
public :: omega_vs_time
public :: omega_vs_time_short
public :: navg
private
integer :: ntg_out
integer :: nwrite, nsave, navg
real :: halfnavg
integer :: stdout_unit, fluxes_unit, omega_unit
logical :: save_for_restart
logical :: write_omega
logical :: write_moments
logical :: write_phi_vs_time
logical :: write_gvmus
logical :: write_gzvs
logical :: write_kspectra
logical :: write_radial_fluxes
logical :: write_fluxes_kxkyz
logical :: flux_norm
! Arrays needed for averaging in x,y,z
real, dimension (:), allocatable :: pflux_avg, vflux_avg, qflux_avg, heat_avg
real, dimension (:,:,:), allocatable :: pflux, vflux, qflux, exchange
! Needed for calculating growth rates and frequencies
complex, dimension (:,:,:), allocatable :: omega_vs_time, omega_vs_time_short
! Initialized variables
integer :: nout = 1
logical :: diagnostics_initialized = .false.
! Debugging
logical :: debug = .false.
contains
!==============================================
!====== INITIATE STELLA DIAGNOSTICS ===========
!==============================================
! Broadcast the parameters from the namelist "stella_diagnostics_knobs"
! and open/append the netcdf file and the ascii files.
subroutine init_stella_diagnostics (restart, tstart)
use zgrid, only: init_zgrid
use kt_grids, only: init_kt_grids
use physics_parameters, only: init_physics_parameters
use run_parameters, only: init_run_parameters
use species, only: init_species
use dist_fn, only: init_dist_fn
use init_g, only: init_init_g
use stella_io, only: init_stella_io, get_nout
use mp, only: broadcast, proc0
implicit none
logical, intent (in) :: restart
real, intent (in) :: tstart
! Only initialize the diagnostics once
if (diagnostics_initialized) return
diagnostics_initialized = .true.
! Only debug on the first processor
debug = debug .and. proc0
! Make sure the other routines are intialized
call init_zgrid
call init_physics_parameters
call init_kt_grids
call init_run_parameters
call init_species
call init_init_g
call init_dist_fn
! Read the namelist "stella_diagnostics_knobs" in the input file
call read_parameters
call allocate_arrays
! Broadcast the variables to all processors
call broadcast (nwrite)
call broadcast (navg)
call broadcast (nsave)
call broadcast (save_for_restart)
call broadcast (write_omega)
call broadcast (write_kspectra)
call broadcast (write_moments)
call broadcast (write_phi_vs_time)
call broadcast (write_gvmus)
call broadcast (write_gzvs)
call broadcast (write_radial_fluxes)
call broadcast (write_fluxes_kxkyz)
call broadcast (flux_norm)
! Initiate the netcdf file with extension '.out.nc'
call init_stella_io (restart, write_phi_vs_time, write_kspectra, &
write_gvmus, write_gzvs, write_moments, write_radial_fluxes, &
write_fluxes_kxkyz)
! Open the '.out', '.fluxes' and '.omega' file
if (proc0) call open_loop_ascii_files (restart)
! Get the final position [[nout]] of the time axis in the netcdf file
if (proc0) call get_nout (tstart,nout)
call broadcast (nout)
end subroutine init_stella_diagnostics
!==============================================
!============== READ PARAMETERS ===============
!==============================================
subroutine read_parameters
use mp, only: proc0
use file_utils, only: input_unit_exist
use zgrid, only: nperiod, nzed
use physics_flags, only: radial_variation
implicit none
logical :: exist
integer :: in_file
namelist /stella_diagnostics_knobs/ nwrite, navg, nsave, &
save_for_restart, write_phi_vs_time, write_gvmus, write_gzvs, &
write_omega, write_kspectra, write_moments, write_radial_fluxes, &
write_fluxes_kxkyz, flux_norm
if (proc0) then
nwrite = 50
navg = 50
halfnavg = navg/2
nsave = -1
save_for_restart = .false.
write_omega = .false.
write_phi_vs_time = .false.
write_gvmus = .false.
write_gzvs = .false.
write_kspectra = .false.
write_moments = .false.
write_radial_fluxes = radial_variation
write_fluxes_kxkyz = .false.
flux_norm = .true.
in_file = input_unit_exist ("stella_diagnostics_knobs", exist)
if (exist) read (unit=in_file, nml=stella_diagnostics_knobs)
if (.not. save_for_restart) nsave = -1
end if
ntg_out = nzed/2 + (nperiod-1)*nzed
end subroutine read_parameters
!==============================================
!============== ALLOCATE ARRAYS ===============
!==============================================
subroutine allocate_arrays
use species, only: nspec
use kt_grids, only: nakx, naky
implicit none
if (.not.allocated(pflux)) allocate(pflux (nakx,naky,nspec)) ; pflux = 0.
if (.not.allocated(qflux)) allocate(qflux (nakx,naky,nspec)) ; qflux = 0.
if (.not.allocated(vflux)) allocate(vflux (nakx,naky,nspec)) ; vflux = 0.
if (.not.allocated(exchange)) allocate(exchange (nakx,naky,nspec)) ; exchange = 0.
if (.not.allocated(pflux_avg)) allocate(pflux_avg(nspec)) ; pflux_avg = 0.
if (.not.allocated(qflux_avg)) allocate(qflux_avg(nspec)) ; qflux_avg = 0.
if (.not.allocated(vflux_avg)) allocate(vflux_avg(nspec)) ; vflux_avg = 0.
if (.not.allocated(heat_avg)) allocate(heat_avg(nspec)) ; heat_avg = 0.
if (.not.allocated(omega_vs_time)) then
if (write_omega) then
allocate (omega_vs_time(navg,naky,nakx))
omega_vs_time = 0.
else
allocate (omega_vs_time(1,1,1))
navg=1
end if
end if
if(.not.allocated(omega_vs_time_short)) then
if (write_omega) then
allocate(omega_vs_time_short(nint(halfnavg),naky,nakx))
omega_vs_time_short = 0.
else
allocate (omega_vs_time_short(1,1,1))
end if
end if
end subroutine allocate_arrays
!==============================================
!============= OPEN ASCII FILES ===============
!==============================================
! Open the '.out' and the '.fluxes' file.
! When running a new simulation, create a new file or replace an old file.
! When restarting a simulation, append the old files.
subroutine open_loop_ascii_files(restart)
use file_utils, only: open_output_file
use species, only: nspec
implicit none
logical, intent (in) :: restart
character (3) :: nspec_str
character (100) :: str
logical :: overwrite
! Do not overwrite, but append files, when we restart the simulation.
overwrite = .not.restart
! Open the '.out' and the '.fluxes' files.
call open_output_file (stdout_unit,'.out',overwrite)
call open_output_file (fluxes_unit,'.fluxes',overwrite)
! Create the header for the .fluxes file.
! Every column is made up of 12 spaces, so make sure the headers
! are placed correctly since we have the following columns for nspec=3:
! #time pflx1 pflx2 pflx3 vflx1 vflx2 vflx32 qflx1 qflx2 qflx3
if (.not.restart) then
write (nspec_str,'(i3)') nspec*12
str = trim('(2a12,2a'//trim(nspec_str)//')')
write (fluxes_unit,str) '#time', 'pflx', 'vflx', 'qflx'
end if
! Open the '.omega' file and create its header.
if (write_omega) then
call open_output_file (omega_unit,'.omega',overwrite)
if (.not.restart) then
write (omega_unit,'(7a12)') '#time', 'ky', 'kx', &
'Re[om]', 'Im[om]', 'Re[omavg]', 'Im[omavg]'
end if
end if
end subroutine open_loop_ascii_files
!==============================================
!============= CLOSE ASCII FILES ==============
!==============================================
subroutine close_loop_ascii_files
use file_utils, only: close_output_file
implicit none
call close_output_file (stdout_unit)
call close_output_file (fluxes_unit)
if (write_omega) call close_output_file (omega_unit)
end subroutine close_loop_ascii_files
!==============================================
!============== DIAGNOSE STELLA ===============
!==============================================
subroutine diagnose_stella (istep)
use mp, only: proc0
use constants, only: zi
use redistribute, only: scatter
use fields_arrays, only: phi, apar
use fields_arrays, only: phi_old, phi_corr_QN
use fields, only: fields_updated, advance_fields
use dist_fn_arrays, only: gvmu, gnew
use g_tofrom_h, only: g_to_h
use stella_io, only: write_time_nc
use stella_io, only: write_phi2_nc
use stella_io, only: write_phi_nc
use stella_io, only: write_gvmus_nc
use stella_io, only: write_gzvs_nc
use stella_io, only: write_kspectra_nc
use stella_io, only: write_moments_nc
use stella_io, only: write_radial_fluxes_nc
use stella_io, only: write_fluxes_kxkyz_nc
use stella_io, only: sync_nc
! use stella_io, only: write_symmetry_nc
use stella_time, only: code_time, code_dt
use zgrid, only: nztot, nzgrid, ntubes
use vpamu_grids, only: nmu, nvpa
use species, only: nspec
use kt_grids, only: naky, nakx
use dist_redistribute, only: kxkyz2vmu
use physics_flags, only: radial_variation
use volume_averages, only: volume_average, fieldline_average
use fields, only: advance_fields, fields_updated
use dist_fn_arrays, only: gnew
implicit none
integer, intent (in) :: istep
real :: phi2, apar2
real :: zero
real, dimension (:,:,:), allocatable :: gvmus
real, dimension (:,:,:,:), allocatable :: gzvs
! real, dimension (:,:,:), allocatable :: pflx_zvpa, vflx_zvpa, qflx_zvpa
real, dimension (:), allocatable :: part_flux, mom_flux, heat_flux
real, dimension (:,:), allocatable :: part_flux_x, mom_flux_x, heat_flux_x
real, dimension (:,:), allocatable :: phi2_vs_kxky
real, dimension (:,:,:,:,:), allocatable :: pflx_kxkyz, vflx_kxkyz, qflx_kxkyz
complex, dimension (:,:,:,:,:), allocatable :: density, upar, temperature
complex, dimension (:,:), allocatable :: omega_avg
complex, dimension (:,:), allocatable :: phiavg, phioldavg
complex, dimension (:,:,:,:), allocatable :: phi_out
! calculation of omega requires computation of omega more
! frequently than every nwrite time steps
if (write_omega .and. proc0) then
zero = 100.*epsilon(0.)
if (istep > 0) then
allocate (phiavg(naky,nakx))
allocate (phioldavg(naky,nakx))
call fieldline_average (phi, phiavg)
call fieldline_average (phi_old, phioldavg)
where (abs(phiavg) < zero .or. abs(phioldavg) < zero)
omega_vs_time(mod(istep,navg)+1,:,:) = 0.0
omega_vs_time_short(mod(istep,nint(halfnavg))+1,:,:) = 0.0
elsewhere
omega_vs_time(mod(istep,navg)+1,:,:) = log(phiavg/phioldavg)*zi/code_dt
omega_vs_time_short(mod(istep,nint(halfnavg))+1,:,:) = log(phiavg/phioldavg)*zi/code_dt
end where
deallocate (phiavg, phioldavg)
write(*,*) 'omega vs time', omega_vs_time(mod(istep,navg)+1,1,1)
end if
end if
! only write data to file every nwrite time steps
if (mod(istep,nwrite) /= 0) return
if(radial_variation) fields_updated = .false.
call advance_fields(gnew, phi, apar, dist='gbar')
allocate(phi_out(naky,nakx,-nzgrid:nzgrid,ntubes))
phi_out = phi
if(radial_variation) then
phi_out = phi_out + phi_corr_QN
endif
allocate (part_flux(nspec))
allocate (mom_flux(nspec))
allocate (heat_flux(nspec))
allocate (pflx_kxkyz(naky,nakx,nztot,ntubes,nspec))
allocate (vflx_kxkyz(naky,nakx,nztot,ntubes,nspec))
allocate (qflx_kxkyz(naky,nakx,nztot,ntubes,nspec))
if(write_radial_fluxes) then
allocate (part_flux_x(nakx,nspec))
allocate (mom_flux_x(nakx,nspec))
allocate (heat_flux_x(nakx,nspec))
endif
! obtain turbulent fluxes
if(radial_variation.or.write_radial_fluxes) then
! handle g_to_h in get_fluxes_vmulo to eliminate x^2 terms
! call g_to_h (gnew, phi, fphi, phi_corr_QN)
call get_fluxes_vmulo (gnew, phi_out, part_flux, mom_flux, heat_flux, &
part_flux_x,mom_flux_x, heat_flux_x)
! call g_to_h (gnew, phi, -fphi, phi_corr_QN)
else
call scatter (kxkyz2vmu, gnew, gvmu)
! call g_to_h (gvmu, phi, fphi)
call get_fluxes (gvmu, part_flux, mom_flux, heat_flux, &
pflx_kxkyz, vflx_kxkyz, qflx_kxkyz)
! call g_to_h (gvmu, phi, -fphi)
endif
if (proc0) then
if(write_omega) then
allocate (omega_avg(naky,nakx))
omega_avg = sum(omega_vs_time,dim=1)/real(navg)
else
allocate (omega_avg(1,1))
endif
call volume_average (phi_out, phi2)
call volume_average (apar, apar2)
! Print information to stella.out, the header is printed in stella.f90
write (*,'(A2,I7,A2,ES12.4,A2,ES12.4,A2,ES12.4)') " ",istep," ",code_time," ",code_dt," ", phi2
call write_loop_ascii_files (istep, phi2, apar2, part_flux, mom_flux, heat_flux, &
omega_vs_time(mod(istep,navg)+1,:,:), omega_avg)
! do not need omega_avg again this time step
deallocate (omega_avg)
end if
if (proc0) then
if (debug) write (*,*) 'stella_diagnostics::write_time_nc'
call write_time_nc (nout, code_time)
call write_phi2_nc (nout, phi2)
if (write_phi_vs_time) then
if (debug) write (*,*) 'stella_diagnostics::diagnose_stella::write_phi_nc'
call write_phi_nc (nout, phi_out)
end if
if (write_kspectra) then
if (debug) write (*,*) 'stella_diagnostics::diagnose_stella::write_kspectra'
allocate (phi2_vs_kxky(naky,nakx))
call fieldline_average (real(phi_out*conjg(phi_out)),phi2_vs_kxky)
call write_kspectra_nc (nout, phi2_vs_kxky)
deallocate (phi2_vs_kxky)
end if
if (write_radial_fluxes) then
call write_radial_fluxes_nc(nout,part_flux_x,mom_flux_x,heat_flux_x)
endif
end if
if (write_moments) then
if (debug) write (*,*) 'stella_diagnostics::diagnose_stella::write_moments'
allocate (density(naky,nakx,nztot,ntubes,nspec))
allocate (upar(naky,nakx,nztot,ntubes,nspec))
allocate (temperature(naky,nakx,nztot,ntubes,nspec))
call get_moments (gnew, density, upar, temperature)
if (proc0) call write_moments_nc (nout, density, upar, temperature)
deallocate (density, upar, temperature)
end if
IF (write_fluxes_kxkyz) then
IF (debug) write (*,*) 'stella_diagnostics::diagnose_stella::write_fluxes_kxkyz'
IF (proc0) call write_fluxes_kxkyz_nc (nout, pflx_kxkyz, vflx_kxkyz, qflx_kxkyz)
ENDIF
if (write_gvmus) then
allocate (gvmus(nvpa,nmu,nspec))
if (debug) write (*,*) 'stella_diagnostics::diagnose_stella::get_gvmus'
! note that gvmus is h at this point
call get_gvmus (gvmu, gvmus)
if (debug) write (*,*) 'stella_diagnostics::diagnose_stella::write_gvmus_nc'
if (proc0) call write_gvmus_nc (nout, gvmus)
deallocate (gvmus)
end if
if (write_gzvs) then
allocate (gzvs(ntubes,nztot,nvpa,nspec))
if (debug) write (*,*) 'stella_diagnostics::diagnose_stella::get_gzvs'
call get_gzvs (gnew, gzvs)
if (debug) write (*,*) 'stella_diagnostics::diagnose_stella::write_gzvs_nc'
if (proc0) call write_gzvs_nc (nout, gzvs)
deallocate (gzvs)
end if
! if (write_symmetry) then
! allocate (pflx_zvpa(nztot,nvpa,nspec))
! allocate (vflx_zvpa(nztot,nvpa,nspec))
! allocate (qflx_zvpa(nztot,nvpa,nspec))
! call get_fluxes_vs_zvpa (gnew, pflx_zvpa, vflx_zvpa, qflx_zvpa)
! deallocate (pflx_zvpa, vflx_zvpa, qflx_zvpa)
! end if
if (proc0) call sync_nc
deallocate (part_flux, mom_flux, heat_flux)
if(allocated(part_flux_x)) deallocate (part_flux_x)
if(allocated(mom_flux_x)) deallocate (mom_flux_x)
if(allocated(heat_flux_x)) deallocate (heat_flux_x)
deallocate(phi_out)
nout = nout + 1
end subroutine diagnose_stella
!==============================================
!=============== GET FLUXES ===================
!==============================================
! assumes that the non-Boltzmann part of df is passed in (aka h)
subroutine get_fluxes (g, pflx, vflx, qflx,&
pflx_vs_kxkyz, vflx_vs_kxkyz, qflx_vs_kxkyz)
use mp, only: sum_reduce
use constants, only: zi
use fields_arrays, only: phi, apar
use stella_layouts, only: kxkyz_lo
use stella_layouts, only: iky_idx, ikx_idx, iz_idx, it_idx, is_idx
use species, only: spec, nspec
use stella_geometry, only: jacob, grho, bmag, btor
use stella_geometry, only: gds21, gds22
use stella_geometry, only: geo_surf
use zgrid, only: delzed, nzgrid, ntubes
use vpamu_grids, only: nvpa, nmu
use vpamu_grids, only: vperp2, vpa
use run_parameters, only: fphi, fapar
use kt_grids, ONLY: nakx, naky
use kt_grids, only: aky, theta0
use gyro_averages, only: gyro_average, gyro_average_j1
implicit none
complex, dimension (:,:,kxkyz_lo%llim_proc:), intent (in) :: g
real, dimension (:), intent (out) :: pflx, vflx, qflx
real, dimension (:,:,-nzgrid:,:,:), intent (out) :: pflx_vs_kxkyz, vflx_vs_kxkyz, qflx_vs_kxkyz
integer :: ikxkyz, iky, ikx, iz, it, is, ia
real, dimension (:), allocatable :: flx_norm
real :: flx_norm_partial
complex, dimension (:,:), allocatable :: gtmp1, gtmp2, gtmp3
allocate (flx_norm(-nzgrid:nzgrid))
allocate (gtmp1(nvpa,nmu), gtmp2(nvpa,nmu), gtmp3(nvpa,nmu))
pflx = 0. ; vflx = 0. ; qflx = 0.
pflx_vs_kxkyz = 0. ; vflx_vs_kxkyz = 0. ; qflx_vs_kxkyz = 0.
flx_norm = jacob(1,:)*delzed
flx_norm(-nzgrid) = 0.5*flx_norm(-nzgrid)
flx_norm( nzgrid) = 0.5*flx_norm( nzgrid)
if (flux_norm) then
! Flux definition with an extra factor 1/<\nabla\rho> in front.
flx_norm_partial = sum(flx_norm)/sum(flx_norm*grho(1,:))
flx_norm = flx_norm/sum(flx_norm*grho(1,:))
else
! Flux definition witou the extra factor.
flx_norm_partial = 1.0
flx_norm = flx_norm/sum(flx_norm)
endif
ia = 1
! get electrostatic contributions to fluxes
if (fphi > epsilon(0.0)) then
do ikxkyz = kxkyz_lo%llim_proc, kxkyz_lo%ulim_proc
iky = iky_idx(kxkyz_lo,ikxkyz)
ikx = ikx_idx(kxkyz_lo,ikxkyz)
iz = iz_idx(kxkyz_lo,ikxkyz)
it = it_idx(kxkyz_lo,ikxkyz)
is = is_idx(kxkyz_lo,ikxkyz)
! get particle flux
call gyro_average (g(:,:,ikxkyz), ikxkyz, gtmp1)
call get_one_flux (iky, iz, flx_norm(iz), gtmp1, phi(iky,ikx,iz,it), pflx(is))
call get_one_flux (iky, iz, flx_norm_partial, gtmp1, phi(iky,ikx,iz,it), pflx_vs_kxkyz(iky,ikx,iz,it,is))
! get heat flux
! NEEDS TO BE MODIFIED TO TREAT ENERGY = ENERGY(ALPHA)
gtmp1 = gtmp1*(spread(vpa**2,2,nmu)+spread(vperp2(1,iz,:),1,nvpa))
call get_one_flux (iky, iz, flx_norm(iz), gtmp1, phi(iky,ikx,iz,it), qflx(is))
call get_one_flux (iky, iz, flx_norm_partial, gtmp1, phi(iky,ikx,iz,it), qflx_vs_kxkyz(iky,ikx,iz,it,is))
! get momentum flux
! parallel component
gtmp1 = g(:,:,ikxkyz)*spread(vpa,2,nmu)*geo_surf%rmaj*btor(iz)/bmag(ia,iz)
call gyro_average (gtmp1, ikxkyz, gtmp2)
gtmp1 = -g(:,:,ikxkyz)*zi*aky(iky)*spread(vperp2(ia,iz,:),1,nvpa)*geo_surf%rhoc &
* (gds21(ia,iz)+theta0(iky,ikx)*gds22(ia,iz))*spec(is)%smz &
/ (geo_surf%qinp*geo_surf%shat*bmag(ia,iz)**2)
call gyro_average_j1 (gtmp1, ikxkyz, gtmp3)
gtmp1 = gtmp2 + gtmp3
call get_one_flux (iky, iz, flx_norm(iz), gtmp1, phi(iky,ikx,iz,it), vflx(is))
call get_one_flux (iky, iz, flx_norm_partial, gtmp1, phi(iky,ikx,iz,it), vflx_vs_kxkyz(iky,ikx,iz,it,is))
end do
end if
if (fapar > epsilon(0.0)) then
! particle flux
do ikxkyz = kxkyz_lo%llim_proc, kxkyz_lo%ulim_proc
iky = iky_idx(kxkyz_lo,ikxkyz)
ikx = ikx_idx(kxkyz_lo,ikxkyz)
iz = iz_idx(kxkyz_lo,ikxkyz)
it = it_idx(kxkyz_lo,ikxkyz)
is = is_idx(kxkyz_lo,ikxkyz)
! Apar contribution to particle flux
gtmp1 = -g(:,:,ikxkyz)*spec(is)%stm*spread(vpa,2,nmu)
call gyro_average (gtmp1, ikxkyz, gtmp2)
call get_one_flux (iky, iz, flx_norm(iz), gtmp2, apar(iky,ikx,iz,it), pflx(is))
! Apar contribution to heat flux
gtmp2 = gtmp2*(spread(vpa**2,2,nmu)+spread(vperp2(ia,iz,:),1,nvpa))
call get_one_flux (iky, iz, flx_norm(iz), gtmp2, apar(iky,ikx,iz,it), qflx(is))
! Apar contribution to momentum flux
! parallel component
gtmp1 = -spread(vpa**2,2,nmu)*spec(is)%stm*g(:,:,ikxkyz) &
* geo_surf%rmaj*btor(iz)/bmag(1,iz)
call gyro_average (gtmp1, ikxkyz, gtmp2)
! perp component
gtmp1 = spread(vpa,2,nmu)*spec(is)%stm*g(:,:,ikxkyz) &
* zi*aky(iky)*spread(vperp2(ia,iz,:),1,nvpa)*geo_surf%rhoc &
* (gds21(ia,iz)+theta0(iky,ikx)*gds22(ia,iz))*spec(is)%smz &
/ (geo_surf%qinp*geo_surf%shat*bmag(ia,iz)**2)
call gyro_average_j1 (gtmp1, ikxkyz, gtmp3)
! FLAG -- NEED TO ADD IN CONTRIBUTION FROM BOLTZMANN PIECE !!
gtmp1 = gtmp2 + gtmp3
call get_one_flux (iky, iz, flx_norm(iz), gtmp1, apar(iky,ikx,iz,it), vflx(is))
end do
end if
call sum_reduce (pflx, 0) ; pflx = pflx*spec%dens_psi0
call sum_reduce (qflx, 0) ; qflx = qflx*spec%dens_psi0*spec%temp_psi0
call sum_reduce (vflx, 0) ; vflx = vflx*spec%dens_psi0*sqrt(spec%mass*spec%temp_psi0)
! normalise to account for contributions from multiple flux tubes
! in flux tube train
pflx = pflx/real(ntubes)
qflx = qflx/real(ntubes)
vflx = vflx/real(ntubes)
call sum_reduce (pflx_vs_kxkyz, 0)
call sum_reduce (qflx_vs_kxkyz, 0)
call sum_reduce (vflx_vs_kxkyz, 0)
do is = 1, nspec
pflx_vs_kxkyz(:,:,:,:,is) = pflx_vs_kxkyz(:,:,:,:,is)*spec(is)%dens_psi0
qflx_vs_kxkyz(:,:,:,:,is) = qflx_vs_kxkyz(:,:,:,:,is)*spec(is)%dens_psi0*spec(is)%temp_psi0
vflx_vs_kxkyz(:,:,:,:,is) = vflx_vs_kxkyz(:,:,:,:,is)*spec(is)%dens_psi0*sqrt(spec(is)%mass*spec(is)%temp_psi0)
enddo
deallocate (gtmp1, gtmp2, gtmp3)
deallocate (flx_norm)
end subroutine get_fluxes
!==============================================
!============ GET FLUXES VMULO ================
!==============================================
subroutine get_fluxes_vmulo (g, phi, pflx, vflx, qflx, pflx_x, vflx_x, qflx_x)
use mp, only: sum_reduce
use constants, only: zi
use dist_fn_arrays, only: g1, g2, kperp2, dkperp2dr
use stella_layouts, only: vmu_lo
use stella_layouts, only: iv_idx, imu_idx, is_idx
use species, only: spec
use stella_geometry, only: jacob, grho, bmag, btor
use stella_geometry, only: drhodpsi
use stella_geometry, only: gds21, gds22
use stella_geometry, only: dgds21dr, dgds22dr
use stella_geometry, only: geo_surf
use stella_geometry, only: dBdrho, dIdrho
use stella_geometry, only: dl_over_b, d_dl_over_b_drho
use zgrid, only: delzed, nzgrid, ntubes
use vpamu_grids, only: vperp2, vpa, mu
use vpamu_grids, only: maxwell_vpa, maxwell_mu, maxwell_fac
use run_parameters, only: fphi
use kt_grids, only: aky, theta0, naky, nakx, multiply_by_rho
use physics_flags, only: radial_variation
use gyro_averages, only: gyro_average, gyro_average_j1, aj0x, aj1x
implicit none
complex, dimension (:,:,-nzgrid:,:,vmu_lo%llim_proc:), intent (in) :: g
complex, dimension (:,:,-nzgrid:,:), intent (in) :: phi
real, dimension (:), intent (out) :: pflx, vflx, qflx
real, dimension (:,:), intent (out) :: pflx_x, vflx_x, qflx_x
integer :: ivmu, imu, iv, iz, it, is, ia
real, dimension (:), allocatable :: flx_norm, dflx_norm
complex, dimension (:,:), allocatable :: g0k,g1k
allocate (flx_norm(-nzgrid:nzgrid))
allocate (dflx_norm(-nzgrid:nzgrid))
pflx = 0. ; vflx = 0. ; qflx = 0.
pflx_x = 0. ; vflx_x = 0. ; qflx_x = 0.
ia = 1
flx_norm = jacob(1,:)*delzed
flx_norm(-nzgrid) = 0.5*flx_norm(-nzgrid)
flx_norm( nzgrid) = 0.5*flx_norm( nzgrid)
flx_norm = flx_norm/sum(flx_norm*grho(1,:))
allocate (g0k(naky,nakx))
allocate (g1k(naky,nakx))
if(radial_variation) then
where (dl_over_b(ia,:) .gt. epsilon(0.0))
dflx_norm = d_dl_over_b_drho(ia,:)/dl_over_b(ia,:)
elsewhere
dflx_norm = 0.
endwhere
endif
! FLAG - electrostatic for now
! get electrostatic contributions to fluxes
if (fphi > epsilon(0.0)) then
ia = 1
!get particle flux
do ivmu = vmu_lo%llim_proc, vmu_lo%ulim_proc
iv = iv_idx(vmu_lo,ivmu)
imu = imu_idx(vmu_lo,ivmu)
is = is_idx(vmu_lo,ivmu)
call gyro_average (g(:,:,:,:,ivmu), ivmu, g1(:,:,:,:,ivmu))
do it = 1, ntubes
do iz= -nzgrid, nzgrid
if(radial_variation) then
g0k = g1(:,:,iz,it,ivmu) &
* (-0.5*aj1x(:,:,iz,ivmu)/aj0x(:,:,iz,ivmu)*(spec(is)%smz)**2 &
* (kperp2(:,:,ia,iz)*vperp2(ia,iz,imu)/bmag(ia,iz)**2) &
* (dkperp2dr(:,:,ia,iz) - dBdrho(iz)/bmag(ia,iz)) &
+ dBdrho(iz)/bmag(ia,iz) + dflx_norm(iz))
call multiply_by_rho(g0k)
g1(:,:,iz,it,ivmu) = g1(:,:,iz,it,ivmu) + g0k
endif
!subtract adiabatic contribution part of g
g0k = spec(is)%zt*fphi*phi(:,:,iz,it)*aj0x(:,:,iz,ivmu)**2 &
*maxwell_vpa(iv,is)*maxwell_mu(ia,iz,imu,is)*maxwell_fac(is)
if(radial_variation) then
g1k = g0k*( -spec(is)%tprim*(vpa(iv)**2+vperp2(ia,iz,imu) - 2.5) &
-spec(is)%fprim - 2.0*dBdrho(iz)*mu(imu) &
-aj1x(:,:,iz,ivmu)/aj0x(:,:,iz,ivmu)*(spec(is)%smz)**2 &
* (kperp2(:,:,ia,iz)*vperp2(ia,iz,imu)/bmag(ia,iz)**2) &
* (dkperp2dr(:,:,ia,iz) - dBdrho(iz)/bmag(ia,iz)) &
+ dBdrho(iz)/bmag(ia,iz)+dflx_norm(iz))
call multiply_by_rho(g1k)
g0k = g0k + g1k
endif
g1(:,:,iz,it,ivmu) = g1(:,:,iz,it,ivmu) + g0k
enddo
enddo
enddo
call get_one_flux_vmulo (flx_norm, spec%dens_psi0, g1, phi, pflx)
if(write_radial_fluxes) then
call get_one_flux_radial (flx_norm, spec%dens_psi0, g1, phi, pflx_x)
endif
!get heat flux
do ivmu = vmu_lo%llim_proc, vmu_lo%ulim_proc
iv = iv_idx(vmu_lo,ivmu)
imu = imu_idx(vmu_lo,ivmu)
is = is_idx(vmu_lo,ivmu)
call gyro_average (g(:,:,:,:,ivmu), ivmu, g1(:,:,:,:,ivmu))
do it = 1, ntubes
do iz= -nzgrid, nzgrid
g1(:,:,iz,it,ivmu) = g1(:,:,iz,it,ivmu)*(vpa(iv)**2+vperp2(ia,iz,imu))
if(radial_variation) then
g0k = g1(:,:,iz,it,ivmu) &
* (-0.5*aj1x(:,:,iz,ivmu)/aj0x(:,:,iz,ivmu)*(spec(is)%smz)**2 &
* (kperp2(:,:,ia,iz)*vperp2(ia,iz,imu)/bmag(ia,iz)**2) &
* (dkperp2dr(:,:,ia,iz) - dBdrho(iz)/bmag(ia,iz)) &
+ dBdrho(iz)/bmag(ia,iz) &
+ 2.0*mu(imu)*dBdrho(iz)/(vpa(iv)**2+vperp2(ia,iz,imu)) &
+ dflx_norm(iz))
call multiply_by_rho(g0k)
g1(:,:,iz,it,ivmu) = g1(:,:,iz,it,ivmu) + g0k
endif
!subtract adiabatic contribution part of g
g0k = spec(is)%zt*fphi*phi(:,:,iz,it)*aj0x(:,:,iz,ivmu)**2 &
*(vpa(iv)**2+vperp2(ia,iz,imu)) &
*maxwell_vpa(iv,is)*maxwell_mu(ia,iz,imu,is)*maxwell_fac(is)
if(radial_variation) then
g1k = g0k*( -spec(is)%tprim*(vpa(iv)**2+vperp2(ia,iz,imu) - 2.5) &
-spec(is)%fprim - 2.0*dBdrho(iz)*mu(imu) &
-aj1x(:,:,iz,ivmu)/aj0x(:,:,iz,ivmu)*(spec(is)%smz)**2 &
* (kperp2(:,:,ia,iz)*vperp2(ia,iz,imu)/bmag(ia,iz)**2) &
* (dkperp2dr(:,:,ia,iz) - dBdrho(iz)/bmag(ia,iz)) &
+ dBdrho(iz)/bmag(ia,iz)+dflx_norm(iz) &
+ 2.0*mu(imu)*dBdrho(iz)/(vpa(iv)**2+vperp2(ia,iz,imu)))
call multiply_by_rho(g1k)
g0k = g0k + g1k
endif
g1(:,:,iz,it,ivmu) = g1(:,:,iz,it,ivmu) + g0k
enddo
enddo
enddo
call get_one_flux_vmulo (flx_norm,spec%dens_psi0*spec%temp_psi0, g1, phi, qflx)
if(write_radial_fluxes) then
call get_one_flux_radial (flx_norm,spec%dens_psi0*spec%temp_psi0, g1, phi, qflx_x)
endif
! get momentum flux
do ivmu = vmu_lo%llim_proc, vmu_lo%ulim_proc
iv = iv_idx(vmu_lo,ivmu)
imu = imu_idx(vmu_lo,ivmu)
is = is_idx(vmu_lo,ivmu)
do it = 1, ntubes
do iz= -nzgrid, nzgrid
! parallel component
g0k = g(:,:,iz,it,ivmu)*vpa(iv)*geo_surf%rmaj*btor(iz)/bmag(ia,iz)
call gyro_average (g0k, iz, ivmu, g1(:,:,iz,it,ivmu))
if(radial_variation) then
g0k = g1(:,:,iz,it,ivmu) &
* (-0.5*aj1x(:,:,iz,ivmu)/aj0x(:,:,iz,ivmu)*(spec(is)%smz)**2 &
* (kperp2(:,:,ia,iz)*vperp2(ia,iz,imu)/bmag(ia,iz)**2) &
* (dkperp2dr(:,:,ia,iz) - dBdrho(iz)/bmag(ia,iz)) &
+ dIdrho/(geo_surf%rmaj*btor(iz)) &
+ dflx_norm(iz))
call multiply_by_rho(g0k)
g1(:,:,iz,it,ivmu) = g1(:,:,iz,it,ivmu) + g0k
endif
!subtract adiabatic contribution part of g
g0k = spec(is)%zt*fphi*phi(:,:,iz,it)*aj0x(:,:,iz,ivmu)**2 &
*vpa(iv)*geo_surf%rmaj*btor(iz)/bmag(ia,iz) &
*maxwell_vpa(iv,is)*maxwell_mu(ia,iz,imu,is)*maxwell_fac(is)
if(radial_variation) then
g1k = g0k*( -spec(is)%tprim*(vpa(iv)**2+vperp2(ia,iz,imu) - 2.5) &
-spec(is)%fprim - 2.0*dBdrho(iz)*mu(imu) &
-aj1x(:,:,iz,ivmu)/aj0x(:,:,iz,ivmu)*(spec(is)%smz)**2 &
* (kperp2(:,:,ia,iz)*vperp2(ia,iz,imu)/bmag(ia,iz)**2) &
* (dkperp2dr(:,:,ia,iz) - dBdrho(iz)/bmag(ia,iz)) &
+ dflx_norm(iz) + dIdrho/(geo_surf%rmaj*btor(iz)))
call multiply_by_rho(g1k)
g0k = g0k + g1k
endif
g1(:,:,iz,it,ivmu) = g1(:,:,iz,it,ivmu) + g0k
! perpendicular component
g0k = -g(:,:,iz,it,ivmu)*zi*spread(aky,2,nakx)*vperp2(ia,iz,imu)*geo_surf%rhoc &
* (gds21(ia,iz)+theta0*gds22(ia,iz))*spec(is)%smz &
/ (geo_surf%qinp*geo_surf%shat*bmag(ia,iz)**2)
call gyro_average_j1 (g0k, iz, ivmu, g2(:,:,iz,it,ivmu))
if(radial_variation) then
g0k = - g(:,:,iz,it,ivmu)*zi*spread(aky,2,nakx)*vperp2(ia,iz,imu)*geo_surf%rhoc &
* (dgds21dr(ia,iz)+theta0*dgds22dr(ia,iz))*aj1x(:,:,iz,ivmu)*spec(is)%smz &
/ (geo_surf%qinp*geo_surf%shat*bmag(ia,iz)**2)
g0k = g0k - g(:,:,iz,it,ivmu)*zi*spread(aky,2,nakx)*vperp2(ia,iz,imu)*geo_surf%rhoc &
* (gds21(ia,iz)+theta0*gds22(ia,iz))*spec(is)%smz &
/ (geo_surf%qinp*geo_surf%shat*bmag(ia,iz)**2) &
* (0.5*aj0x(:,:,iz,ivmu) - aj1x(:,:,iz,ivmu)) &
* (dkperp2dr(:,:,ia,iz) - dBdrho(iz)/bmag(ia,iz))
g0k = g0k + g2(:,:,iz,it,ivmu) &
* (- geo_surf%d2qdr2*geo_surf%rhoc/(geo_surf%shat*geo_surf%qinp) &
- geo_surf%d2psidr2*drhodpsi + dflx_norm(iz))
call multiply_by_rho(g0k)
g2(:,:,iz,it,ivmu) = g2(:,:,iz,it,ivmu) + g0k
endif
!subtract adiabatic contribution part of g
g0k = -spec(is)%zt*fphi*phi(:,:,iz,it)*aj0x(:,:,iz,ivmu)*aj1x(:,:,iz,ivmu) &
*maxwell_vpa(iv,is)*maxwell_mu(ia,iz,imu,is)*maxwell_fac(is) &
*zi*spread(aky,2,nakx)*vperp2(ia,iz,imu)*geo_surf%rhoc &
* (gds21(ia,iz)+theta0*gds22(ia,iz))*spec(is)%smz &
/ (geo_surf%qinp*geo_surf%shat*bmag(ia,iz)**2)
if(radial_variation) then
g1k = -spec(is)%zt*fphi*phi(:,:,iz,it)*aj0x(:,:,iz,ivmu)*aj1x(:,:,iz,ivmu) &
*maxwell_vpa(iv,is)*maxwell_mu(ia,iz,imu,is)*maxwell_fac(is) &
*zi*spread(aky,2,nakx)*vperp2(ia,iz,imu)*geo_surf%rhoc &
* (dgds21dr(ia,iz)+theta0*dgds22dr(ia,iz))*spec(is)%smz &
/ (geo_surf%qinp*geo_surf%shat*bmag(ia,iz)**2)
g1k = g1k -spec(is)%zt*fphi*phi(:,:,iz,it)*aj0x(:,:,iz,ivmu) &
*maxwell_vpa(iv,is)*maxwell_mu(ia,iz,imu,is)*maxwell_fac(is) &
*zi*spread(aky,2,nakx)*vperp2(ia,iz,imu)*geo_surf%rhoc &
* (gds21(ia,iz)+theta0*gds22(ia,iz))*spec(is)%smz &
/ (geo_surf%qinp*geo_surf%shat*bmag(ia,iz)**2) &
* (0.5*aj0x(:,:,iz,ivmu) - aj1x(:,:,iz,ivmu)) &
* (dkperp2dr(:,:,ia,iz) - dBdrho(iz)/bmag(ia,iz))
g1k = g1k + &
g0k*(-spec(is)%tprim*(vpa(iv)**2+vperp2(ia,iz,imu) - 2.5) &
-spec(is)%fprim - 2.0*dBdrho(iz)*mu(imu) &
-0.5*aj1x(:,:,iz,ivmu)/aj0x(:,:,iz,ivmu)*(spec(is)%smz)**2 &
* (kperp2(:,:,ia,iz)*vperp2(ia,iz,imu)/bmag(ia,iz)**2) &
* (dkperp2dr(:,:,ia,iz) - dBdrho(iz)/bmag(ia,iz)) &
- geo_surf%d2qdr2*geo_surf%rhoc/(geo_surf%shat*geo_surf%qinp) &
- geo_surf%d2psidr2*drhodpsi + dflx_norm(iz))
call multiply_by_rho(g1k)
g0k = g0k + g1k
endif
g2(:,:,iz,it,ivmu) = g2(:,:,iz,it,ivmu) + g0k
enddo
enddo
enddo
g1 = g1 + g2
call get_one_flux_vmulo (flx_norm,spec%dens_psi0*sqrt(spec%mass*spec%temp_psi0), g1, phi, vflx)
if(write_radial_fluxes) then
call get_one_flux_radial (flx_norm, spec%dens_psi0*sqrt(spec%mass*spec%temp_psi0), g1, phi, vflx_x)
endif
end if
! normalise to account for contributions from multiple flux tubes
! in flux tube train
pflx = pflx/real(ntubes)
qflx = qflx/real(ntubes)
vflx = vflx/real(ntubes)
deallocate (flx_norm)
if(allocated(dflx_norm)) deallocate(dflx_norm)
if(allocated(g0k)) deallocate(g0k)
if(allocated(g1k)) deallocate(g1k)
end subroutine get_fluxes_vmulo
!==============================================
!=============== GET ONE FLUX =================
!==============================================
subroutine get_one_flux (iky, iz, norm, gin, fld, flxout)
use vpamu_grids, only: integrate_vmu
use volume_averages, only: mode_fac
use kt_grids, only: aky
implicit none
integer, intent (in) :: iky, iz
real, intent (in) :: norm
complex, dimension (:,:), intent (in) :: gin
complex, intent (in) :: fld
real, intent (in out) :: flxout
complex :: flx
call integrate_vmu (gin,iz,flx)
flxout = flxout &
+ 0.5*mode_fac(iky)*aky(iky)*aimag(flx*conjg(fld))*norm
end subroutine get_one_flux
!==============================================
!============ GET ONE FLUX VMULO ==============
!==============================================
subroutine get_one_flux_vmulo (norm, weights, gin, fld, flxout)
use vpamu_grids, only: integrate_vmu
use stella_layouts, only: vmu_lo
use kt_grids, only: aky, nakx, naky
use zgrid, only: nzgrid, ntubes
use species, only: nspec
use volume_averages, only: mode_fac
implicit none
real, dimension (-nzgrid:), intent (in) :: norm
real, dimension (:), intent (in) :: weights
complex, dimension (:,:,-nzgrid:,:,vmu_lo%llim_proc:), intent (in) :: gin
complex, dimension (:,:,-nzgrid:,:), intent (in) :: fld
real, dimension (:), intent (in out) :: flxout
complex, dimension (:,:,:,:,:), allocatable :: totals
integer :: is, it, iz, ikx
allocate (totals(naky,nakx,-nzgrid:nzgrid,ntubes,nspec))
call integrate_vmu (gin,weights,totals)
do is = 1, nspec
do it = 1, ntubes
do iz= -nzgrid, nzgrid
do ikx = 1, nakx
flxout(is) = flxout(is) &
+ sum(0.5*mode_fac*aky*aimag(totals(:,ikx,iz,it,is)*conjg(fld(:,ikx,iz,it)))*norm(iz))
enddo
enddo
enddo
enddo
deallocate (totals)
end subroutine get_one_flux_vmulo
!==============================================
!=========== GET ONE FLUX RADIAL ==============
!==============================================
subroutine get_one_flux_radial (norm, weights, gin, fld, flxout)
use vpamu_grids, only: integrate_vmu
use stella_layouts, only: vmu_lo
use kt_grids, only: aky, nakx, naky
use zgrid, only: nzgrid, ntubes
use species, only: nspec
use volume_averages, only: mode_fac
use stella_transforms, only: transform_kx2x_unpadded